lsst.jointcal  14.0-26-gc4bc114+1
Eigenstuff.h
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 #ifndef LSST_JOINTCAL_EIGENSTUFF_H
3 #define LSST_JOINTCAL_EIGENSTUFF_H
4 
5 #include "Eigen/CholmodSupport" // to switch to cholmod
6 
7 #include "Eigen/Core"
8 
9 typedef Eigen::Matrix<double, Eigen::Dynamic, 2> MatrixX2d;
10 
11 typedef Eigen::SparseMatrix<double> SpMat;
12 
13 /* Cholesky factorization class using cholmod, with the small-rank update capability.
14  *
15  * Class derived from Eigen's CholmodBase, to add the factorization
16  * update capability to the interface. Besides this addition, it
17  * behaves the same way as Eigen's native Cholesky factorization
18  * classes. It relies on the simplicial LDLt factorization.
19  *
20  * @Seealso Eigen::CholmodSimplicialLDLT
21  */
22 template <typename _MatrixType, int _UpLo = Eigen::Lower>
24  : public Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT2<_MatrixType, _UpLo>> {
25  typedef Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT2> Base;
26  using Base::m_cholmod;
27 
28 public:
29  typedef _MatrixType MatrixType;
30  typedef typename MatrixType::Index Index;
31  typedef typename MatrixType::RealScalar RealScalar;
32 
33  CholmodSimplicialLDLT2() : Base() { init(); }
34 
35  CholmodSimplicialLDLT2(MatrixType const &matrix) : Base() {
36  init();
37  this->compute(matrix);
38  }
39 
40  // this routine is the one we added
41  int update(SpMat const &H, bool UpOrDown) {
42  // check size
43  Index const size = Base::m_cholmodFactor->n;
44  EIGEN_UNUSED_VARIABLE(size);
45  eigen_assert(size == H.rows());
46 
47  cholmod_sparse C_cs = viewAsCholmod(H);
48  /* We have to apply the magic permutation to the update matrix,
49  read page 117 of Cholmod UserGuide.pdf */
50  cholmod_sparse *C_cs_perm =
51  cholmod_submatrix(&C_cs, (int *)Base::m_cholmodFactor->Perm, Base::m_cholmodFactor->n,
52  nullptr, -1, true, true, &this->cholmod());
53  assert(C_cs_perm);
54  int ret = cholmod_updown(UpOrDown, C_cs_perm, Base::m_cholmodFactor, &this->cholmod());
55  cholmod_free_sparse(&C_cs_perm, &this->cholmod());
56  assert(ret != 0);
57  return ret;
58  }
59 
60 protected:
61  void init() {
62  m_cholmod.final_asis = 1;
63  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
64  // In CholmodBase::CholmodBase(), the following statement is missing in
65  // SuiteSparse 3.2.0.8. Fixed in 3.2.7
66  Base::m_shiftOffset[0] = Base::m_shiftOffset[1] = RealScalar(0.0);
67  }
68 };
69 
70 #endif // LSST_JOINTCAL_EIGENSTUFF_H
int update(SpMat const &H, bool UpOrDown)
Definition: Eigenstuff.h:41
_MatrixType MatrixType
Definition: Eigenstuff.h:29
CholmodSimplicialLDLT2(MatrixType const &matrix)
Definition: Eigenstuff.h:35
MatrixType::RealScalar RealScalar
Definition: Eigenstuff.h:31
Eigen::SparseMatrix< double > SpMat
Definition: Eigenstuff.h:11
MatrixType::Index Index
Definition: Eigenstuff.h:30
Eigen::Matrix< double, Eigen::Dynamic, 2 > MatrixX2d
Definition: Eigenstuff.h:9