mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
Calculation of the modal covariance in the Fourier basis. More...
Calculation of the modal covariance in the Fourier basis.
Definition in file fourierCovariance.hpp.
#include <gsl/gsl_integration.h>
#include <gsl/gsl_errno.h>
#include "../../math/constants.hpp"
#include "../../math/func/jinc.hpp"
#include "../../sigproc/fourierModes.hpp"
#include "../../ioutils/fits/fitsFile.hpp"
#include "../../improc/eigenImage.hpp"
#include "../../improc/eigenCube.hpp"
#include "../../mxlib.hpp"
#include "../../ipc/ompLoopWatcher.hpp"
#include "../../sys/timeUtils.hpp"
#include "../../math/eigenLapack.hpp"
#include "../../math/func/airyPattern.hpp"
#include "aoAtmosphere.hpp"
#include "aoPSDs.hpp"
#include "aoSystem.hpp"
#include "varmapToImage.hpp"
Go to the source code of this file.
Classes | |
struct | mx::AO::analysis::fourierCovariance< realT, aosysT > |
Structure to manage the Fourier mode covariance calculation, passed to integration functions. More... | |
Namespaces | |
namespace | mx |
The mxlib c++ namespace. | |
Macros | |
#define | WSZ 100000 |
The size of the gsl_integration workspaces. | |
Functions | |
template<typename realT , typename aosysT > | |
realT | mx::AO::analysis::phiInt_basic (realT phi, void *params) |
Worker for the azimuthal integral (in phi) for the basic Fourier mode covariance. | |
template<typename realT , typename aosysT > | |
realT | mx::AO::analysis::phiInt_mod (realT phi, void *params) |
Worker for the azimuthal integral (in phi) for the modified Fourier mode covariance. | |
template<typename realT , typename aosysT > | |
realT | mx::AO::analysis::kInt (realT k, void *params) |
Worker function for the radial integral in the covariance calculation. | |
template<typename realT , typename aosysT > | |
int | mx::AO::analysis::fourierVarVec (const std::string &fname, int N, aosysT &aosys, realT absTol, realT relTol, bool modifed=true) |
Calculate a vector of Fourier mode variances. | |
template<typename realT , typename aosysT > | |
int | mx::AO::analysis::fourierPSDMap (improc::eigenImage< realT > &var, improc::eigenImage< realT > &psd, int N, int overSample, aosysT &aosys) |
Calculate a map of Fourier variances by convolution with the PSD. | |
template<typename realT > | |
int | mx::AO::analysis::fourierCovarMap (const std::string &fname, int N, realT D, realT L_0, bool subPist, bool subTilt, realT absTol, realT relTol, bool modified=true) |
#define WSZ 100000 |
The size of the gsl_integration workspaces.
Definition at line 44 of file fourierCovariance.hpp.
int mx::AO::analysis::fourierCovarMap | ( | const std::string & | fname, |
int | N, | ||
realT | D, | ||
realT | L_0, | ||
bool | subPist, | ||
bool | subTilt, | ||
realT | absTol, | ||
realT | relTol, | ||
bool | modified = true |
||
) |
[out] | fname | the path where the output FITS file will be written |
[in] | N | the linear number of Fourier modes across the aperture. The Nyquist frequency is set by N/2. |
Definition at line 463 of file fourierCovariance.hpp.
References mx::AO::analysis::fourierCovariance< realT, aosysT >::absTol, mx::AO::analysis::fourierCovariance< realT, aosysT >::aosys, mx::fits::fitsHeader::append(), mx::AO::analysis::aoSystem< _realT, _inputSpectT, iosT >::D(), mx::error, mx::AO::analysis::fourierCovarMap(), mx::AO::analysis::fourierCovariance< realT, aosysT >::getVariance(), mx::ipc::ompLoopWatcher< _outputT, _printPretty, _printLoops, _printPercent, _printNLine, _time >::incrementAndOutputStatus(), mx::AO::analysis::aoAtmosphere< _realT >::L_0(), mx::AO::analysis::aoSystem< _realT, _inputSpectT, iosT >::loadMagAOX(), mx::AO::analysis::fourierCovariance< realT, aosysT >::m, mx::sigproc::makeFourierModeFreqs_Rect(), mx::AO::analysis::modified, mx::AO::analysis::fourierCovariance< realT, aosysT >::mp, mx::AO::analysis::fourierCovariance< realT, aosysT >::n, mx::AO::analysis::fourierCovariance< realT, aosysT >::np, mx::AO::analysis::fourierCovariance< realT, aosysT >::p, mx::AO::analysis::fourierCovariance< realT, aosysT >::pp, mx::AO::analysis::aoAtmosphere< _realT >::r_0(), mx::AO::analysis::fourierCovariance< realT, aosysT >::relTol, mx::AO::analysis::fourierCovariance< realT, aosysT >::useBasic, and mx::fits::fitsFile< dataT >::write().
Referenced by mx::AO::analysis::fourierCovarMap().
int mx::AO::analysis::fourierPSDMap | ( | improc::eigenImage< realT > & | var, |
improc::eigenImage< realT > & | psd, | ||
int | N, | ||
int | overSample, | ||
aosysT & | aosys | ||
) |
Calculate a map of Fourier variances by convolution with the PSD.
Uses the Airy pattern for the circularly unobstructed aperture.
[out] | var | The variance estimated by convolution with the PSD |
[out] | psd | the PSD map |
[in] | N | the number of components to analyze |
aosys | [in[ the AO system defining the PSD characteristics. |
Definition at line 405 of file fourierCovariance.hpp.
References mx::math::func::airyPattern(), mx::AO::analysis::fourierPSDMap(), mx::math::six_fifths(), and mx::AO::analysis::varmapToImage().
Referenced by mx::AO::analysis::fourierPSDMap().
int mx::AO::analysis::fourierVarVec | ( | const std::string & | fname, |
int | N, | ||
aosysT & | aosys, | ||
realT | absTol, | ||
realT | relTol, | ||
bool | modifed = true |
||
) |
Calculate a vector of Fourier mode variances.
Definition at line 321 of file fourierCovariance.hpp.
References mx::AO::analysis::fourierCovariance< realT, aosysT >::absTol, mx::AO::analysis::fourierCovariance< realT, aosysT >::aosys, mx::error, mx::AO::analysis::fourierVarVec(), mx::AO::analysis::fourierCovariance< realT, aosysT >::getVariance(), mx::ipc::ompLoopWatcher< _outputT, _printPretty, _printLoops, _printPercent, _printNLine, _time >::incrementAndOutputStatus(), mx::AO::analysis::fourierCovariance< realT, aosysT >::m, mx::AO::analysis::fourierCovariance< realT, aosysT >::mnCon, mx::AO::analysis::fourierCovariance< realT, aosysT >::mp, mx::AO::analysis::fourierCovariance< realT, aosysT >::n, mx::AO::analysis::fourierCovariance< realT, aosysT >::np, mx::AO::analysis::fourierCovariance< realT, aosysT >::p, mx::AO::analysis::fourierCovariance< realT, aosysT >::pp, mx::AO::analysis::fourierCovariance< realT, aosysT >::relTol, and mx::math::six_fifths().
Referenced by mx::AO::analysis::fourierVarVec().
realT mx::AO::analysis::kInt | ( | realT | k, |
void * | params | ||
) |
Worker function for the radial integral in the covariance calculation.
k | the spatial frequency at which to evaluate the integrand |
params | a pointer to a object of type fourierCovariance<realT, aosyT> |
realT | a floating point type used for all calculations. As of Nov 2016 must be double due to gsl_integration. |
aosysT | the type of the AO system structure |
Definition at line 205 of file fourierCovariance.hpp.
References mx::error, mx::AO::analysis::fourierCovariance< realT, aosysT >::k, mx::AO::analysis::kInt(), mx::AO::analysis::fourierCovariance< realT, aosysT >::phi_w, mx::AO::analysis::fourierCovariance< realT, aosysT >::useBasic, and WSZ.
Referenced by mx::AO::analysis::kInt().
realT mx::AO::analysis::phiInt_mod | ( | realT | phi, |
void * | params | ||
) |
Worker for the azimuthal integral (in phi) for the modified Fourier mode covariance.
phi | the angle at which to evaluate the integrand |
params | a pointer to a object of type fourierCovariance<realT, aosyT> |
realT | a floating point type used for all calculations |
aosysT | the type of the AO system structure |
Definition at line 127 of file fourierCovariance.hpp.
References mx::AO::analysis::fourierCovariance< realT, aosysT >::aosys, mx::math::func::jinc(), mx::AO::analysis::fourierCovariance< realT, aosysT >::k, mx::AO::analysis::fourierCovariance< realT, aosysT >::m, mx::AO::analysis::fourierCovariance< realT, aosysT >::mnCon, mx::AO::analysis::fourierCovariance< realT, aosysT >::mp, mx::AO::analysis::fourierCovariance< realT, aosysT >::n, mx::AO::analysis::fourierCovariance< realT, aosysT >::np, mx::AO::analysis::fourierCovariance< realT, aosysT >::p, mx::AO::analysis::phiInt_mod(), mx::AO::analysis::fourierCovariance< realT, aosysT >::pp, and mx::math::six_fifths().
Referenced by mx::AO::analysis::phiInt_mod().