ANIMA
4.0
|
Cholesky decomposition: decomposes a symmetric matrix A in the form L D L^T, where L is lower triangular and D diagonal. May be used to solve efficiently any linear system using the solve methods. Refer to Gill, Golub et al. Methods for modifying matrix factorizations. 1974. More...
#include <animaCholeskyDecomposition.h>
Public Types | |
typedef vnl_diag_matrix< double > | DiagonalType |
typedef vnl_matrix< double > | MatrixType |
typedef vnl_vector< double > | VectorType |
Public Member Functions | |
CholeskyDecomposition () | |
CholeskyDecomposition (const unsigned int matrixDimension, const double epsilon) | |
CholeskyDecomposition (const MatrixType &val) | |
double | GetConditionNumber () |
Get condition number from decomposition. More... | |
DiagonalType & | GetDMatrix () |
MatrixType & | GetInputMatrix () |
Get input matrix, useful if decomposition updated through rank one method. More... | |
MatrixType & | GetLMatrix () |
void | PerformDecomposition () |
Performs LDL decomposition from input matrix. More... | |
void | Recompose () |
Recompose input matrix from LDL decomposition. More... | |
void | SetDMatrix (const DiagonalType &val) |
void | SetInputMatrix (const MatrixType &matrix) |
Set input matrix to decompose. More... | |
void | SetLMatrix (const MatrixType &val) |
VectorType & | SolveLinearSystem (const VectorType &b) |
Solves linear system Ax=b and outputs result in new variable. More... | |
void | SolveLinearSystemInPlace (VectorType &b) |
Solves linear system Ax=b and outputs result in input b variable. More... | |
void | Update (const VectorType &x) |
Update decomposition with x so that LDL matches A + x x^T. More... | |
Cholesky decomposition: decomposes a symmetric matrix A in the form L D L^T, where L is lower triangular and D diagonal. May be used to solve efficiently any linear system using the solve methods. Refer to Gill, Golub et al. Methods for modifying matrix factorizations. 1974.
Definition at line 17 of file animaCholeskyDecomposition.h.
typedef vnl_diag_matrix<double> anima::CholeskyDecomposition::DiagonalType |
Definition at line 21 of file animaCholeskyDecomposition.h.
typedef vnl_matrix<double> anima::CholeskyDecomposition::MatrixType |
Definition at line 20 of file animaCholeskyDecomposition.h.
typedef vnl_vector<double> anima::CholeskyDecomposition::VectorType |
Definition at line 22 of file animaCholeskyDecomposition.h.
|
inline |
Definition at line 24 of file animaCholeskyDecomposition.h.
|
inline |
Definition at line 26 of file animaCholeskyDecomposition.h.
|
inline |
Definition at line 38 of file animaCholeskyDecomposition.h.
double anima::CholeskyDecomposition::GetConditionNumber | ( | ) |
Get condition number from decomposition.
Definition at line 41 of file animaCholeskyDecomposition.cxx.
|
inline |
Definition at line 70 of file animaCholeskyDecomposition.h.
|
inline |
Get input matrix, useful if decomposition updated through rank one method.
Definition at line 65 of file animaCholeskyDecomposition.h.
|
inline |
Definition at line 68 of file animaCholeskyDecomposition.h.
void anima::CholeskyDecomposition::PerformDecomposition | ( | ) |
Performs LDL decomposition from input matrix.
Definition at line 14 of file animaCholeskyDecomposition.cxx.
Referenced by anima::GaussianMCMVariableProjectionCost::PrepareDataForLLS(), anima::B1GMMRelaxometryCostFunction::PrepareDataForLLS(), anima::B1GammaMixtureT2RelaxometryCostFunction::PrepareDataForLLS(), and anima::NNLSOptimizer::StartOptimization().
void anima::CholeskyDecomposition::Recompose | ( | ) |
Recompose input matrix from LDL decomposition.
Definition at line 108 of file animaCholeskyDecomposition.cxx.
|
inline |
Definition at line 69 of file animaCholeskyDecomposition.h.
void anima::CholeskyDecomposition::SetInputMatrix | ( | const MatrixType & | matrix | ) |
Set input matrix to decompose.
Definition at line 8 of file animaCholeskyDecomposition.cxx.
Referenced by anima::GaussianMCMVariableProjectionCost::PrepareDataForLLS(), anima::B1GMMRelaxometryCostFunction::PrepareDataForLLS(), anima::B1GammaMixtureT2RelaxometryCostFunction::PrepareDataForLLS(), and anima::NNLSOptimizer::StartOptimization().
|
inline |
Definition at line 67 of file animaCholeskyDecomposition.h.
CholeskyDecomposition::VectorType & anima::CholeskyDecomposition::SolveLinearSystem | ( | const VectorType & | b | ) |
Solves linear system Ax=b and outputs result in new variable.
Definition at line 62 of file animaCholeskyDecomposition.cxx.
References SolveLinearSystemInPlace().
Referenced by anima::GaussianMCMVariableProjectionCost::PrepareDataForLLS(), anima::B1GMMRelaxometryCostFunction::PrepareDataForLLS(), and anima::B1GammaMixtureT2RelaxometryCostFunction::PrepareDataForLLS().
void anima::CholeskyDecomposition::SolveLinearSystemInPlace | ( | VectorType & | b | ) |
Solves linear system Ax=b and outputs result in input b variable.
Definition at line 69 of file animaCholeskyDecomposition.cxx.
Referenced by SolveLinearSystem(), and anima::NNLSOptimizer::StartOptimization().
void anima::CholeskyDecomposition::Update | ( | const VectorType & | x | ) |
Update decomposition with x so that LDL matches A + x x^T.
Definition at line 87 of file animaCholeskyDecomposition.cxx.