SimplicialCholesky.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2012 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_SIMPLICIAL_CHOLESKY_H
11 #define EIGEN_SIMPLICIAL_CHOLESKY_H
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
19 
20 namespace internal {
21 template <typename CholMatrixType, typename InputMatrixType>
23  typedef CholMatrixType const* ConstCholMatrixPtr;
24  static void run(const InputMatrixType& input, ConstCholMatrixPtr& pmat, CholMatrixType& tmp) {
25  tmp = input;
26  pmat = &tmp;
27  }
28 };
29 
30 template <typename MatrixType>
32  typedef MatrixType const* ConstMatrixPtr;
33  static void run(const MatrixType& input, ConstMatrixPtr& pmat, MatrixType& /*tmp*/) { pmat = &input; }
34 };
35 } // end namespace internal
36 
50 template <typename Derived>
51 class SimplicialCholeskyBase : public SparseSolverBase<Derived> {
54 
55  public:
59  typedef typename MatrixType::Scalar Scalar;
62  typedef typename MatrixType::StorageIndex StorageIndex;
67 
68  enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
69 
70  public:
71  using Base::derived;
72 
76 
79  derived().compute(matrix);
80  }
81 
83 
84  Derived& derived() { return *static_cast<Derived*>(this); }
85  const Derived& derived() const { return *static_cast<const Derived*>(this); }
86 
87  inline Index cols() const { return m_matrix.cols(); }
88  inline Index rows() const { return m_matrix.rows(); }
89 
96  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
97  return m_info;
98  }
99 
103 
107 
118  Derived& setShift(const DiagonalScalar& offset, const DiagonalScalar& scale = 1) {
119  m_shiftOffset = offset;
120  m_shiftScale = scale;
121  return derived();
122  }
123 
124 #ifndef EIGEN_PARSED_BY_DOXYGEN
126  template <typename Stream>
127  void dumpMemory(Stream& s) {
128  int total = 0;
129  s << " L: "
130  << ((total += (m_matrix.cols() + 1) * sizeof(int) + m_matrix.nonZeros() * (sizeof(int) + sizeof(Scalar))) >> 20)
131  << "Mb"
132  << "\n";
133  s << " diag: " << ((total += m_diag.size() * sizeof(Scalar)) >> 20) << "Mb"
134  << "\n";
135  s << " tree: " << ((total += m_parent.size() * sizeof(int)) >> 20) << "Mb"
136  << "\n";
137  s << " nonzeros: " << ((total += m_workSpace.size() * sizeof(int)) >> 20) << "Mb"
138  << "\n";
139  s << " perm: " << ((total += m_P.size() * sizeof(int)) >> 20) << "Mb"
140  << "\n";
141  s << " perm^-1: " << ((total += m_Pinv.size() * sizeof(int)) >> 20) << "Mb"
142  << "\n";
143  s << " TOTAL: " << (total >> 20) << "Mb"
144  << "\n";
145  }
146 
148  template <typename Rhs, typename Dest>
149  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
151  "The decomposition is not in a valid state for solving, you must first call either compute() or "
152  "symbolic()/numeric()");
153  eigen_assert(m_matrix.rows() == b.rows());
154 
155  if (m_info != Success) return;
156 
157  if (m_P.size() > 0)
158  dest = m_P * b;
159  else
160  dest = b;
161 
162  if (m_matrix.nonZeros() > 0) // otherwise L==I
163  derived().matrixL().solveInPlace(dest);
164 
165  if (m_diag.size() > 0) dest = m_diag.asDiagonal().inverse() * dest;
166 
167  if (m_matrix.nonZeros() > 0) // otherwise U==I
168  derived().matrixU().solveInPlace(dest);
169 
170  if (m_P.size() > 0) dest = m_Pinv * dest;
171  }
172 
173  template <typename Rhs, typename Dest>
176  }
177 
178 #endif // EIGEN_PARSED_BY_DOXYGEN
179 
180  protected:
182  template <bool DoLDLT, bool NonHermitian>
183  void compute(const MatrixType& matrix) {
184  eigen_assert(matrix.rows() == matrix.cols());
185  Index size = matrix.cols();
187  ConstCholMatrixPtr pmat;
188  ordering<NonHermitian>(matrix, pmat, tmp);
189  analyzePattern_preordered(*pmat, DoLDLT);
190  factorize_preordered<DoLDLT, NonHermitian>(*pmat);
191  }
192 
193  template <bool DoLDLT, bool NonHermitian>
194  void factorize(const MatrixType& a) {
195  eigen_assert(a.rows() == a.cols());
196  Index size = a.cols();
198  ConstCholMatrixPtr pmat;
199 
200  if (m_P.size() == 0 && (int(UpLo) & int(Upper)) == Upper) {
201  // If there is no ordering, try to directly use the input matrix without any copy
203  } else {
204  internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, tmp, m_P.indices().data());
205  pmat = &tmp;
206  }
207 
208  factorize_preordered<DoLDLT, NonHermitian>(*pmat);
209  }
210 
211  template <bool DoLDLT, bool NonHermitian>
213 
214  template <bool DoLDLT, bool NonHermitian>
215  void analyzePattern(const MatrixType& a) {
216  eigen_assert(a.rows() == a.cols());
217  Index size = a.cols();
219  ConstCholMatrixPtr pmat;
220  ordering<NonHermitian>(a, pmat, tmp);
221  analyzePattern_preordered(*pmat, DoLDLT);
222  }
223  void analyzePattern_preordered(const CholMatrixType& a, bool doLDLT);
224 
225  template <bool NonHermitian>
227 
230 
232  struct keep_diag {
233  inline bool operator()(const Index& row, const Index& col, const Scalar&) const { return row != col; }
234  };
235 
239 
241  VectorType m_diag; // the diagonal coefficients (LDLT mode)
242  VectorI m_parent; // elimination tree
246 
249 };
250 
251 template <typename MatrixType_, int UpLo_ = Lower,
253 class SimplicialLLT;
254 template <typename MatrixType_, int UpLo_ = Lower,
256 class SimplicialLDLT;
257 template <typename MatrixType_, int UpLo_ = Lower,
260 template <typename MatrixType_, int UpLo_ = Lower,
263 template <typename MatrixType_, int UpLo_ = Lower,
265 class SimplicialCholesky;
266 
267 namespace internal {
268 
269 template <typename MatrixType_, int UpLo_, typename Ordering_>
270 struct traits<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
271  typedef MatrixType_ MatrixType;
272  typedef Ordering_ OrderingType;
273  enum { UpLo = UpLo_ };
274  typedef typename MatrixType::Scalar Scalar;
276  typedef typename MatrixType::StorageIndex StorageIndex;
280  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
281  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.adjoint()); }
282  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
283  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
284 };
285 
286 template <typename MatrixType_, int UpLo_, typename Ordering_>
287 struct traits<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
288  typedef MatrixType_ MatrixType;
289  typedef Ordering_ OrderingType;
290  enum { UpLo = UpLo_ };
291  typedef typename MatrixType::Scalar Scalar;
293  typedef typename MatrixType::StorageIndex StorageIndex;
297  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
298  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.adjoint()); }
299  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
300  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
301 };
302 
303 template <typename MatrixType_, int UpLo_, typename Ordering_>
304 struct traits<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
305  typedef MatrixType_ MatrixType;
306  typedef Ordering_ OrderingType;
307  enum { UpLo = UpLo_ };
308  typedef typename MatrixType::Scalar Scalar;
310  typedef typename MatrixType::StorageIndex StorageIndex;
314  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
315  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.transpose()); }
316  static inline DiagonalScalar getDiag(Scalar x) { return x; }
317  static inline Scalar getSymm(Scalar x) { return x; }
318 };
319 
320 template <typename MatrixType_, int UpLo_, typename Ordering_>
321 struct traits<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
322  typedef MatrixType_ MatrixType;
323  typedef Ordering_ OrderingType;
324  enum { UpLo = UpLo_ };
325  typedef typename MatrixType::Scalar Scalar;
327  typedef typename MatrixType::StorageIndex StorageIndex;
331  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
332  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.transpose()); }
333  static inline DiagonalScalar getDiag(Scalar x) { return x; }
334  static inline Scalar getSymm(Scalar x) { return x; }
335 };
336 
337 template <typename MatrixType_, int UpLo_, typename Ordering_>
338 struct traits<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
339  typedef MatrixType_ MatrixType;
340  typedef Ordering_ OrderingType;
341  enum { UpLo = UpLo_ };
342  typedef typename MatrixType::Scalar Scalar;
344  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
345  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
346 };
347 
348 } // namespace internal
349 
370 template <typename MatrixType_, int UpLo_, typename Ordering_>
371 class SimplicialLLT : public SimplicialCholeskyBase<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
372  public:
373  typedef MatrixType_ MatrixType;
374  enum { UpLo = UpLo_ };
376  typedef typename MatrixType::Scalar Scalar;
378  typedef typename MatrixType::StorageIndex StorageIndex;
382  typedef typename Traits::MatrixL MatrixL;
383  typedef typename Traits::MatrixU MatrixU;
384 
385  public:
389  explicit SimplicialLLT(const MatrixType& matrix) : Base(matrix) {}
390 
392  inline const MatrixL matrixL() const {
393  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
394  return Traits::getL(Base::m_matrix);
395  }
396 
398  inline const MatrixU matrixU() const {
399  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
400  return Traits::getU(Base::m_matrix);
401  }
402 
405  Base::template compute<false, false>(matrix);
406  return *this;
407  }
408 
415  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<false, false>(a); }
416 
423  void factorize(const MatrixType& a) { Base::template factorize<false, false>(a); }
424 
426  Scalar determinant() const {
427  Scalar detL = Base::m_matrix.diagonal().prod();
428  return numext::abs2(detL);
429  }
430 };
431 
452 template <typename MatrixType_, int UpLo_, typename Ordering_>
453 class SimplicialLDLT : public SimplicialCholeskyBase<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
454  public:
455  typedef MatrixType_ MatrixType;
456  enum { UpLo = UpLo_ };
458  typedef typename MatrixType::Scalar Scalar;
460  typedef typename MatrixType::StorageIndex StorageIndex;
464  typedef typename Traits::MatrixL MatrixL;
465  typedef typename Traits::MatrixU MatrixU;
466 
467  public:
470 
472  explicit SimplicialLDLT(const MatrixType& matrix) : Base(matrix) {}
473 
475  inline const VectorType vectorD() const {
476  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
477  return Base::m_diag;
478  }
480  inline const MatrixL matrixL() const {
481  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
482  return Traits::getL(Base::m_matrix);
483  }
484 
486  inline const MatrixU matrixU() const {
487  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
488  return Traits::getU(Base::m_matrix);
489  }
490 
493  Base::template compute<true, false>(matrix);
494  return *this;
495  }
496 
503  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, false>(a); }
504 
511  void factorize(const MatrixType& a) { Base::template factorize<true, false>(a); }
512 
514  Scalar determinant() const { return Base::m_diag.prod(); }
515 };
516 
537 template <typename MatrixType_, int UpLo_, typename Ordering_>
539  : public SimplicialCholeskyBase<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
540  public:
541  typedef MatrixType_ MatrixType;
542  enum { UpLo = UpLo_ };
544  typedef typename MatrixType::Scalar Scalar;
546  typedef typename MatrixType::StorageIndex StorageIndex;
550  typedef typename Traits::MatrixL MatrixL;
551  typedef typename Traits::MatrixU MatrixU;
552 
553  public:
556 
559 
561  inline const MatrixL matrixL() const {
562  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
563  return Traits::getL(Base::m_matrix);
564  }
565 
567  inline const MatrixU matrixU() const {
568  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
569  return Traits::getU(Base::m_matrix);
570  }
571 
574  Base::template compute<false, true>(matrix);
575  return *this;
576  }
577 
584  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<false, true>(a); }
585 
592  void factorize(const MatrixType& a) { Base::template factorize<false, true>(a); }
593 
595  Scalar determinant() const {
596  Scalar detL = Base::m_matrix.diagonal().prod();
597  return detL * detL;
598  }
599 };
600 
621 template <typename MatrixType_, int UpLo_, typename Ordering_>
623  : public SimplicialCholeskyBase<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
624  public:
625  typedef MatrixType_ MatrixType;
626  enum { UpLo = UpLo_ };
628  typedef typename MatrixType::Scalar Scalar;
630  typedef typename MatrixType::StorageIndex StorageIndex;
634  typedef typename Traits::MatrixL MatrixL;
635  typedef typename Traits::MatrixU MatrixU;
636 
637  public:
640 
643 
645  inline const VectorType vectorD() const {
646  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
647  return Base::m_diag;
648  }
650  inline const MatrixL matrixL() const {
651  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
652  return Traits::getL(Base::m_matrix);
653  }
654 
656  inline const MatrixU matrixU() const {
657  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
658  return Traits::getU(Base::m_matrix);
659  }
660 
663  Base::template compute<true, true>(matrix);
664  return *this;
665  }
666 
673  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, true>(a); }
674 
681  void factorize(const MatrixType& a) { Base::template factorize<true, true>(a); }
682 
684  Scalar determinant() const { return Base::m_diag.prod(); }
685 };
686 
693 template <typename MatrixType_, int UpLo_, typename Ordering_>
694 class SimplicialCholesky : public SimplicialCholeskyBase<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
695  public:
696  typedef MatrixType_ MatrixType;
697  enum { UpLo = UpLo_ };
699  typedef typename MatrixType::Scalar Scalar;
701  typedef typename MatrixType::StorageIndex StorageIndex;
706 
707  public:
709 
710  explicit SimplicialCholesky(const MatrixType& matrix) : Base(), m_LDLT(true) { compute(matrix); }
711 
713  switch (mode) {
715  m_LDLT = false;
716  break;
718  m_LDLT = true;
719  break;
720  default:
721  break;
722  }
723 
724  return *this;
725  }
726 
727  inline const VectorType vectorD() const {
728  eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
729  return Base::m_diag;
730  }
731  inline const CholMatrixType rawMatrix() const {
732  eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
733  return Base::m_matrix;
734  }
735 
738  if (m_LDLT)
739  Base::template compute<true, false>(matrix);
740  else
741  Base::template compute<false, false>(matrix);
742  return *this;
743  }
744 
751  void analyzePattern(const MatrixType& a) {
752  if (m_LDLT)
753  Base::template analyzePattern<true, false>(a);
754  else
755  Base::template analyzePattern<false, false>(a);
756  }
757 
764  void factorize(const MatrixType& a) {
765  if (m_LDLT)
766  Base::template factorize<true, false>(a);
767  else
768  Base::template factorize<false, false>(a);
769  }
770 
772  template <typename Rhs, typename Dest>
773  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
775  "The decomposition is not in a valid state for solving, you must first call either compute() or "
776  "symbolic()/numeric()");
777  eigen_assert(Base::m_matrix.rows() == b.rows());
778 
779  if (Base::m_info != Success) return;
780 
781  if (Base::m_P.size() > 0)
782  dest = Base::m_P * b;
783  else
784  dest = b;
785 
786  if (Base::m_matrix.nonZeros() > 0) // otherwise L==I
787  {
788  if (m_LDLT)
789  LDLTTraits::getL(Base::m_matrix).solveInPlace(dest);
790  else
791  LLTTraits::getL(Base::m_matrix).solveInPlace(dest);
792  }
793 
794  if (Base::m_diag.size() > 0) dest = Base::m_diag.real().asDiagonal().inverse() * dest;
795 
796  if (Base::m_matrix.nonZeros() > 0) // otherwise I==I
797  {
798  if (m_LDLT)
799  LDLTTraits::getU(Base::m_matrix).solveInPlace(dest);
800  else
801  LLTTraits::getU(Base::m_matrix).solveInPlace(dest);
802  }
803 
804  if (Base::m_P.size() > 0) dest = Base::m_Pinv * dest;
805  }
806 
808  template <typename Rhs, typename Dest>
811  }
812 
813  Scalar determinant() const {
814  if (m_LDLT) {
815  return Base::m_diag.prod();
816  } else {
818  return numext::abs2(detL);
819  }
820  }
821 
822  protected:
823  bool m_LDLT;
824 };
825 
826 template <typename Derived>
827 template <bool NonHermitian>
829  eigen_assert(a.rows() == a.cols());
830  const Index size = a.rows();
831  pmat = &ap;
832  // Note that ordering methods compute the inverse permutation
834  {
836  internal::permute_symm_to_fullsymm<UpLo, NonHermitian>(a, C, NULL);
837 
838  OrderingType ordering;
839  ordering(C, m_Pinv);
840  }
841 
842  if (m_Pinv.size() > 0)
843  m_P = m_Pinv.inverse();
844  else
845  m_P.resize(0);
846 
847  ap.resize(size, size);
848  internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, m_P.indices().data());
849  } else {
850  m_Pinv.resize(0);
851  m_P.resize(0);
852  if (int(UpLo) == int(Lower) || MatrixType::IsRowMajor) {
853  // we have to transpose the lower part to to the upper one
854  ap.resize(size, size);
855  internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, NULL);
856  } else
858  }
859 }
860 
861 } // end namespace Eigen
862 
863 #endif // EIGEN_SIMPLICIAL_CHOLESKY_H
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
#define eigen_assert(x)
Definition: Macros.h:910
m col(1)
m row(1)
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Definition: Ordering.h:50
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:68
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Definition: Ordering.h:89
EIGEN_DEVICE_FUNC Index size() const
Definition: PermutationMatrix.h:96
const IndicesType & indices() const
Definition: PermutationMatrix.h:334
A base class for direct sparse Cholesky factorizations.
Definition: SimplicialCholesky.h:51
Matrix< StorageIndex, Dynamic, 1 > VectorI
Definition: SimplicialCholesky.h:66
SparseSolverBase< Derived > Base
Definition: SimplicialCholesky.h:52
SimplicialCholeskyBase()
Definition: SimplicialCholesky.h:74
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:194
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:60
VectorI m_parent
Definition: SimplicialCholesky.h:242
Derived & derived()
Definition: SimplicialCholesky.h:84
CholMatrixType const * ConstCholMatrixPtr
Definition: SimplicialCholesky.h:64
void compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:183
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition: SimplicialCholesky.h:149
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SimplicialCholesky.h:95
bool m_analysisIsOk
Definition: SimplicialCholesky.h:238
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:59
void dumpMemory(Stream &s)
Definition: SimplicialCholesky.h:127
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:63
CholMatrixType m_matrix
Definition: SimplicialCholesky.h:240
void _solve_impl(const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
Definition: SimplicialCholesky.h:174
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:215
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv() const
Definition: SimplicialCholesky.h:106
@ UpLo
Definition: SimplicialCholesky.h:58
DiagonalScalar m_shiftScale
Definition: SimplicialCholesky.h:248
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:62
SimplicialCholeskyBase(const MatrixType &matrix)
Definition: SimplicialCholesky.h:77
DiagonalScalar m_shiftOffset
Definition: SimplicialCholesky.h:247
Index cols() const
Definition: SimplicialCholesky.h:87
Derived & setShift(const DiagonalScalar &offset, const DiagonalScalar &scale=1)
Definition: SimplicialCholesky.h:118
@ ColsAtCompileTime
Definition: SimplicialCholesky.h:68
@ MaxColsAtCompileTime
Definition: SimplicialCholesky.h:68
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP() const
Definition: SimplicialCholesky.h:102
internal::traits< Derived >::DiagonalScalar DiagonalScalar
Definition: SimplicialCholesky.h:61
VectorType m_diag
Definition: SimplicialCholesky.h:241
Index rows() const
Definition: SimplicialCholesky.h:88
const Derived & derived() const
Definition: SimplicialCholesky.h:85
Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:229
void analyzePattern_preordered(const CholMatrixType &a, bool doLDLT)
Definition: SimplicialCholesky_impl.h:280
DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:228
internal::traits< Derived >::MatrixType MatrixType
Definition: SimplicialCholesky.h:56
bool m_isInitialized
Definition: SparseSolverBase.h:110
void ordering(const MatrixType &a, ConstCholMatrixPtr &pmat, CholMatrixType &ap)
Definition: SimplicialCholesky.h:828
~SimplicialCholeskyBase()
Definition: SimplicialCholesky.h:82
bool m_factorizationIsOk
Definition: SimplicialCholesky.h:237
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:65
void factorize_preordered(const CholMatrixType &a)
Definition: SimplicialCholesky_impl.h:296
VectorI m_workSpace
Definition: SimplicialCholesky.h:243
ComputationInfo m_info
Definition: SimplicialCholesky.h:236
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition: SimplicialCholesky.h:245
internal::traits< Derived >::OrderingType OrderingType
Definition: SimplicialCholesky.h:57
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
Definition: SimplicialCholesky.h:244
Definition: SimplicialCholesky.h:694
SimplicialCholesky()
Definition: SimplicialCholesky.h:708
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition: SimplicialCholesky.h:773
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:699
bool m_LDLT
Definition: SimplicialCholesky.h:823
internal::traits< SimplicialLLT< MatrixType, UpLo > > LLTTraits
Definition: SimplicialCholesky.h:705
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:700
SimplicialCholesky & compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:737
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:764
SimplicialCholesky & setMode(SimplicialCholeskyMode mode)
Definition: SimplicialCholesky.h:712
void _solve_impl(const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
Definition: SimplicialCholesky.h:809
const VectorType vectorD() const
Definition: SimplicialCholesky.h:727
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:703
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:751
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:696
internal::traits< SimplicialLDLT< MatrixType, UpLo > > LDLTTraits
Definition: SimplicialCholesky.h:704
Scalar determinant() const
Definition: SimplicialCholesky.h:813
@ UpLo
Definition: SimplicialCholesky.h:697
SimplicialCholesky(const MatrixType &matrix)
Definition: SimplicialCholesky.h:710
SimplicialCholeskyBase< SimplicialCholesky > Base
Definition: SimplicialCholesky.h:698
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:701
const CholMatrixType rawMatrix() const
Definition: SimplicialCholesky.h:731
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:702
A direct sparse LDLT Cholesky factorizations without square root.
Definition: SimplicialCholesky.h:453
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:459
Traits::MatrixL MatrixL
Definition: SimplicialCholesky.h:464
Scalar determinant() const
Definition: SimplicialCholesky.h:514
const MatrixL matrixL() const
Definition: SimplicialCholesky.h:480
SimplicialLDLT(const MatrixType &matrix)
Definition: SimplicialCholesky.h:472
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:458
const MatrixU matrixU() const
Definition: SimplicialCholesky.h:486
@ UpLo
Definition: SimplicialCholesky.h:456
SimplicialCholeskyBase< SimplicialLDLT > Base
Definition: SimplicialCholesky.h:457
SimplicialLDLT & compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:492
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:461
internal::traits< SimplicialLDLT > Traits
Definition: SimplicialCholesky.h:463
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:503
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:511
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:455
const VectorType vectorD() const
Definition: SimplicialCholesky.h:475
Traits::MatrixU MatrixU
Definition: SimplicialCholesky.h:465
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:462
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:460
SimplicialLDLT()
Definition: SimplicialCholesky.h:469
A direct sparse LLT Cholesky factorizations.
Definition: SimplicialCholesky.h:371
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:376
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:378
const MatrixL matrixL() const
Definition: SimplicialCholesky.h:392
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:380
Scalar determinant() const
Definition: SimplicialCholesky.h:426
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:373
SimplicialCholeskyBase< SimplicialLLT > Base
Definition: SimplicialCholesky.h:375
Traits::MatrixL MatrixL
Definition: SimplicialCholesky.h:382
@ UpLo
Definition: SimplicialCholesky.h:374
internal::traits< SimplicialLLT > Traits
Definition: SimplicialCholesky.h:381
SimplicialLLT & compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:404
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:377
const MatrixU matrixU() const
Definition: SimplicialCholesky.h:398
SimplicialLLT()
Definition: SimplicialCholesky.h:387
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:415
SparseMatrix< Scalar, ColMajor, Index > CholMatrixType
Definition: SimplicialCholesky.h:379
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:423
Traits::MatrixU MatrixU
Definition: SimplicialCholesky.h:383
SimplicialLLT(const MatrixType &matrix)
Definition: SimplicialCholesky.h:389
A direct sparse LDLT Cholesky factorizations without square root, for symmetric non-hermitian matrice...
Definition: SimplicialCholesky.h:623
SimplicialNonHermitianLDLT()
Definition: SimplicialCholesky.h:639
Scalar determinant() const
Definition: SimplicialCholesky.h:684
SimplicialCholeskyBase< SimplicialNonHermitianLDLT > Base
Definition: SimplicialCholesky.h:627
const MatrixU matrixU() const
Definition: SimplicialCholesky.h:656
Traits::MatrixU MatrixU
Definition: SimplicialCholesky.h:635
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:625
@ UpLo
Definition: SimplicialCholesky.h:626
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:632
Traits::MatrixL MatrixL
Definition: SimplicialCholesky.h:634
internal::traits< SimplicialNonHermitianLDLT > Traits
Definition: SimplicialCholesky.h:633
const MatrixL matrixL() const
Definition: SimplicialCholesky.h:650
const VectorType vectorD() const
Definition: SimplicialCholesky.h:645
SimplicialNonHermitianLDLT(const MatrixType &matrix)
Definition: SimplicialCholesky.h:642
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:681
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:629
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:630
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:631
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:673
SimplicialNonHermitianLDLT & compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:662
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:628
A direct sparse LLT Cholesky factorizations, for symmetric non-hermitian matrices.
Definition: SimplicialCholesky.h:539
void factorize(const MatrixType &a)
Definition: SimplicialCholesky.h:592
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:546
const MatrixU matrixU() const
Definition: SimplicialCholesky.h:567
Matrix< Scalar, Dynamic, 1 > VectorType
Definition: SimplicialCholesky.h:548
internal::traits< SimplicialNonHermitianLLT > Traits
Definition: SimplicialCholesky.h:549
SimplicialNonHermitianLLT()
Definition: SimplicialCholesky.h:555
Traits::MatrixU MatrixU
Definition: SimplicialCholesky.h:551
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:544
Scalar determinant() const
Definition: SimplicialCholesky.h:595
SimplicialCholeskyBase< SimplicialNonHermitianLLT > Base
Definition: SimplicialCholesky.h:543
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:541
const MatrixL matrixL() const
Definition: SimplicialCholesky.h:561
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:547
Traits::MatrixL MatrixL
Definition: SimplicialCholesky.h:550
MatrixType::RealScalar RealScalar
Definition: SimplicialCholesky.h:545
void analyzePattern(const MatrixType &a)
Definition: SimplicialCholesky.h:584
SimplicialNonHermitianLLT(const MatrixType &matrix)
Definition: SimplicialCholesky.h:558
SimplicialNonHermitianLLT & compute(const MatrixType &matrix)
Definition: SimplicialCholesky.h:573
@ UpLo
Definition: SimplicialCholesky.h:542
Base class of any sparse matrices or sparse expressions.
Definition: SparseMatrixBase.h:30
Index nonZeros() const
Definition: SparseCompressedBase.h:64
const ConstDiagonalReturnType diagonal() const
Definition: SparseMatrix.h:757
Index cols() const
Definition: SparseMatrix.h:161
Index rows() const
Definition: SparseMatrix.h:159
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
bool m_isInitialized
Definition: SparseSolverBase.h:110
Derived & derived()
Definition: SparseSolverBase.h:76
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:167
Definition: matrices.h:74
Concept for reading and writing characters.
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
float real
Definition: datatypes.h:10
ComputationInfo
Definition: Constants.h:438
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ Success
Definition: Constants.h:440
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
Scalar * ap
Definition: level2_cplx_impl.h:161
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
std::enable_if_t< Rhs::ColsAtCompileTime !=1 &&Dest::ColsAtCompileTime !=1 > solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition: SparseSolverBase.h:25
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
SimplicialCholeskyMode
Definition: SimplicialCholesky.h:18
@ SimplicialCholeskyLDLT
Definition: SimplicialCholesky.h:18
@ SimplicialCholeskyLLT
Definition: SimplicialCholesky.h:18
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: SimplicialCholesky.h:232
bool operator()(const Index &row, const Index &col, const Scalar &) const
Definition: SimplicialCholesky.h:233
Definition: Meta.h:205
static void run(const MatrixType &input, ConstMatrixPtr &pmat, MatrixType &)
Definition: SimplicialCholesky.h:33
MatrixType const * ConstMatrixPtr
Definition: SimplicialCholesky.h:32
Definition: SimplicialCholesky.h:22
CholMatrixType const * ConstCholMatrixPtr
Definition: SimplicialCholesky.h:23
static void run(const InputMatrixType &input, ConstCholMatrixPtr &pmat, CholMatrixType &tmp)
Definition: SimplicialCholesky.h:24
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:342
static Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:345
static DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:344
MatrixType::RealScalar DiagonalScalar
Definition: SimplicialCholesky.h:343
MatrixType::RealScalar DiagonalScalar
Definition: SimplicialCholesky.h:292
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:288
TriangularView< const CholMatrixType, Eigen::UnitLower > MatrixL
Definition: SimplicialCholesky.h:295
TriangularView< const typename CholMatrixType::AdjointReturnType, Eigen::UnitUpper > MatrixU
Definition: SimplicialCholesky.h:296
static MatrixU getU(const CholMatrixType &m)
Definition: SimplicialCholesky.h:298
static DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:299
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:291
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:294
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:293
static Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:300
static MatrixL getL(const CholMatrixType &m)
Definition: SimplicialCholesky.h:297
static MatrixU getU(const CholMatrixType &m)
Definition: SimplicialCholesky.h:281
TriangularView< const typename CholMatrixType::AdjointReturnType, Eigen::Upper > MatrixU
Definition: SimplicialCholesky.h:279
static DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:282
MatrixType_ MatrixType
Definition: SimplicialCholesky.h:271
static Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:283
Ordering_ OrderingType
Definition: SimplicialCholesky.h:272
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:276
static MatrixL getL(const CholMatrixType &m)
Definition: SimplicialCholesky.h:280
MatrixType::RealScalar DiagonalScalar
Definition: SimplicialCholesky.h:275
TriangularView< const CholMatrixType, Eigen::Lower > MatrixL
Definition: SimplicialCholesky.h:278
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:274
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:277
static MatrixU getU(const CholMatrixType &m)
Definition: SimplicialCholesky.h:332
TriangularView< const CholMatrixType, Eigen::UnitLower > MatrixL
Definition: SimplicialCholesky.h:329
static Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:334
TriangularView< const typename CholMatrixType::ConstTransposeReturnType, Eigen::UnitUpper > MatrixU
Definition: SimplicialCholesky.h:330
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:327
MatrixType::Scalar DiagonalScalar
Definition: SimplicialCholesky.h:326
static MatrixL getL(const CholMatrixType &m)
Definition: SimplicialCholesky.h:331
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:328
static DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:333
static DiagonalScalar getDiag(Scalar x)
Definition: SimplicialCholesky.h:316
MatrixType::Scalar DiagonalScalar
Definition: SimplicialCholesky.h:309
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:310
static Scalar getSymm(Scalar x)
Definition: SimplicialCholesky.h:317
static MatrixL getL(const CholMatrixType &m)
Definition: SimplicialCholesky.h:314
static MatrixU getU(const CholMatrixType &m)
Definition: SimplicialCholesky.h:315
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition: SimplicialCholesky.h:311
TriangularView< const CholMatrixType, Eigen::Lower > MatrixL
Definition: SimplicialCholesky.h:312
TriangularView< const typename CholMatrixType::ConstTransposeReturnType, Eigen::Upper > MatrixU
Definition: SimplicialCholesky.h:313
Definition: ForwardDeclarations.h:21