mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT > Struct Template Reference

template<typename _realT, typename aosysT>
struct mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >

Class to manage the calculation of temporal PSDs of the Fourier modes in atmospheric turbulence.

Works with both basic (sines/cosines) and modified Fourier modes.

Template Parameters
realTis a real floating point type for calculations. Currently must be double due to gsl_integration.
aosysTis an AO system type, usually of type ao_system.
Todo:

Split off the integration parameters in a separate structure.

once integration parameters are in a separate structure, make this a class with protected members.

GSL error handling

Definition at line 101 of file zernikeTemporalPSD.hpp.

#include <ao/analysis/zernikeTemporalPSD.hpp>

Public Member Functions

 zernikeTemporalPSD ()
 Default c'tor. More...
 
 zernikeTemporalPSD (bool allocate)
 Constructor with workspace allocation. More...
 
 ~zernikeTemporalPSD ()
 Destructor. More...
 
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 Fourier mode for a single layer. More...
 
template<bool parallel = true>
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. More...
 
GSL Integration Tolerances

For good results it seems that absolute tolerance (absTol) needs to be 1e-10. Lower tolerances cause some frequencies to drop out, etc. Relative tolerance (relTol) seems to be less sensitive, and 1e-4 works on cases tested as of 1 Jan, 2017.

See the documentation for the GSL Library integrators at (https://www.gnu.org/software/gsl/manual/htmlm_node/QAGI-adaptive-integration-on-infinite-intervals.html)

void absTol (realT at)
 Set absolute tolerance. More...
 
realT absTol ()
 Get the current absolute tolerance. More...
 
void relTol (realT rt)
 Set relative tolerance. More...
 
realT relTol ()
 Get the current relative tolerance. More...
 

Public Attributes

aosysT * m_aosys {nullptr}
 Pointer to an AO system structure. More...
 
realT m_f {0}
 the current temporal frequency0 More...
 
realT m_zern_j {0}
 the current mode number More...
 
int m_zern_m {0}
 The current mode m. More...
 
int m_zern_n {0}
 The current mode n. More...
 
realT m_cq {0}
 The cosine of the wind direction. More...
 
realT m_sq {0}
 The sine of the wind direction. More...
 
realT m_spatialFilter {false}
 Flag indicating if a spatial filter is applied. More...
 
int _layer_i
 The index of the current layer. More...
 
gsl_integration_workspace * _w
 Worskspace for the gsl integrators, allocated to WSZ if constructed as worker (with allocate == true). More...
 
realT _absTol
 The absolute tolerance to use in the GSL integrator. More...
 
realT _relTol
 The relative tolerance to use in the GSL integrator. More...
 

Protected Member Functions

void initialize ()
 Initialize parameters to default values. More...
 

Constructor & Destructor Documentation

◆ zernikeTemporalPSD() [1/2]

template<typename realT , typename aosysT >
mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::zernikeTemporalPSD

Default c'tor.

Definition at line 252 of file zernikeTemporalPSD.hpp.

◆ zernikeTemporalPSD() [2/2]

template<typename realT , typename aosysT >
mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::zernikeTemporalPSD ( bool  allocate)
explicit

Constructor with workspace allocation.

Parameters
allocateif true, then the workspace for GSL integrators is allocated.

Definition at line 259 of file zernikeTemporalPSD.hpp.

◆ ~zernikeTemporalPSD()

template<typename realT , typename aosysT >
mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::~zernikeTemporalPSD

Destructor.

Frees GSL workspace if it was allocated.

Definition at line 271 of file zernikeTemporalPSD.hpp.

Member Function Documentation

◆ absTol() [1/2]

template<typename realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::absTol

Get the current absolute tolerance.

Returns
_absTol

Definition at line 295 of file zernikeTemporalPSD.hpp.

◆ absTol() [2/2]

template<typename realT , typename aosysT >
void mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::absTol ( realT  at)

Set absolute tolerance.

Parameters
atis the new absolute tolerance.

Definition at line 289 of file zernikeTemporalPSD.hpp.

◆ initialize()

template<typename realT , typename aosysT >
void mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::initialize
protected

Initialize parameters to default values.

Definition at line 280 of file zernikeTemporalPSD.hpp.

◆ multiLayerPSD()

template<typename realT , typename aosysT >
template<bool parallel>
int mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::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.

Template Parameters
parallelcontrols whether layers are calculated in parallel. Default is true. Set to false if this is called inside a parallelized loop, as in makePSDGrid.
Todo:

implement error checking

handle reports of convergence failures form singleLayerPSD when implemented.

Parameters
[out]PSDthe calculated PSD
[in]freqthe populated temporal frequency grid defining at which frequencies the PSD is calculated
[in]fmax[optional] set the maximum temporal frequency for the calculation. The PSD is filled in with a -17/3 power law past this frequency. If 0, then it is taken to be 150 Hz + 2*fastestPeak(m,n).

Definition at line 464 of file zernikeTemporalPSD.hpp.

◆ relTol() [1/2]

template<typename realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::relTol

Get the current relative tolerance.

Returns
_relTol

Definition at line 307 of file zernikeTemporalPSD.hpp.

◆ relTol() [2/2]

template<typename realT , typename aosysT >
void mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::relTol ( realT  rt)

Set relative tolerance.

Parameters
rtis the new relative tolerance.

Definition at line 301 of file zernikeTemporalPSD.hpp.

◆ singleLayerPSD()

template<typename realT , typename aosysT >
int mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::singleLayerPSD ( std::vector< realT > &  PSD,
std::vector< realT > &  freq,
int  zern_j,
int  layer_i,
realT  fmax = 0 
)

Calculate the temporal PSD for a Fourier mode for a single layer.

Todo:

implement error checking.

need a way to track convergence failures in integral without throwing an error.

need better handling of averaging for the -17/3 extension.

Parameters
[out]PSDthe calculated PSD
[in]freqthe populated temporal frequency grid defining the frequencies at which the PSD is calculated
[in]zern_j
[in]layer_ithe index of the layer, for accessing the atmosphere parameters
[in]fmax[optional] set the maximum temporal frequency for the calculation. The PSD is filled in with a -17/3 power law past this frequency.

Definition at line 314 of file zernikeTemporalPSD.hpp.

References mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_layer_i, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_w, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_aosys, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_cq, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_f, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_spatialFilter, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_sq, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_j, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_m, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_n, and mx::sigproc::noll_nm().

Member Data Documentation

◆ _absTol

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_absTol

The absolute tolerance to use in the GSL integrator.

Definition at line 123 of file zernikeTemporalPSD.hpp.

◆ _layer_i

template<typename _realT , typename aosysT >
int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_layer_i

The index of the current layer.

Definition at line 118 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ _relTol

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_relTol

The relative tolerance to use in the GSL integrator.

Definition at line 124 of file zernikeTemporalPSD.hpp.

◆ _w

template<typename _realT , typename aosysT >
gsl_integration_workspace* mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_w

Worskspace for the gsl integrators, allocated to WSZ if constructed as worker (with allocate == true).

Definition at line 121 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_aosys

template<typename _realT , typename aosysT >
aosysT* mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_aosys {nullptr}

Pointer to an AO system structure.

Definition at line 107 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_cq

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_cq {0}

The cosine of the wind direction.

Definition at line 114 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_f

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_f {0}

the current temporal frequency0

Definition at line 109 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_spatialFilter

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_spatialFilter {false}

Flag indicating if a spatial filter is applied.

Definition at line 116 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_sq

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_sq {0}

The sine of the wind direction.

Definition at line 115 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_zern_j

template<typename _realT , typename aosysT >
realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_j {0}

the current mode number

Definition at line 110 of file zernikeTemporalPSD.hpp.

Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().

◆ m_zern_m

template<typename _realT , typename aosysT >
int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_m {0}

◆ m_zern_n

template<typename _realT , typename aosysT >
int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_n {0}

The documentation for this struct was generated from the following file: