ANIMA  4.0
animaCholeskyDecomposition.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <vnl/vnl_matrix.h>
4 #include <vnl/vnl_diag_matrix.h>
5 #include <vnl/vnl_vector.h>
6 
7 #include "AnimaOptimizersExport.h"
8 
9 namespace anima
10 {
11 
17 class ANIMAOPTIMIZERS_EXPORT CholeskyDecomposition
18 {
19 public:
20  typedef vnl_matrix<double> MatrixType;
21  typedef vnl_diag_matrix<double> DiagonalType;
22  typedef vnl_vector<double> VectorType;
23 
25 
26  CholeskyDecomposition(const unsigned int matrixDimension, const double epsilon)
27  {
28  m_InputMatrix.set_size(matrixDimension, matrixDimension);
29  m_InputMatrix.fill(0.0);
30  m_InputMatrix.fill_diagonal(epsilon);
31  m_MatrixSize = matrixDimension;
32  m_LMatrix.set_size(matrixDimension, matrixDimension);
33  m_LMatrix.set_identity();
34  m_DMatrix.set_size(matrixDimension);
35  m_DMatrix.fill_diagonal(epsilon);
36  }
37 
38  CholeskyDecomposition(const MatrixType &val)
39  {
40  this->SetInputMatrix(val);
41  }
42 
44  void PerformDecomposition();
45 
47  void Recompose();
48 
50  VectorType &SolveLinearSystem(const VectorType &b);
51 
53  void SolveLinearSystemInPlace(VectorType &b);
54 
56  void Update(const VectorType &x);
57 
59  void SetInputMatrix(const MatrixType &matrix);
60 
62  double GetConditionNumber();
63 
65  MatrixType &GetInputMatrix() {return m_InputMatrix;}
66 
67  void SetLMatrix(const MatrixType &val) {m_LMatrix = val;}
68  MatrixType &GetLMatrix() {return m_LMatrix;}
69  void SetDMatrix(const DiagonalType &val) {m_DMatrix = val;}
70  DiagonalType &GetDMatrix() {return m_DMatrix;}
71 
72 private:
73  MatrixType m_InputMatrix, m_LMatrix;
74  DiagonalType m_DMatrix;
75  unsigned int m_MatrixSize;
76 
77  VectorType m_ProblemSolution;
78  VectorType m_WorkVector;
79 };
80 
81 } // end of namespace anima
void SetLMatrix(const MatrixType &val)
CholeskyDecomposition(const MatrixType &val)
MatrixType & GetInputMatrix()
Get input matrix, useful if decomposition updated through rank one method.
vnl_diag_matrix< double > DiagonalType
void SetDMatrix(const DiagonalType &val)
CholeskyDecomposition(const unsigned int matrixDimension, const double epsilon)
Cholesky decomposition: decomposes a symmetric matrix A in the form L D L^T, where L is lower triangu...