mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
Structure to manage the zernike mode covariance calculation, passed to integration functions.
realT | a floating point type used for all calculations. As of Dec 2023 must be double due to gsl_integration. |
aosysT | the type of the AO system structure |
Definition at line 45 of file zernikeCovariance.hpp.
#include <ao/analysis/zernikeCovariance.hpp>
Public Member Functions | |
zernikeCovariance () | |
Constructor. | |
~zernikeCovariance () | |
Destructor. | |
realT | getCovariance (realT &error) |
Calculate the covariance between the two modes. | |
Static Public Member Functions | |
static realT | kInt (realT k, void *params) |
Worker function for the radial integral in the covariance calculation. | |
static realT | phiInt (realT phi, void *params) |
Worker for the azimuthal integral (in phi) for the zernike mode covariance. | |
Public Attributes | |
aosysT * | m_aosys { nullptr } |
Pointer to an AO system, which contains the relevant spatial PSD of turbulence. | |
realT | m_n |
The n-index of the unprimed mode. | |
realT | m_m |
The m-index of the unprimed mode. | |
realT | m_np |
The n-indexof the primed mode, corresponding to the \( k_v = n/D \) component of spatial frequency. | |
realT | m_mp |
The m-index of the primed mode, corresponding to the \( k_u = m/D \) component of spatial frequency. | |
realT | m_k |
Spatial frequency being calculated, passed for use in the integrand worker functions. | |
realT | m_kIntEpsAbs { 1e-10 } |
Absolute tolerance for the radial integral. Default is 1e-10. | |
realT | m_kIntEpsRel { 0 } |
Relative tolerance for the radial integral. Default is 0, meaning absolute is used. | |
realT | m_phiIntEpsAbs { 1e-10 } |
Absolute tolerance for the azimuthal integral. Default is 1e-10. | |
realT | m_phiIntEpsRel { 0 } |
Relative tolerance for the azimuthal integral. Default is 0, meaning absolute is used. | |
Protected Attributes | |
gsl_integration_workspace * | m_workspacePhi { nullptr } |
Working memory for the azimuthal integral. | |
gsl_integration_workspace * | m_workspaceK { nullptr } |
Working memory for the radial integral. | |
|
inline |
Constructor.
Definition at line 88 of file zernikeCovariance.hpp.
References mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_workspaceK, and mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_workspacePhi.
|
inline |
Destructor.
Definition at line 95 of file zernikeCovariance.hpp.
References mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_workspaceK, and mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_workspacePhi.
|
inline |
Calculate the covariance between the two modes.
document me
handle gsl errors
Definition at line 124 of file zernikeCovariance.hpp.
References mx::AO::analysis::zernikeCovariance< realT, aosysT >::kInt(), mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_aosys, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_kIntEpsAbs, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_kIntEpsRel, and mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_workspaceK.
|
static |
Worker function for the radial integral in the covariance calculation.
k
. [in] | k | the spatial frequency at which to evaluate the integrand, meters |
[in] | params | a pointer to a object of type zernikeCovariance<realT, aosyT> |
Definition at line 164 of file zernikeCovariance.hpp.
References mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_k, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_phiIntEpsAbs, and mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_phiIntEpsRel.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::getCovariance().
|
static |
Worker for the azimuthal integral (in phi) for the zernike mode covariance.
phi
[in] | phi | the angle at which to evaluate the integrand, radians |
[in] | params | a pointer to a object of type zernikeCovariance<realT, aosyT> |
Definition at line 193 of file zernikeCovariance.hpp.
References mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_aosys, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_k, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_m, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_mp, mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_n, and mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_np.
aosysT* mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_aosys { nullptr } |
Pointer to an AO system, which contains the relevant spatial PSD of turbulence.
Definition at line 48 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::getCovariance(), and mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_k |
Spatial frequency being calculated, passed for use in the integrand worker functions.
Definition at line 63 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::kInt(), and mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_kIntEpsAbs { 1e-10 } |
Absolute tolerance for the radial integral. Default is 1e-10.
Definition at line 66 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::getCovariance().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_kIntEpsRel { 0 } |
Relative tolerance for the radial integral. Default is 0, meaning absolute is used.
Definition at line 69 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::getCovariance().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_m |
The m-index of the unprimed mode.
Definition at line 54 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_mp |
The m-index of the primed mode, corresponding to the \( k_u = m/D \) component of spatial frequency.
Definition at line 60 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_n |
The n-index of the unprimed mode.
Definition at line 51 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_np |
The n-indexof the primed mode, corresponding to the \( k_v = n/D \) component of spatial frequency.
Definition at line 57 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::phiInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_phiIntEpsAbs { 1e-10 } |
Absolute tolerance for the azimuthal integral. Default is 1e-10.
Definition at line 72 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::kInt().
realT mx::AO::analysis::zernikeCovariance< realT, aosysT >::m_phiIntEpsRel { 0 } |
Relative tolerance for the azimuthal integral. Default is 0, meaning absolute is used.
Definition at line 75 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::kInt().
|
protected |
Working memory for the radial integral.
Definition at line 84 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::zernikeCovariance(), mx::AO::analysis::zernikeCovariance< realT, aosysT >::~zernikeCovariance(), and mx::AO::analysis::zernikeCovariance< realT, aosysT >::getCovariance().
|
protected |
Working memory for the azimuthal integral.
Definition at line 81 of file zernikeCovariance.hpp.
Referenced by mx::AO::analysis::zernikeCovariance< realT, aosysT >::zernikeCovariance(), and mx::AO::analysis::zernikeCovariance< realT, aosysT >::~zernikeCovariance().