27#ifndef zernikeTemporalPSD_hpp
28#define zernikeTemporalPSD_hpp
36#include <gsl/gsl_integration.h>
37#include <gsl/gsl_errno.h>
41#include "../../math/constants.hpp"
42#include "../../math/func/jinc.hpp"
43#include "../../math/func/airyPattern.hpp"
44#include "../../math/vectorUtils.hpp"
45#include "../../ioutils/fits/fitsFile.hpp"
46#include "../../sigproc/zernike.hpp"
47#include "../../ioutils/stringUtils.hpp"
48#include "../../ioutils/readColumns.hpp"
49#include "../../ioutils/binVector.hpp"
50#include "../../ioutils/fileUtils.hpp"
52#include "../../ipc/ompLoopWatcher.hpp"
53#include "../../mxError.hpp"
82template <
typename realT,
typename aosysT>
97template <
typename _realT,
typename aosysT>
100 typedef _realT realT;
101 typedef std::complex<realT> complexT;
107 realT m_apertureRatio{ 1 };
113 m_apertureRatio = ar;
117 realT apertureRatio4()
134 gsl_integration_workspace *
_w;
139 std::vector<realT> Jps;
140 std::vector<realT> Jms;
142 std::vector<realT> ms;
143 std::vector<realT> ns;
210 std::vector<realT> &freq,
223 template <
bool m_parallel>
230 std::vector<realT> &PSD, std::vector<realT> &freq,
int zern_j, realT fmax, isParallel<true> parallel );
234 std::vector<realT> &PSD, std::vector<realT> &freq,
int zern_j, realT fmax, isParallel<false> parallel );
249 template <
bool parallel = true>
251 std::vector<realT> &PSD,
261template <
typename realT,
typename aosysT>
268template <
typename realT,
typename aosysT>
276 _w = gsl_integration_workspace_alloc( WSZ );
280template <
typename realT,
typename aosysT>
285 gsl_integration_workspace_free( _w );
289template <
typename realT,
typename aosysT>
298template <
typename realT,
typename aosysT>
304template <
typename realT,
typename aosysT>
310template <
typename realT,
typename aosysT>
316template <
typename realT,
typename aosysT>
322template <
typename realT,
typename aosysT>
324 std::vector<realT> &PSD, std::vector<realT> &freq,
int zern_j,
int layer_i, realT fmax )
327 fmax = freq[freq.size() - 1];
329 realT v_wind = m_aosys->atm.layer_v_wind( layer_i );
330 realT q_wind = m_aosys->atm.layer_dir( layer_i );
333 realT cq = cos( q_wind );
334 realT sq = sin( q_wind );
336 realT scale = 2 * ( 1 / v_wind );
339 gsl_set_error_handler_off();
345 params.m_apertureRatio = m_apertureRatio;
352 if( m_aosys->spatialFilter_ku() < std::numeric_limits<realT>::max() ||
353 m_aosys->spatialFilter_kv() < std::numeric_limits<realT>::max() )
360 func.function = &F_zernike<realT, aosysT>;
362 func.params = ¶ms;
366 while( freq[i] <= fmax )
368 params.
m_f = freq[i];
370 int ec = gsl_integration_qagi( &func, _absTol, _relTol, WSZ, params.
_w, &result, &error );
372 if( ec == GSL_EDIVERGE )
374 std::cerr <<
"GSL_EDIVERGE:" <<
" " << freq[i] <<
" " << v_wind <<
" " << zern_j <<
"\n";
375 std::cerr <<
"ignoring . . .\n";
378 PSD[i] = scale * result;
381 if( i >= freq.size() )
410template <
typename realT,
typename aosysT>
412 std::vector<realT> &PSD, std::vector<realT> &freq,
int zern_j, realT fmax, isParallel<true> parallel )
414 static_cast<void>( parallel );
419 std::vector<realT> single_PSD( freq.size() );
422 for(
size_t i = 0; i < m_aosys->atm.n_layers(); ++i )
424 singleLayerPSD( single_PSD, freq, zern_j, i, fmax );
428 for(
size_t j = 0; j < freq.size(); ++j )
430 PSD[j] += m_aosys->atm.layer_Cn2( i ) * single_PSD[j];
438template <
typename realT,
typename aosysT>
439int zernikeTemporalPSD<realT, aosysT>::m_multiLayerPSD(
440 std::vector<realT> &PSD, std::vector<realT> &freq,
int zern_j, realT fmax, isParallel<false> parallel )
442 static_cast<void>( parallel );
445 std::vector<realT> single_PSD( freq.size() );
447 for(
size_t i = 0; i < m_aosys->atm.n_layers(); ++i )
449 singleLayerPSD( single_PSD, freq, zern_j, i, fmax );
452 for(
size_t j = 0; j < freq.size(); ++j )
454 PSD[j] += m_aosys->atm.layer_Cn2( i ) * single_PSD[j];
461template <
typename realT,
typename aosysT>
462template <
bool parallel>
464 std::vector<realT> &freq,
469 for(
size_t j = 0; j < PSD.size(); ++j )
477 return m_multiLayerPSD( PSD, freq, zern_j, fmax, isParallel<parallel>() );
483template <
typename realT,
typename aosysT>
492 realT apertureRatio4 = Fp->apertureRatio4();
494 realT ku = f / v_wind;
496 realT phi = atan( kv / ku );
498 realT k = sqrt( pow( ku, 2 ) + pow( kv, 2 ) );
Calculate and provide constants related to adaptive optics.
Spatial power spectra used in adaptive optics.
Declares and defines an analytical AO system.
Provides a class to manage closed loop gain linear predictor determination.
Provides a class to manage closed loop gain optimization.
realT F_zernike(realT kv, void *params)
Worker function for GSL Integration for the Zernike modes.
int noll_nm(int &n, int &m, int j)
Get the Zernike coefficients n,m corrresponding the Noll index j.
realT zernikeQNorm(realT k, realT phi, int n, int m)
Calculate the square-normed Fourier transform of a Zernike polynomial at position (k,...
Calculates the PSD of speckle intensity given a modified Fourier mode amplitude PSD.
Class to manage the calculation of temporal PSDs of the Fourier modes in atmospheric turbulence.
void apertureRatio(realT ar)
Aperture ratio ^4 used for calculations.
realT m_apertureRatio4
Aperture ratio, < 1, used for segments.
int singleLayerPSD(std::vector< realT > &PSD, std::vector< realT > &freq, int zern_j, int layer_i, realT fmax=0)
Calculate the temporal PSD for a Zernike mode for a single layer.
realT m_spatialFilter
Flag indicating if a spatial filter is applied.
realT m_zern_j
the current mode number
void initialize()
Initialize parameters to default values.
int m_zern_m
The current mode m.
gsl_integration_workspace * _w
Worskspace for the gsl integrators, allocated to WSZ if constructed as worker (with allocate == true)...
realT _relTol
The relative tolerance to use in the GSL integrator.
realT m_cq
The cosine of the wind direction.
int m_zern_n
The current mode n.
aosysT * m_aosys
Pointer to an AO system structure.
int _layer_i
The index of the current layer.
zernikeTemporalPSD()
Default c'tor.
realT m_f
the current temporal frequency
realT relTol()
Get the current relative tolerance.
~zernikeTemporalPSD()
Destructor.
int multiLayerPSD(std::vector< realT > &PSD, std::vector< realT > &freq, int zern_j, realT fmax=0)
Calculate the temporal PSD for a Fourier mode in a multi-layer model.
realT m_sq
The sine of the wind direction.
realT absTol()
Get the current absolute tolerance.
realT _absTol
The absolute tolerance to use in the GSL integrator.
A utility to convert a wavefront variance map to an intensity image.
Declares and defines a function to calculate the measurement noise PSD.