mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
Interfaces to Lapack and BLAS for Eigen-like arrays.
Functions | |
template<typename eigenT1 , typename eigenT2 > | |
void | mx::math::eigenSYRK (eigenT1 &cv, const eigenT2 &ims) |
Calculates the lower triangular part of the covariance matrix of ims. More... | |
template<typename cvT , typename calcT > | |
MXLAPACK_INT | mx::math::eigenSYEVR (Eigen::Array< calcT, Eigen::Dynamic, Eigen::Dynamic > &eigvec, Eigen::Array< calcT, Eigen::Dynamic, Eigen::Dynamic > &eigval, Eigen::Array< cvT, Eigen::Dynamic, Eigen::Dynamic > &X, int ev0=0, int ev1=-1, char UPLO='L', syevrMem< calcT > *mem=0) |
Calculate select eigenvalues and eigenvectors of an Eigen Array. More... | |
template<typename _evCalcT = double, typename eigenT , typename eigenT1 > | |
MXLAPACK_INT | mx::math::calcKLModes (eigenT &klModes, eigenT &cv, const eigenT1 &Rims, int n_modes=0, syevrMem< _evCalcT > *mem=0, double *t_eigenv=nullptr, double *t_klim=nullptr) |
Calculate the K-L modes, or principle components, given a covariance matrix. More... | |
template<typename dataT > | |
MXLAPACK_INT | mx::math::eigenGESDD (Eigen::Array< dataT,-1,-1 > &U, Eigen::Array< dataT,-1,-1 > &S, Eigen::Array< dataT,-1,-1 > &VT, Eigen::Array< dataT,-1,-1 > &A) |
Compute the SVD of an Eigen::Array using LAPACK's xgesdd. More... | |
template<typename dataT > | |
int | mx::math::eigenPseudoInverse (Eigen::Array< dataT, -1, -1 > &PInv, dataT &condition, int &nRejected, Eigen::Array< dataT, -1, -1 > &U, Eigen::Array< dataT, -1, -1 > &S, Eigen::Array< dataT, -1, -1 > &VT, Eigen::Array< dataT, -1, -1 > &A, dataT &maxCondition, dataT alpha=0, int interact=MX_PINV_NO_INTERACT) |
Calculate the pseudo-inverse of a patrix using the SVD. More... | |
template<typename dataT > | |
int | mx::math::eigenPseudoInverse (Eigen::Array< dataT, -1, -1 > &PInv, dataT &condition, int &nRejected, Eigen::Array< dataT, -1, -1 > &A, dataT &maxCondition, dataT alpha=0, int interact=MX_PINV_NO_INTERACT) |
Calculate the pseudo-inverse of a matrix using the SVD. More... | |
MXLAPACK_INT mx::math::calcKLModes | ( | eigenT & | klModes, |
eigenT & | cv, | ||
const eigenT1 & | Rims, | ||
int | n_modes = 0 , |
||
syevrMem< _evCalcT > * | mem = 0 , |
||
double * | t_eigenv = nullptr , |
||
double * | t_klim = nullptr |
||
) |
Calculate the K-L modes, or principle components, given a covariance matrix.
Eigen-decomposition of the covariance matrix is performed using eigenSYEVR().
evCalcT | is the type in which to perform eigen-decomposition. |
eigenT | is a 2D Eigen-like type |
eigenT1 | is a 2D Eigen-like type. |
[out] | klModes | on exit contains the K-L modes (or P.C.s) |
[in] | cv | a lower-triangle (in the Lapack sense) square covariance matrix. |
[in] | Rims | The reference data. cv.rows() == Rims.cols(). |
[in] | n_modes | [optional] Tbe maximum number of modes to solve for. If 0 all modes are solved for. |
[in] | mem | [optional] A memory structure which can be re-used by SYEVR for efficiency. |
[out] | t_eigenv | [optional] if not null, will be filled in with the time taken to calculate eigenvalues. |
[out] | t_klim | [optional] if not null, will be filled in with the time taken to calculate the KL modes. |
Definition at line 264 of file eigenLapack.hpp.
Referenced by mx::improc::KLIPreduction< _realT, _derotFunctObj, _evCalcT >::worker().
MXLAPACK_INT mx::math::eigenGESDD | ( | Eigen::Array< dataT,-1,-1 > & | U, |
Eigen::Array< dataT,-1,-1 > & | S, | ||
Eigen::Array< dataT,-1,-1 > & | VT, | ||
Eigen::Array< dataT,-1,-1 > & | A | ||
) |
Compute the SVD of an Eigen::Array using LAPACK's xgesdd.
Computes the SVD of A, \( A = U S V^T \).
[out] | U | the A.rows() x A.rows() left matrix |
[out] | S | the A.cols() x 1 matrix of singular values |
[out] | VT | the A.cols() x A.cols() right matrix, note this is the transpose. |
[in] | A | the input matrix to be decomposed |
dataT | is either float or double. |
Definition at line 388 of file eigenLapack.hpp.
Referenced by mx::math::eigenPseudoInverse().
int mx::math::eigenPseudoInverse | ( | Eigen::Array< dataT, -1, -1 > & | PInv, |
dataT & | condition, | ||
int & | nRejected, | ||
Eigen::Array< dataT, -1, -1 > & | A, | ||
dataT & | maxCondition, | ||
dataT | alpha = 0 , |
||
int | interact = MX_PINV_NO_INTERACT |
||
) |
Calculate the pseudo-inverse of a matrix using the SVD.
First computes the SVD of A, \( A = U S V^T \), using eigenGESDD. Then the psuedo-inverse is calculated as \( A^+ = V S^+ U^T\). This interface does not provide access to U, S and VT.
The parameter interact
is intepreted as a bitmask. The values can be
interact
is 0 then no interaction is used and maxCondition controls the inversion.dataT | is either float or double. |
[out] | PInv | The pseudo-inverse of A |
[out] | condition | The final condition number. |
[out] | nRejected | The number of eigenvectors rejected |
[in] | A | The matrix to invert, will be altered! |
[in] | maxCondition | If > 0, the maximum condition number desired. If <0 the number of modes to keep. Used to threshold the singular values. Set to 0 to include all eigenvalues/vectors. Ignored if interactive. |
[in] | alpha | [optional] the Tikhonov regularization value, as a fraction of the highest singular value. If alpha < 0, then it is treated as a (positive) floor (as a fraction of highest singular value) for the singular values, which is not the same as Tikhonov (alpha > 0). |
[in] | interact | [optional] a bitmask controlling interaction. See above. |
Definition at line 591 of file eigenLapack.hpp.
int mx::math::eigenPseudoInverse | ( | Eigen::Array< dataT, -1, -1 > & | PInv, |
dataT & | condition, | ||
int & | nRejected, | ||
Eigen::Array< dataT, -1, -1 > & | U, | ||
Eigen::Array< dataT, -1, -1 > & | S, | ||
Eigen::Array< dataT, -1, -1 > & | VT, | ||
Eigen::Array< dataT, -1, -1 > & | A, | ||
dataT & | maxCondition, | ||
dataT | alpha = 0 , |
||
int | interact = MX_PINV_NO_INTERACT |
||
) |
Calculate the pseudo-inverse of a patrix using the SVD.
First computes the SVD of A, \( A = U S V^T \), using eigenGESDD. Then the psuedo-inverse is calculated as \( A^+ = V S^+ U^T\).
The parameter interact
is intepreted as a bitmask. The values can be
interact
is 0 then no interaction is used and maxCondition
controls the inversion.dataT | is either float or double. |
[out] | PInv | The pseudo-inverse of A |
[out] | condition | The final condition number. |
[out] | nRejected | The number of eigenvectors rejected |
[out] | U | the A.rows() x A.rows() left matrix |
[out] | S | the A.cols() x 1 matrix of singular values |
[out] | VT | the A.cols() x A.cols() right matrix, note this is the transpose. |
[in] | A | The matrix to invert. This will be modified! |
[in] | maxCondition | If > 0, the maximum condition number desired. If <0 the number of modes to keep. Used to threshold the singular values. Set to 0 to include all eigenvalues/vectors. Ignored if interactive. |
[in] | alpha | [optional] the Tikhonov regularization value, as a fraction of the highest singular value. If alpha < 0, then it is treated as a (positive) floor (as a fraction of highest singular value) for the singular values, which is not the same as Tikhonov (alpha > 0). |
[in] | interact | [optional] a bitmask controlling interaction. See above. |
Definition at line 439 of file eigenLapack.hpp.
References mx::math::gnuPlot::command(), mx::math::eigenGESDD(), mx::math::gnuPlot::logy(), mx::math::gnuPlot::plot(), and mx::astro::constants::sigma().
Referenced by mx::sigproc::linearPredictor< _realT >::calcCoefficients().
MXLAPACK_INT mx::math::eigenSYEVR | ( | Eigen::Array< calcT, Eigen::Dynamic, Eigen::Dynamic > & | eigvec, |
Eigen::Array< calcT, Eigen::Dynamic, Eigen::Dynamic > & | eigval, | ||
Eigen::Array< cvT, Eigen::Dynamic, Eigen::Dynamic > & | X, | ||
int | ev0 = 0 , |
||
int | ev1 = -1 , |
||
char | UPLO = 'L' , |
||
syevrMem< calcT > * | mem = 0 |
||
) |
Calculate select eigenvalues and eigenvectors of an Eigen Array.
Uses the templateLapack wrapper for syevr.
cvT | is the scalar type of X (a.k.a. the covariance matrix) |
calcT | is the type in which to calculate the eigenvectors/eigenvalues |
[out] | eigvec | will contain the eigenvectors as columns |
[out] | eigval | will contain the eigenvalues |
[in] | X | is a square matrix which is either upper or lower (default) triangular |
[in] | ev0 | [optional] is the first desired eigenvalue (default = 0) |
[in] | ev1 | [optional] if >= ev0 then this is the last desired eigenvalue. If -1 all eigenvalues are returned. |
[in] | UPLO | [optional] specifies whether X is upper ('U') or lower ('L') triangular. Default is ('L'). |
[in] | mem | [optional] holds the working memory arrays, can be re-passed to avoid unnecessary re-allocations |
Definition at line 141 of file eigenLapack.hpp.
void mx::math::eigenSYRK | ( | eigenT1 & | cv, |
const eigenT2 & | ims | ||
) |
Calculates the lower triangular part of the covariance matrix of ims.
Uses cblas_ssyrk. cv is resized to ims.cols() X ims.cols(). Calculates \( cv = A^T*A \).
eigenT1 | is the eigen matrix/array type of cv. |
eigenT2 | is the eigen matrix/array type of ims |
[out] | cv | is the eigen matrix/array where to store the result |
[in] | ims | is the eigen matrix/array (images as columns) to calculate the covariance of |
Definition at line 65 of file eigenLapack.hpp.
Referenced by mx::improc::KLIPreduction< _realT, _derotFunctObj, _evCalcT >::worker().