mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
Namespaces | |
namespace | analysis |
namespace | constants |
namespace | path |
Namespace for paths. | |
namespace | sim |
Classes | |
struct | influenceFunctionGaussianSpec |
Functions | |
template<typename realT > | |
void | applyPupil2Basis (eigenCube< realT > &modes, const std::string &basisName, const std::string &pupilName, realT fwhm=0) |
Multiply a raw modal basis by a pupil mask. | |
template<typename realT , typename spectRealT > | |
int | orthogonalizeBasis (eigenCube< realT > &ortho, Eigen::Array< spectRealT, -1, -1 > &spect, eigenCube< realT > &modes, int method, realT psum=1.0) |
template<typename realT > | |
void | orthogonalizeBasis (const std::string &orthoName, const std::string &basisName, const std::string &pupilName, int method) |
Calculate an orthogonal projection of a basis on a pupil. | |
template<typename spectRealT , typename realT > | |
void | normalizeSpectrum (mx::improc::eigenImage< spectRealT > &spectrum, mx::improc::eigenCube< realT > &modes, mx::improc::eigenCube< realT > &ortho, mx::improc::eigenImage< realT > &pupil) |
template<typename realT > | |
int | slaveBasis (const std::string &outputBasisN, const std::string &inputBasisN, const std::string &pupilN, const std::string &dmN, realT fwhm, realT fsmooth, int firstMode=0) |
template<typename realT > | |
int | apodizeBasis (const std::string &outputBasisN, const std::string &inputBasisN, realT tukeyAlpha, realT centralObs, realT overScan, int firstMode) |
template<typename realT > | |
int | subtractBasis (const std::string &basisName, const std::string &basisName1, const std::string &basisName2, const std::string &pupilName, int method) |
template<typename realT > | |
void | makeModfBasis (const std::string &basisName, int dim, int N, realT ang, int nZern=0) |
Make the modified Fourier basis. | |
template<typename realT > | |
void | influenceFunctionsGaussian (const std::string &dmName, realT dmSz, int linNAct, realT diameter, realT coupling, realT couplingRange, realT pupilSz=0, bool offsetOdd=false) |
Create a set of Gaussian influence functions (IFs) on a square grid. | |
template<typename realT > | |
void | ifPInv (const std::string &dmName, realT maxCondition=-1) |
Calculate the pseudo-inverse and mirror modes for a set of influence functions. | |
template<typename realT > | |
void | m2cMatrix (Eigen::Array< realT, -1, -1 > &M2c, Eigen::Array< realT, -1, -1 > &Ap, mx::improc::eigenCube< realT > &M) |
Calculate the modes-to-commands matrix for a set of modes. | |
template<typename realT > | |
void | m2cMatrix (const std::string &dmName, const std::string &basisName) |
Calculate the modes-to-commands matrix for a set of modes. | |
template<typename realT > | |
void | modalDMM2cMatrix (const std::string &basisName) |
Generate the modes-to-commands matrix for a set of modes on the modal DM. | |
template<typename realT > | |
void | modesOnDM (const std::string &dmName, const std::string &basisName) |
Calculate the basis set as it will be reproduced by the DM. | |
template<typename realT > | |
void | circularPupil (const std::string &pupilName, realT pupilDiamPixels, realT pupilDiamMeters, realT centralObs=0, realT overscan=0) |
Generate a circular pupil and saves it to disk. | |
template<typename realT > | |
void | circularApodizedPupil (const std::string &pupilName, int pupilDiamPixels, realT pupilDiamMeters, realT tukeyAlpha, realT centralObs=0, realT overScan=0) |
Generates a circular apodized pupil and saves it to disk. | |
template<typename realT > | |
void | makeZernikeBasis (const std::string &basisName, const std::string &pupilName, int dim, int N) |
Make the Zernike basis. | |
int mx::AO::apodizeBasis | ( | const std::string & | outputBasisN, |
const std::string & | inputBasisN, | ||
realT | tukeyAlpha, | ||
realT | centralObs, | ||
realT | overScan, | ||
int | firstMode | ||
) |
Definition at line 371 of file basis.hpp.
References mx::AO::path::basis::modes(), and mx::math::six_fifths().
void mx::AO::applyPupil2Basis | ( | eigenCube< realT > & | modes, |
const std::string & | basisName, | ||
const std::string & | pupilName, | ||
realT | fwhm = 0 |
||
) |
Multiply a raw modal basis by a pupil mask.
Definition at line 38 of file basis.hpp.
References mx::improc::filterImage(), mx::AO::path::basis::modes(), mx::improc::padImage(), mx::AO::path::pupil::pupilFile(), and mx::math::six_fifths().
Referenced by orthogonalizeBasis().
void mx::AO::circularApodizedPupil | ( | const std::string & | pupilName, |
int | pupilDiamPixels, | ||
realT | pupilDiamMeters, | ||
realT | tukeyAlpha, | ||
realT | centralObs = 0 , |
||
realT | overScan = 0 |
||
) |
Generates a circular apodized pupil and saves it to disk.
Apodization is with a Tukey window.
Definition at line 62 of file pupil.hpp.
References mx::fits::fitsHeader::append(), mx::AO::path::pupil::pupilFile(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::circularPupil | ( | const std::string & | pupilName, |
realT | pupilDiamPixels, | ||
realT | pupilDiamMeters, | ||
realT | centralObs = 0 , |
||
realT | overscan = 0 |
||
) |
Generate a circular pupil and saves it to disk.
Definition at line 28 of file pupil.hpp.
References mx::fits::fitsHeader::append(), mx::wfp::circularPupil(), mx::AO::path::pupil::pupilFile(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::ifPInv | ( | const std::string & | dmName, |
realT | maxCondition = -1 |
||
) |
Calculate the pseudo-inverse and mirror modes for a set of influence functions.
The influence functions are a cube of 2D maps of mirror deformations, and the matrix \( A \) is formed by vectorizing the influence functions. The pseudo inverse is then calculated from the SVD as follows: \( A = U S V^T \) which gives the Moore-Penrose pseudo-inverse: \( A^+ = V S^+ U^T\). The columns of U contain the orthogonal mirror modes, and S contains the singular values. These are both written to disk along with the pseudo-inverse.
[in] | dmName | is the name of the deformable mirror. |
[in] | maxCondition | [optional] the maximum condition number to accept for the pseudo-inverse. If < 0 [default] then the user is interactively asked what to use. |
realT | is the real floating point type in which to do calculations. |
Definition at line 366 of file influenceFunctions.hpp.
References mx::fits::fitsHeader::append(), mx::improc::eigenCube< dataT >::asVectors(), mx::improc::eigenCube< dataT >::cols(), mx::math::eigenPseudoInverse(), mx::AO::path::dm::influenceFunctions(), mx::AO::path::dm::mirrorModes(), MX_PINV_ASK, MX_PINV_NO_INTERACT, MX_PINV_PLOT, mx::improc::eigenCube< dataT >::planes(), mx::AO::path::dm::pseudoInverse(), mx::fits::fitsFile< dataT >::read(), mx::improc::eigenCube< dataT >::resize(), mx::improc::eigenCube< dataT >::rows(), mx::AO::path::dm::singularValues(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::influenceFunctionsGaussian | ( | const std::string & | dmName, |
realT | dmSz, | ||
int | linNAct, | ||
realT | diameter, | ||
realT | coupling, | ||
realT | couplingRange, | ||
realT | pupilSz = 0 , |
||
bool | offsetOdd = false |
||
) |
Create a set of Gaussian influence functions (IFs) on a square grid.
The width of the IFs is specified by the coupling, which is the height of the IF at 1 actuator separation. For a FWHM = 1 actuator, coupling = 0.0625. The IFs are normalized to a height of 1.
[in] | dmName | the name of the deformable mirror (DM) (the mx::AO path will be constructed from this). |
[in] | dmSz | the linear size of the DM in pixels. |
[in] | linNAct | the linear number of actuators across the DM. |
[in] | diameter | is the diameter of the DM, in actuators. If 0 then the DM is a square. |
[in] | coupling | is the height of the IF at 1 actuator separation. E.G., for a FWHM = 1 actuator, set this to 0.0625. |
[in] | pupilSz | is the size of the pupil, and therefore the size of the maps (pupilSz x pupilSz), Is set to dmSz if 0, can be pupilSz < dmSx for an oversized DM. The pupil is assumed to be centered. |
realT | is the real floating point type in which to do calculations. |
< Uniform deviate, used in shiftRandom.
couplingRange | Uniformly distributed uncertainy in coupling, in fraction of the coupling. |
Definition at line 230 of file influenceFunctions.hpp.
References mx::AO::path::dm::actuatorPositions(), mx::fits::fitsHeader::append(), mx::math::func::fwhm2sigma(), mx::improc::eigenCube< dataT >::image(), mx::AO::path::dm::influenceFunctions(), mx::math::randomT< typeT, _ranengT, _randistT >::seed(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::m2cMatrix | ( | const std::string & | dmName, |
const std::string & | basisName | ||
) |
Calculate the modes-to-commands matrix for a set of modes.
A wrapper for m2cMatrix, where the various matrices are here specified with names, which in turn generate mx::AO paths
[in] | dmName | is the name of the deformable mirror |
[in] | basisName | is the name of the modal basis |
[in] | pupilName | is the name of the pupil |
realT | is the real floating point type in which to do calculations. |
Definition at line 470 of file influenceFunctions.hpp.
References mx::AO::path::dm::M2c(), m2cMatrix(), mx::AO::path::basis::modes(), mx::AO::path::dm::pseudoInverse(), mx::fits::fitsFile< dataT >::read(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::m2cMatrix | ( | Eigen::Array< realT, -1, -1 > & | M2c, |
Eigen::Array< realT, -1, -1 > & | Ap, | ||
mx::improc::eigenCube< realT > & | M | ||
) |
Calculate the modes-to-commands matrix for a set of modes.
Given the pseudo-inverse \( A^+ \) of the influence functions, the commands to the actuators to take a mirror shape \( \vec{s} \) are calculated as
\[ \vec{c} = \mathbf{A^+} \vec{s}. \]
Now given a basis set \( \mathbf{M} \), the mirror shape is specified as
\[ \vec{s}= \sum_i m_i M_i \]
where \( \vec{m} = m_0, m_1, \cdot m_i \cdot \) are the modal coefficients. If the mirror-to-commands matrix, \( \mathbf{M2c} \) gives commands as
\[ \vec{c} = \mathbf{M2c} \vec{m} \]
then we can calculate \( \mathbf{M2c} \) as:
\[ \mathbf{M2c} = \mathbf{A}^{+T} \mathbf{M} \]
[out] | M2c | is the M2c matrix, allocated during the calculationg |
[in] | Ap | is the pseudo-inverse of the influence function basis set |
[in] | M | is the modal basis set, a cube of shapes |
realT | is the real floating point type used for all calculations. |
Definition at line 453 of file influenceFunctions.hpp.
References mx::improc::eigenCube< dataT >::asVectors().
Referenced by m2cMatrix().
void mx::AO::makeModfBasis | ( | const std::string & | basisName, |
int | dim, | ||
int | N, | ||
realT | ang, | ||
int | nZern = 0 |
||
) |
Make the modified Fourier basis.
[in] | basisName | the name of the basis (not including the mx::AO path) |
[in] | dim | the linear size of the maps, that is they will be dimxdim in size. |
[in] | N | is the number of degrees of freedom. Number of modes will be (N+1)(N+1) - 1. |
realT | the real numeric type for calculations |
Definition at line 32 of file fourierBasis.hpp.
References mx::improc::eigenCube< dataT >::image(), mx::sigproc::makeFourierBasis_Rect(), mx::AO::path::basis::modes(), MX_FOURIER_MODIFIED, mx::improc::eigenCube< dataT >::resize(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::makeZernikeBasis | ( | const std::string & | basisName, |
const std::string & | pupilName, | ||
int | dim, | ||
int | N | ||
) |
Make the Zernike basis.
[in] | basisName | the name of the basis (not including the mx::AO path) |
[in] | dim | the linear size of the maps, that is they will be dimxdim in size. |
[in] | N | is the number of degrees of freedom. Number of modes will be (N+1)(N+1) - 1. |
realT | the real numeric type for calculations |
Definition at line 29 of file zernikeBasis.hpp.
References mx::improc::eigenCube< dataT >::image(), mx::AO::path::basis::modes(), mx::improc::eigenCube< dataT >::planes(), mx::AO::path::pupil::pupilFile(), and mx::improc::eigenCube< dataT >::resize().
void mx::AO::modalDMM2cMatrix | ( | const std::string & | basisName | ) |
Generate the modes-to-commands matrix for a set of modes on the modal DM.
Generates an Identity matrix of the appropriate size.
[in] | basisName | is the name of the modal basis |
realT | is the real floating point type in which to do calculations. |
Definition at line 502 of file influenceFunctions.hpp.
References mx::AO::path::dm::M2c(), mx::AO::path::basis::modes(), mx::improc::eigenCube< dataT >::planes(), mx::fits::fitsFile< dataT >::read(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::modesOnDM | ( | const std::string & | dmName, |
const std::string & | basisName | ||
) |
Calculate the basis set as it will be reproduced by the DM.
The projected modes are written to the path specified by mx::AO::path::dm::projectedModes().
[in] | dmName | is the name of the DM. |
[in] | basisName | is the name of the basis. |
realT | is the real floating point type in which to do calculations. |
Definition at line 540 of file influenceFunctions.hpp.
References mx::improc::eigenCube< dataT >::cols(), mx::improc::eigenCube< dataT >::image(), mx::AO::path::dm::influenceFunctions(), mx::AO::path::dm::M2c(), mx::AO::path::dm::projectedModes(), mx::fits::fitsFile< dataT >::read(), mx::improc::eigenCube< dataT >::resize(), mx::improc::eigenCube< dataT >::rows(), mx::improc::eigenCube< dataT >::setZero(), and mx::fits::fitsFile< dataT >::write().
void mx::AO::normalizeSpectrum | ( | mx::improc::eigenImage< spectRealT > & | spectrum, |
mx::improc::eigenCube< realT > & | modes, | ||
mx::improc::eigenCube< realT > & | ortho, | ||
mx::improc::eigenImage< realT > & | pupil | ||
) |
Definition at line 202 of file basis.hpp.
References mx::improc::eigenCube< dataT >::image(), mx::math::six_fifths(), and mx::math::vectorMedianInPlace().
void mx::AO::orthogonalizeBasis | ( | const std::string & | orthoName, |
const std::string & | basisName, | ||
const std::string & | pupilName, | ||
int | method | ||
) |
Calculate an orthogonal projection of a basis on a pupil.
Can use either the stabilized Gramm Schmidt (SGS) procedure, or Singular Value Decomposition (SVD). This calls applyPupil2Basis as a first step..
[in] | basisName | the name of the basis to orthogonalize |
[in] | pupilName | the name of the pupil on which to orthogonalize. |
[in] | method | either MXAO_ORTHO_METHOD_SGS (for SGS) or MXAO_ORTHO_METHOD_SVD (for SVD) |
realT |
Definition at line 143 of file basis.hpp.
References applyPupil2Basis(), mx::improc::eigenCube< dataT >::image(), mx::AO::path::basis::modes(), orthogonalizeBasis(), mx::improc::eigenCube< dataT >::planes(), mx::AO::path::pupil::pupilFile(), mx::math::six_fifths(), and mx::AO::path::basis::spectrum().
int mx::AO::orthogonalizeBasis | ( | eigenCube< realT > & | ortho, |
Eigen::Array< spectRealT, -1, -1 > & | spect, | ||
eigenCube< realT > & | modes, | ||
int | method, | ||
realT | psum = 1.0 |
||
) |
Definition at line 89 of file basis.hpp.
References mx::improc::eigenCube< dataT >::cols(), mx::improc::eigenCube< dataT >::data(), mx::math::eigenGESDD(), mx::improc::eigenCube< dataT >::planes(), mx::improc::eigenCube< dataT >::resize(), mx::improc::eigenCube< dataT >::rows(), and mx::math::six_fifths().
Referenced by orthogonalizeBasis(), and subtractBasis().
int mx::AO::slaveBasis | ( | const std::string & | outputBasisN, |
const std::string & | inputBasisN, | ||
const std::string & | pupilN, | ||
const std::string & | dmN, | ||
realT | fwhm, | ||
realT | fsmooth, | ||
int | firstMode = 0 |
||
) |
Definition at line 266 of file basis.hpp.
References mx::improc::cutPaddedImage(), mx::improc::filterImage(), mx::AO::path::dm::influenceFunctions(), mx::AO::path::basis::modes(), mx::improc::padImage(), mx::AO::path::pupil::pupilFile(), mx::improc::eigenCube< dataT >::resize(), mx::improc::eigenCube< dataT >::rows(), and mx::math::six_fifths().
int mx::AO::subtractBasis | ( | const std::string & | basisName, |
const std::string & | basisName1, | ||
const std::string & | basisName2, | ||
const std::string & | pupilName, | ||
int | method | ||
) |
Definition at line 412 of file basis.hpp.
References mx::improc::eigenCube< dataT >::cols(), mx::improc::cutPaddedImage(), mx::improc::filterImage(), mx::improc::eigenCube< dataT >::image(), mx::AO::path::basis::modes(), orthogonalizeBasis(), mx::improc::padImage(), mx::AO::path::pupil::pupilFile(), mx::improc::eigenCube< dataT >::rows(), and mx::math::six_fifths().