|
mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
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.
| realT | is a real floating point type for calculations. Currently must be double due to gsl_integration. |
| aosysT | is an AO system type, usually of type ao_system. |
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 98 of file zernikeTemporalPSD.hpp.
#include <ao/analysis/zernikeTemporalPSD.hpp>
Public Member Functions | |
| void | apertureRatio (realT ar) |
| Aperture ratio ^4 used for calculations. | |
| zernikeTemporalPSD () | |
| Default c'tor. | |
| zernikeTemporalPSD (bool allocate) | |
| Constructor with workspace allocation. | |
| ~zernikeTemporalPSD () | |
| Destructor. | |
| 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. | |
| 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. | |
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. | |
| realT | absTol () |
| Get the current absolute tolerance. | |
| void | relTol (realT rt) |
| Set relative tolerance. | |
| realT | relTol () |
| Get the current relative tolerance. | |
Public Attributes | |
| aosysT * | m_aosys { nullptr } |
| Pointer to an AO system structure. | |
| realT | m_f { 0 } |
| the current temporal frequency | |
| realT | m_zern_j { 0 } |
| the current mode number | |
| int | m_zern_m { 0 } |
| The current mode m. | |
| int | m_zern_n { 0 } |
| The current mode n. | |
| realT | m_cq { 0 } |
| The cosine of the wind direction. | |
| realT | m_sq { 0 } |
| The sine of the wind direction. | |
| realT | m_spatialFilter { false } |
| Flag indicating if a spatial filter is applied. | |
| int | _layer_i |
| The index of the current layer. | |
| gsl_integration_workspace * | _w |
| Worskspace for the gsl integrators, allocated to WSZ if constructed as worker (with allocate == true). | |
| realT | _absTol |
| The absolute tolerance to use in the GSL integrator. | |
| realT | _relTol |
| The relative tolerance to use in the GSL integrator. | |
Protected Member Functions | |
| void | initialize () |
| Initialize parameters to default values. | |
Protected Attributes | |
| realT | m_apertureRatio4 { 1 } |
| Aperture ratio, < 1, used for segments. | |
| mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::zernikeTemporalPSD | ( | ) |
Default c'tor.
Definition at line 262 of file zernikeTemporalPSD.hpp.
|
explicit |
Constructor with workspace allocation.
| allocate | if true, then the workspace for GSL integrators is allocated. |
Definition at line 269 of file zernikeTemporalPSD.hpp.
References WSZ.
| mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::~zernikeTemporalPSD | ( | ) |
Destructor.
Frees GSL workspace if it was allocated.
Definition at line 281 of file zernikeTemporalPSD.hpp.
| realT mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::absTol | ( | ) |
Get the current absolute tolerance.
Definition at line 305 of file zernikeTemporalPSD.hpp.
| void mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::absTol | ( | realT | at | ) |
Set absolute tolerance.
| at | is the new absolute tolerance. |
Definition at line 299 of file zernikeTemporalPSD.hpp.
|
inline |
Aperture ratio ^4 used for calculations.
Definition at line 111 of file zernikeTemporalPSD.hpp.
References mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_apertureRatio4.
|
protected |
Initialize parameters to default values.
Definition at line 290 of file zernikeTemporalPSD.hpp.
| 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.
| parallel | controls whether layers are calculated in parallel. Default is true. Set to false if this is called inside a parallelized loop, as in makePSDGrid. |
implement error checking
handle reports of convergence failures form singleLayerPSD when implemented.
| [out] | PSD | the calculated PSD |
| [in] | freq | the 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 463 of file zernikeTemporalPSD.hpp.
| realT mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::relTol | ( | ) |
Get the current relative tolerance.
Definition at line 317 of file zernikeTemporalPSD.hpp.
| void mx::AO::analysis::zernikeTemporalPSD< realT, aosysT >::relTol | ( | realT | rt | ) |
Set relative tolerance.
| rt | is the new relative tolerance. |
Definition at line 311 of file zernikeTemporalPSD.hpp.
| 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 Zernike mode for a single layer.
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.
| [out] | PSD | the calculated PSD |
| [in] | freq | the populated temporal frequency grid defining the frequencies at which the PSD is calculated |
| [in] | zern_j | |
| [in] | layer_i | the 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 323 of file zernikeTemporalPSD.hpp.
References mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_layer_i, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_w, mx::error, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_aosys, mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_apertureRatio4, 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, mx::sigproc::noll_nm(), and WSZ.
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_absTol |
The absolute tolerance to use in the GSL integrator.
Definition at line 136 of file zernikeTemporalPSD.hpp.
| int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_layer_i |
The index of the current layer.
Definition at line 131 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::_relTol |
The relative tolerance to use in the GSL integrator.
Definition at line 137 of file zernikeTemporalPSD.hpp.
| 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 134 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| aosysT* mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_aosys { nullptr } |
Pointer to an AO system structure.
Definition at line 104 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
|
protected |
Aperture ratio, < 1, used for segments.
Definition at line 108 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::apertureRatio(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_cq { 0 } |
The cosine of the wind direction.
Definition at line 127 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_f { 0 } |
the current temporal frequency
Definition at line 122 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_spatialFilter { false } |
Flag indicating if a spatial filter is applied.
Definition at line 129 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_sq { 0 } |
The sine of the wind direction.
Definition at line 128 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| realT mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_j { 0 } |
the current mode number
Definition at line 123 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_m { 0 } |
The current mode m.
Definition at line 124 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().
| int mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::m_zern_n { 0 } |
The current mode n.
Definition at line 125 of file zernikeTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_zernike(), and mx::AO::analysis::zernikeTemporalPSD< _realT, aosysT >::singleLayerPSD().