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 113 of file fourierTemporalPSD.hpp.
#include <ao/analysis/fourierTemporalPSD.hpp>
Public Types | |
typedef _realT | realT |
The type for arithmetic. More... | |
typedef std::complex< realT > | complexT |
The complex type for arithmetic More... | |
Public Member Functions | |
fourierTemporalPSD () | |
Default c'tor. More... | |
fourierTemporalPSD (bool allocate) | |
Constructor with workspace allocation. More... | |
~fourierTemporalPSD () | |
Destructor. More... | |
realT | fastestPeak (int m, int n) |
Determine the frequency of the highest V-dot-k peak. More... | |
int | singleLayerPSD (std::vector< realT > &PSD, std::vector< realT > &freq, realT m, realT n, int layer_i, int p, 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, realT m, realT n, int p, realT fmax=0) |
Calculate the temporal PSD for a Fourier mode in a multi-layer model. More... | |
void | makePSDGrid (const std::string &dir, int mnMax, realT dFreq, realT maxFreq, realT fmax=0) |
Calculate PSDs over a grid of spatial frequencies. More... | |
int | analyzePSDGrid (const std::string &subDir, const std::string &psdDir, int mnMax, int mnCon, realT gfixed, int lpNc, std::vector< realT > &mags, int lifetimeTrials=0, bool uncontrolledLifetimes=false, bool writePSDs=false, bool writeXfer=false) |
Analyze a PSD grid under closed-loop control. More... | |
int | intensityPSD (const std::string &subDir, const std::string &psdDir, const std::string &CvdPath, int mnMax, int mnCon, const std::string &si_or_lp, std::vector< realT > &mags, int lifetimeTrials, bool writePSDs) |
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... | |
Disk Storage | |
These methods handle writing to and reading from disk. The calculated PSDs are store in the mx::BinVector binary format. A grid of PSDs is specified by its directory name. The directory contains one frequency file (freq.binv), and a set of PSD files, named according to psd_<m>_<n>_.binv. | |
int | getGridFreq (std::vector< realT > &freq, const std::string &dir) |
Get the frequency scale for a PSD grid. More... | |
int | getGridPSD (std::vector< realT > &psd, const std::string &dir, int m, int n) |
Get a single PSD from a PSD grid. More... | |
int | getGridPSD (std::vector< realT > &freq, std::vector< realT > &psd, const std::string &dir, int m, int n) |
Get both the frequency scale and a single PSD from a PSD grid. More... | |
Public Attributes | |
aosysT * | m_aosys {nullptr} |
Pointer to an AO system structure. More... | |
realT | m_f {0} |
the current temporal frequency More... | |
realT | m_m {0} |
the spatial frequency m index More... | |
realT | m_n {0} |
the spatial frequency n index 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... | |
realT | m_f0 {0} |
the Berdja boiling parameter More... | |
int | m_p {1} |
The parity of the mode, +/- 1. If _useBasis==MXAO_FTPSD_BASIS_BASIC then +1 indicates cosine, -1 indicates sine. More... | |
int | _layer_i |
The index of the current layer. More... | |
int | _useBasis |
Set to MXAO_FTPSD_BASIS_BASIC/MODIFIED/PROJECTED_* to use the basic sin/cos modes, the modified Fourier modes, or a projection of them. More... | |
gsl_integration_workspace * | _w |
Workspace 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... | |
int | m_mode_i |
Projected basis mode index. More... | |
Eigen::Array< realT, -1, -1 > | m_modeCoeffs |
Coeeficients of the projection onto the Fourier modes. More... | |
Protected Member Functions | |
void | initialize () |
Initialize parameters to default values. More... | |
typedef std::complex<realT> mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::complexT |
The complex type for arithmetic
Definition at line 119 of file fourierTemporalPSD.hpp.
typedef _realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::realT |
The type for arithmetic.
Definition at line 116 of file fourierTemporalPSD.hpp.
mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::fourierTemporalPSD |
Default c'tor.
Definition at line 402 of file fourierTemporalPSD.hpp.
|
explicit |
Constructor with workspace allocation.
allocate | if true, then the workspace for GSL integrators is allocated. |
Definition at line 409 of file fourierTemporalPSD.hpp.
mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::~fourierTemporalPSD |
Destructor.
Frees GSL workspace if it was allocated.
Definition at line 421 of file fourierTemporalPSD.hpp.
realT mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::absTol |
Get the current absolute tolerance.
Definition at line 446 of file fourierTemporalPSD.hpp.
void mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::absTol | ( | realT | at | ) |
Set absolute tolerance.
at | is the new absolute tolerance. |
Definition at line 440 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::analyzePSDGrid | ( | const std::string & | subDir, |
const std::string & | psdDir, | ||
int | mnMax, | ||
int | mnCon, | ||
realT | gfixed, | ||
int | lpNc, | ||
std::vector< realT > & | mags, | ||
int | lifetimeTrials = 0 , |
||
bool | uncontrolledLifetimes = false , |
||
bool | writePSDs = false , |
||
bool | writeXfer = false |
||
) |
Analyze a PSD grid under closed-loop control.
This always analyzes the simple integrator, and can also analyze the linear predictor controller. Outputs maps of optimum gains, predictor coefficients, variances, and contrasts for a range of guide star magnitudes. Optionally calculates speckle lifetimes. Optionally writes the closed-loop PSDs and transfer functions.
[out] | subDir | the sub-directory of psdDir where to write the results. Is created. |
[in] | psdDir | the directory containing the grid of PSDs. |
[in] | mnMax | the maximum value of m and n in the grid. |
[in] | mnCon | the maximum value of m and n which can be controlled. |
[in] | gfixed | if > 0 then this fixed gain is used in the SI. |
[in] | lpNc | the number of linear predictor coefficients to analyze. If 0 then LP is not analyzed. |
[in] | mags | the guide star magnitudes to analyze for. |
[in] | lifetimeTrials | [optional] number of trials used for calculating speckle lifetimes. If 0, lifetimes are not calculated. |
[in] | uncontrolledLifetimes | [optional] flag controlling whether lifetimes are calculated for uncontrolled modes. |
[in] | writePSDs | [optional] flag controlling if resultant PSDs are saved |
[in] | writeXfer | [optional] flag controlling if resultant Xfer functions are saved |
Definition at line 781 of file fourierTemporalPSD.hpp.
realT mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::fastestPeak | ( | int | m, |
int | n | ||
) |
Determine the frequency of the highest V-dot-k peak.
m | the spatial frequency u index |
n | the spatial frequency v index |
Definition at line 464 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::getGridFreq | ( | std::vector< realT > & | freq, |
const std::string & | dir | ||
) |
Get the frequency scale for a PSD grid.
[out] | freq | the vector to populate with the frequency scale. |
[in] | dir | specifies the directory containing the grid. |
Definition at line 2133 of file fourierTemporalPSD.hpp.
References mx::ioutils::readBinVector().
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::getGridPSD | ( | std::vector< realT > & | freq, |
std::vector< realT > & | psd, | ||
const std::string & | dir, | ||
int | m, | ||
int | n | ||
) |
Get both the frequency scale and a single PSD from a PSD grid.
[out] | freq | the vector to populate with the frequency scale. |
[out] | psd | the vector to populate with the PSD. |
[in] | dir | specifies the directory containing the grid. |
[in] | m | specifies the u component of spatial frequency. |
[in] | n | specifies the v component of spatial frequency. |
Definition at line 2153 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::getGridPSD | ( | std::vector< realT > & | psd, |
const std::string & | dir, | ||
int | m, | ||
int | n | ||
) |
Get a single PSD from a PSD grid.
[out] | psd | the vector to populate with the PSD. |
[in] | dir | specifies the directory containing the grid. |
[in] | m | specifies the u component of spatial frequency. |
[in] | n | specifies the v component of spatial frequency. |
Definition at line 2142 of file fourierTemporalPSD.hpp.
References mx::ioutils::convertToString(), and mx::ioutils::readBinVector().
|
protected |
Initialize parameters to default values.
Definition at line 430 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::intensityPSD | ( | const std::string & | subDir, |
const std::string & | psdDir, | ||
const std::string & | CvdPath, | ||
int | mnMax, | ||
int | mnCon, | ||
const std::string & | si_or_lp, | ||
std::vector< realT > & | mags, | ||
int | lifetimeTrials, | ||
bool | writePSDs | ||
) |
Definition at line 1452 of file fourierTemporalPSD.hpp.
void mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::makePSDGrid | ( | const std::string & | dir, |
int | mnMax, | ||
realT | dFreq, | ||
realT | maxFreq, | ||
realT | fmax = 0 |
||
) |
Calculate PSDs over a grid of spatial frequencies.
The grid of spatial frequencies is square, set by the maximum value of m and n.
The PSDs are written as mx::binVector binary files to a directory. We do not use FITS since this adds overhead and cfitisio handles parallelization poorly due to the limitation on number of created file pointers.
[in] | dir | the directory for output of the PSDs |
[in] | mnMax | the maximum value of m and n in the grid. |
[in] | dFreq | the temporal frequency spacing. |
[in] | maxFreq | the maximum temporal frequency to calculate |
[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 694 of file fourierTemporalPSD.hpp.
References mx::ioutils::convertToString(), mx::ipc::ompLoopWatcher< _outputT, _printPretty, _printLoops, _printPercent, _printNLine, _time >::incrementAndOutputStatus(), mx::math::vectorScale(), and mx::ioutils::writeBinVector().
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::multiLayerPSD | ( | std::vector< realT > & | PSD, |
std::vector< realT > & | freq, | ||
realT | m, | ||
realT | n, | ||
int | p, | ||
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] | m | the first index of the spatial frequency |
[in] | n | the second index of the spatial frequency |
[in] | p | sets which mode is calculated (if basic modes, p = -1 for sine, p = +1 for cosine) |
[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 673 of file fourierTemporalPSD.hpp.
realT mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::relTol |
Get the current relative tolerance.
Definition at line 458 of file fourierTemporalPSD.hpp.
void mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::relTol | ( | realT | rt | ) |
Set relative tolerance.
rt | is the new relative tolerance. |
Definition at line 452 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< realT, aosysT >::singleLayerPSD | ( | std::vector< realT > & | PSD, |
std::vector< realT > & | freq, | ||
realT | m, | ||
realT | n, | ||
int | layer_i, | ||
int | p, | ||
realT | fmax = 0 |
||
) |
Calculate the temporal PSD for a Fourier 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] | m | the first index of the spatial frequency |
[in] | n | the second index of the spatial frequency |
[in] | layer_i | the index of the layer, for accessing the atmosphere parameters |
[in] | p | sets which mode is calculated (if basic modes, p = -1 for sine, p = +1 for cosine) |
[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 487 of file fourierTemporalPSD.hpp.
References mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_layer_i, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_aosys, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_cq, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_m, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_mode_i, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_modeCoeffs, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_n, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_p, mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_spatialFilter, and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_sq.
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_absTol |
The absolute tolerance to use in the GSL integrator.
Definition at line 143 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_layer_i |
The index of the current layer.
Definition at line 136 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), mx::AO::analysis::F_mod(), mx::AO::analysis::Fm_projMod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_relTol |
The relative tolerance to use in the GSL integrator.
Definition at line 144 of file fourierTemporalPSD.hpp.
int mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_useBasis |
Set to MXAO_FTPSD_BASIS_BASIC/MODIFIED/PROJECTED_* to use the basic sin/cos modes, the modified Fourier modes, or a projection of them.
Definition at line 138 of file fourierTemporalPSD.hpp.
gsl_integration_workspace* mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::_w |
Workspace for the gsl integrators, allocated to WSZ if constructed as worker (with allocate == true).
Definition at line 141 of file fourierTemporalPSD.hpp.
aosysT* mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_aosys {nullptr} |
Pointer to an AO system structure.
Definition at line 122 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), mx::AO::analysis::F_mod(), mx::AO::analysis::Fm_projMod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_cq {0} |
The cosine of the wind direction.
Definition at line 127 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_mod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_f {0} |
the current temporal frequency
Definition at line 124 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), mx::AO::analysis::F_mod(), and mx::AO::analysis::Fm_projMod().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_f0 {0} |
the Berdja boiling parameter
Definition at line 133 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_mod().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_m {0} |
the spatial frequency m index
Definition at line 125 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), mx::AO::analysis::F_mod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
int mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_mode_i |
Projected basis mode index.
Definition at line 146 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::Fm_projMod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
Eigen::Array<realT, -1, -1> mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_modeCoeffs |
Coeeficients of the projection onto the Fourier modes.
Definition at line 148 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::Fm_projMod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_n {0} |
the spatial frequency n index
Definition at line 126 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), mx::AO::analysis::F_mod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
int mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_p {1} |
The parity of the mode, +/- 1. If _useBasis==MXAO_FTPSD_BASIS_BASIC then +1 indicates cosine, -1 indicates sine.
Definition at line 135 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_basic(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_spatialFilter {false} |
Flag indicating if a spatial filter is applied.
Definition at line 129 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_mod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().
realT mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::m_sq {0} |
The sine of the wind direction.
Definition at line 128 of file fourierTemporalPSD.hpp.
Referenced by mx::AO::analysis::F_mod(), and mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::singleLayerPSD().