Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
SimplicialCholesky.h
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
16namespace Eigen {
17
18enum SimplicialCholeskyMode { SimplicialCholeskyLLT, SimplicialCholeskyLDLT };
19
20namespace internal {
21template <typename CholMatrixType, typename InputMatrixType>
22struct simplicial_cholesky_grab_input {
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
30template <typename MatrixType>
31struct simplicial_cholesky_grab_input<MatrixType, 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
50template <typename Derived>
52 typedef SparseSolverBase<Derived> Base;
53 using Base::m_isInitialized;
54
55 public:
56 typedef typename internal::traits<Derived>::MatrixType MatrixType;
57 typedef typename internal::traits<Derived>::OrderingType OrderingType;
58 enum { UpLo = internal::traits<Derived>::UpLo };
59 typedef typename MatrixType::Scalar Scalar;
60 typedef typename MatrixType::RealScalar RealScalar;
61 typedef typename internal::traits<Derived>::DiagonalScalar DiagonalScalar;
62 typedef typename MatrixType::StorageIndex StorageIndex;
64 typedef CholMatrixType const* ConstCholMatrixPtr;
65 typedef Matrix<Scalar, Dynamic, 1> VectorType;
67
68 enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
69
70 public:
71 using Base::derived;
72
75 : m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {}
76
77 explicit SimplicialCholeskyBase(const MatrixType& matrix)
78 : m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {
79 derived().compute(matrix);
80 }
81
82 ~SimplicialCholeskyBase() {}
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 {
150 eigen_assert(m_factorizationIsOk &&
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>
174 void _solve_impl(const SparseMatrixBase<Rhs>& b, SparseMatrixBase<Dest>& dest) const {
175 internal::solve_sparse_through_dense_panels(derived(), b, 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();
186 CholMatrixType tmp(size, size);
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();
197 CholMatrixType tmp(size, size);
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
202 internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, tmp);
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>
212 void factorize_preordered(const CholMatrixType& a);
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();
218 CholMatrixType tmp(size, size);
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>
226 void ordering(const MatrixType& a, ConstCholMatrixPtr& pmat, CholMatrixType& ap);
227
228 inline DiagonalScalar getDiag(Scalar x) { return internal::traits<Derived>::getDiag(x); }
229 inline Scalar getSymm(Scalar x) { return internal::traits<Derived>::getSymm(x); }
230
232 struct keep_diag {
233 inline bool operator()(const Index& row, const Index& col, const Scalar&) const { return row != col; }
234 };
235
236 mutable ComputationInfo m_info;
237 bool m_factorizationIsOk;
238 bool m_analysisIsOk;
239
240 CholMatrixType m_matrix;
241 VectorType m_diag; // the diagonal coefficients (LDLT mode)
242 VectorI m_parent; // elimination tree
243 VectorI m_workSpace;
245 PermutationMatrix<Dynamic, Dynamic, StorageIndex> m_Pinv; // the inverse permutation
246
247 DiagonalScalar m_shiftOffset;
248 DiagonalScalar m_shiftScale;
249};
250
251template <typename MatrixType_, int UpLo_ = Lower,
253class SimplicialLLT;
254template <typename MatrixType_, int UpLo_ = Lower,
256class SimplicialLDLT;
257template <typename MatrixType_, int UpLo_ = Lower,
260template <typename MatrixType_, int UpLo_ = Lower,
263template <typename MatrixType_, int UpLo_ = Lower,
266
267namespace internal {
268
269template <typename MatrixType_, int UpLo_, typename Ordering_>
270struct traits<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
271 typedef MatrixType_ MatrixType;
272 typedef Ordering_ OrderingType;
273 enum { UpLo = UpLo_ };
274 typedef typename MatrixType::Scalar Scalar;
275 typedef typename MatrixType::RealScalar DiagonalScalar;
276 typedef typename MatrixType::StorageIndex StorageIndex;
277 typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
278 typedef TriangularView<const CholMatrixType, Eigen::Lower> MatrixL;
279 typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::Upper> MatrixU;
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
286template <typename MatrixType_, int UpLo_, typename Ordering_>
287struct traits<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
288 typedef MatrixType_ MatrixType;
289 typedef Ordering_ OrderingType;
290 enum { UpLo = UpLo_ };
291 typedef typename MatrixType::Scalar Scalar;
292 typedef typename MatrixType::RealScalar DiagonalScalar;
293 typedef typename MatrixType::StorageIndex StorageIndex;
294 typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
295 typedef TriangularView<const CholMatrixType, Eigen::UnitLower> MatrixL;
296 typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::UnitUpper> MatrixU;
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
303template <typename MatrixType_, int UpLo_, typename Ordering_>
304struct traits<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
305 typedef MatrixType_ MatrixType;
306 typedef Ordering_ OrderingType;
307 enum { UpLo = UpLo_ };
308 typedef typename MatrixType::Scalar Scalar;
309 typedef typename MatrixType::Scalar DiagonalScalar;
310 typedef typename MatrixType::StorageIndex StorageIndex;
311 typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
312 typedef TriangularView<const CholMatrixType, Eigen::Lower> MatrixL;
313 typedef TriangularView<const typename CholMatrixType::ConstTransposeReturnType, Eigen::Upper> MatrixU;
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
320template <typename MatrixType_, int UpLo_, typename Ordering_>
321struct traits<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
322 typedef MatrixType_ MatrixType;
323 typedef Ordering_ OrderingType;
324 enum { UpLo = UpLo_ };
325 typedef typename MatrixType::Scalar Scalar;
326 typedef typename MatrixType::Scalar DiagonalScalar;
327 typedef typename MatrixType::StorageIndex StorageIndex;
328 typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
329 typedef TriangularView<const CholMatrixType, Eigen::UnitLower> MatrixL;
330 typedef TriangularView<const typename CholMatrixType::ConstTransposeReturnType, Eigen::UnitUpper> MatrixU;
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
337template <typename MatrixType_, int UpLo_, typename Ordering_>
338struct traits<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
339 typedef MatrixType_ MatrixType;
340 typedef Ordering_ OrderingType;
341 enum { UpLo = UpLo_ };
342 typedef typename MatrixType::Scalar Scalar;
343 typedef typename MatrixType::RealScalar DiagonalScalar;
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
370template <typename MatrixType_, int UpLo_, typename Ordering_>
371class SimplicialLLT : public SimplicialCholeskyBase<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
372 public:
373 typedef MatrixType_ MatrixType;
374 enum { UpLo = UpLo_ };
376 typedef typename MatrixType::Scalar Scalar;
377 typedef typename MatrixType::RealScalar RealScalar;
378 typedef typename MatrixType::StorageIndex StorageIndex;
379 typedef SparseMatrix<Scalar, ColMajor, Index> CholMatrixType;
380 typedef Matrix<Scalar, Dynamic, 1> VectorType;
381 typedef internal::traits<SimplicialLLT> Traits;
382 typedef typename Traits::MatrixL MatrixL;
383 typedef typename Traits::MatrixU MatrixU;
384
385 public:
387 SimplicialLLT() : Base() {}
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
404 SimplicialLLT& compute(const MatrixType& matrix) {
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
424 void factorize(const MatrixType& a) { Base::template factorize<false, false>(a); }
425
427 Scalar determinant() const {
428 Scalar detL = Base::m_matrix.diagonal().prod();
429 return numext::abs2(detL);
430 }
431};
432
453template <typename MatrixType_, int UpLo_, typename Ordering_>
454class SimplicialLDLT : public SimplicialCholeskyBase<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
455 public:
456 typedef MatrixType_ MatrixType;
457 enum { UpLo = UpLo_ };
459 typedef typename MatrixType::Scalar Scalar;
460 typedef typename MatrixType::RealScalar RealScalar;
461 typedef typename MatrixType::StorageIndex StorageIndex;
463 typedef Matrix<Scalar, Dynamic, 1> VectorType;
464 typedef internal::traits<SimplicialLDLT> Traits;
465 typedef typename Traits::MatrixL MatrixL;
466 typedef typename Traits::MatrixU MatrixU;
467
468 public:
470 SimplicialLDLT() : Base() {}
471
473 explicit SimplicialLDLT(const MatrixType& matrix) : Base(matrix) {}
474
476 inline const VectorType vectorD() const {
477 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
478 return Base::m_diag;
479 }
480
481 inline const MatrixL matrixL() const {
482 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
483 return Traits::getL(Base::m_matrix);
484 }
485
487 inline const MatrixU matrixU() const {
488 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
489 return Traits::getU(Base::m_matrix);
490 }
491
493 SimplicialLDLT& compute(const MatrixType& matrix) {
494 Base::template compute<true, false>(matrix);
495 return *this;
496 }
497
504 void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, false>(a); }
505
513 void factorize(const MatrixType& a) { Base::template factorize<true, false>(a); }
514
516 Scalar determinant() const { return Base::m_diag.prod(); }
517};
518
539template <typename MatrixType_, int UpLo_, typename Ordering_>
541 : public SimplicialCholeskyBase<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
542 public:
543 typedef MatrixType_ MatrixType;
544 enum { UpLo = UpLo_ };
546 typedef typename MatrixType::Scalar Scalar;
547 typedef typename MatrixType::RealScalar RealScalar;
548 typedef typename MatrixType::StorageIndex StorageIndex;
550 typedef Matrix<Scalar, Dynamic, 1> VectorType;
551 typedef internal::traits<SimplicialNonHermitianLLT> Traits;
552 typedef typename Traits::MatrixL MatrixL;
553 typedef typename Traits::MatrixU MatrixU;
554
555 public:
558
560 explicit SimplicialNonHermitianLLT(const MatrixType& matrix) : Base(matrix) {}
561
563 inline const MatrixL matrixL() const {
564 eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
565 return Traits::getL(Base::m_matrix);
566 }
567
569 inline const MatrixU matrixU() const {
570 eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
571 return Traits::getU(Base::m_matrix);
572 }
573
575 SimplicialNonHermitianLLT& compute(const MatrixType& matrix) {
576 Base::template compute<false, true>(matrix);
577 return *this;
578 }
579
586 void analyzePattern(const MatrixType& a) { Base::template analyzePattern<false, true>(a); }
587
595 void factorize(const MatrixType& a) { Base::template factorize<false, true>(a); }
596
598 Scalar determinant() const {
599 Scalar detL = Base::m_matrix.diagonal().prod();
600 return detL * detL;
601 }
602};
603
624template <typename MatrixType_, int UpLo_, typename Ordering_>
626 : public SimplicialCholeskyBase<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
627 public:
628 typedef MatrixType_ MatrixType;
629 enum { UpLo = UpLo_ };
631 typedef typename MatrixType::Scalar Scalar;
632 typedef typename MatrixType::RealScalar RealScalar;
633 typedef typename MatrixType::StorageIndex StorageIndex;
635 typedef Matrix<Scalar, Dynamic, 1> VectorType;
636 typedef internal::traits<SimplicialNonHermitianLDLT> Traits;
637 typedef typename Traits::MatrixL MatrixL;
638 typedef typename Traits::MatrixU MatrixU;
639
640 public:
643
645 explicit SimplicialNonHermitianLDLT(const MatrixType& matrix) : Base(matrix) {}
646
648 inline const VectorType vectorD() const {
649 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
650 return Base::m_diag;
651 }
652
653 inline const MatrixL matrixL() const {
654 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
655 return Traits::getL(Base::m_matrix);
656 }
657
659 inline const MatrixU matrixU() const {
660 eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
661 return Traits::getU(Base::m_matrix);
662 }
663
665 SimplicialNonHermitianLDLT& compute(const MatrixType& matrix) {
666 Base::template compute<true, true>(matrix);
667 return *this;
668 }
669
676 void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, true>(a); }
677
685 void factorize(const MatrixType& a) { Base::template factorize<true, true>(a); }
686
688 Scalar determinant() const { return Base::m_diag.prod(); }
689};
690
697template <typename MatrixType_, int UpLo_, typename Ordering_>
698class SimplicialCholesky : public SimplicialCholeskyBase<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
699 public:
700 typedef MatrixType_ MatrixType;
701 enum { UpLo = UpLo_ };
703 typedef typename MatrixType::Scalar Scalar;
704 typedef typename MatrixType::RealScalar RealScalar;
705 typedef typename MatrixType::StorageIndex StorageIndex;
707 typedef Matrix<Scalar, Dynamic, 1> VectorType;
708 typedef internal::traits<SimplicialLDLT<MatrixType, UpLo> > LDLTTraits;
709 typedef internal::traits<SimplicialLLT<MatrixType, UpLo> > LLTTraits;
710
711 public:
712 SimplicialCholesky() : Base(), m_LDLT(true) {}
713
714 explicit SimplicialCholesky(const MatrixType& matrix) : Base(), m_LDLT(true) { compute(matrix); }
715
716 SimplicialCholesky& setMode(SimplicialCholeskyMode mode) {
717 switch (mode) {
718 case SimplicialCholeskyLLT:
719 m_LDLT = false;
720 break;
721 case SimplicialCholeskyLDLT:
722 m_LDLT = true;
723 break;
724 default:
725 break;
726 }
727
728 return *this;
729 }
730
731 inline const VectorType vectorD() const {
732 eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
733 return Base::m_diag;
734 }
735 inline const CholMatrixType rawMatrix() const {
736 eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
737 return Base::m_matrix;
738 }
739
741 SimplicialCholesky& compute(const MatrixType& matrix) {
742 if (m_LDLT)
743 Base::template compute<true, false>(matrix);
744 else
745 Base::template compute<false, false>(matrix);
746 return *this;
747 }
748
755 void analyzePattern(const MatrixType& a) {
756 if (m_LDLT)
757 Base::template analyzePattern<true, false>(a);
758 else
759 Base::template analyzePattern<false, false>(a);
760 }
761
769 void factorize(const MatrixType& a) {
770 if (m_LDLT)
771 Base::template factorize<true, false>(a);
772 else
773 Base::template factorize<false, false>(a);
774 }
775
777 template <typename Rhs, typename Dest>
778 void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
779 eigen_assert(Base::m_factorizationIsOk &&
780 "The decomposition is not in a valid state for solving, you must first call either compute() or "
781 "symbolic()/numeric()");
782 eigen_assert(Base::m_matrix.rows() == b.rows());
783
784 if (Base::m_info != Success) return;
785
786 if (Base::m_P.size() > 0)
787 dest = Base::m_P * b;
788 else
789 dest = b;
790
791 if (Base::m_matrix.nonZeros() > 0) // otherwise L==I
792 {
793 if (m_LDLT)
794 LDLTTraits::getL(Base::m_matrix).solveInPlace(dest);
795 else
796 LLTTraits::getL(Base::m_matrix).solveInPlace(dest);
797 }
798
799 if (Base::m_diag.size() > 0) dest = Base::m_diag.real().asDiagonal().inverse() * dest;
800
801 if (Base::m_matrix.nonZeros() > 0) // otherwise I==I
802 {
803 if (m_LDLT)
804 LDLTTraits::getU(Base::m_matrix).solveInPlace(dest);
805 else
806 LLTTraits::getU(Base::m_matrix).solveInPlace(dest);
807 }
808
809 if (Base::m_P.size() > 0) dest = Base::m_Pinv * dest;
810 }
811
813 template <typename Rhs, typename Dest>
814 void _solve_impl(const SparseMatrixBase<Rhs>& b, SparseMatrixBase<Dest>& dest) const {
815 internal::solve_sparse_through_dense_panels(*this, b, dest);
816 }
817
818 Scalar determinant() const {
819 if (m_LDLT) {
820 return Base::m_diag.prod();
821 } else {
822 Scalar detL = Diagonal<const CholMatrixType>(Base::m_matrix).prod();
823 return numext::abs2(detL);
824 }
825 }
826
827 protected:
828 bool m_LDLT;
829};
830
831template <typename Derived>
832template <bool NonHermitian>
833void SimplicialCholeskyBase<Derived>::ordering(const MatrixType& a, ConstCholMatrixPtr& pmat, CholMatrixType& ap) {
834 eigen_assert(a.rows() == a.cols());
835 const Index size = a.rows();
836 pmat = &ap;
837 // Note that ordering methods compute the inverse permutation
838 if (!internal::is_same<OrderingType, NaturalOrdering<StorageIndex> >::value) {
839 {
840 CholMatrixType C;
841 internal::permute_symm_to_fullsymm<UpLo, NonHermitian>(a, C, NULL);
842
843 OrderingType ordering;
844 ordering(C, m_Pinv);
845 }
846
847 if (m_Pinv.size() > 0)
848 m_P = m_Pinv.inverse();
849 else
850 m_P.resize(0);
851
852 ap.resize(size, size);
853 internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, m_P.indices().data());
854 } else {
855 m_Pinv.resize(0);
856 m_P.resize(0);
857 if (int(UpLo) == int(Lower) || MatrixType::IsRowMajor) {
858 // we have to transpose the lower part to to the upper one
859 ap.resize(size, size);
860 internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, NULL);
861 } else
862 internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, ap);
863 }
864}
865
866} // end namespace Eigen
867
868#endif // EIGEN_SIMPLICIAL_CHOLESKY_H
Definition Ordering.h:48
RealReturnType real() const
Definition DenseBase.h:86
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Definition Ordering.h:87
Index size() const
Definition PermutationMatrix.h:96
Permutation matrix.
Definition PermutationMatrix.h:283
const IndicesType & indices() const
Definition PermutationMatrix.h:337
SimplicialCholeskyBase()
Definition SimplicialCholesky.h:74
void compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:183
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SimplicialCholesky.h:95
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP() const
Definition SimplicialCholesky.h:102
Derived & setShift(const DiagonalScalar &offset, const DiagonalScalar &scale=1)
Definition SimplicialCholesky.h:118
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv() const
Definition SimplicialCholesky.h:106
Definition SimplicialCholesky.h:698
void factorize(const MatrixType &a)
Definition SimplicialCholesky.h:769
void analyzePattern(const MatrixType &a)
Definition SimplicialCholesky.h:755
SimplicialCholesky & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:741
A direct sparse LDLT Cholesky factorizations without square root.
Definition SimplicialCholesky.h:454
Scalar determinant() const
Definition SimplicialCholesky.h:516
const MatrixL matrixL() const
Definition SimplicialCholesky.h:481
SimplicialLDLT(const MatrixType &matrix)
Definition SimplicialCholesky.h:473
const MatrixU matrixU() const
Definition SimplicialCholesky.h:487
void analyzePattern(const MatrixType &a)
Definition SimplicialCholesky.h:504
void factorize(const MatrixType &a)
Definition SimplicialCholesky.h:513
const VectorType vectorD() const
Definition SimplicialCholesky.h:476
SimplicialLDLT & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:493
SimplicialLDLT()
Definition SimplicialCholesky.h:470
A direct sparse LLT Cholesky factorizations.
Definition SimplicialCholesky.h:371
const MatrixL matrixL() const
Definition SimplicialCholesky.h:392
SimplicialLLT & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:404
Scalar determinant() const
Definition SimplicialCholesky.h:427
const MatrixU matrixU() const
Definition SimplicialCholesky.h:398
SimplicialLLT()
Definition SimplicialCholesky.h:387
void analyzePattern(const MatrixType &a)
Definition SimplicialCholesky.h:415
void factorize(const MatrixType &a)
Definition SimplicialCholesky.h:424
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:626
SimplicialNonHermitianLDLT()
Definition SimplicialCholesky.h:642
Scalar determinant() const
Definition SimplicialCholesky.h:688
SimplicialNonHermitianLDLT & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:665
const MatrixU matrixU() const
Definition SimplicialCholesky.h:659
const MatrixL matrixL() const
Definition SimplicialCholesky.h:653
const VectorType vectorD() const
Definition SimplicialCholesky.h:648
SimplicialNonHermitianLDLT(const MatrixType &matrix)
Definition SimplicialCholesky.h:645
void factorize(const MatrixType &a)
Definition SimplicialCholesky.h:685
void analyzePattern(const MatrixType &a)
Definition SimplicialCholesky.h:676
A direct sparse LLT Cholesky factorizations, for symmetric non-hermitian matrices.
Definition SimplicialCholesky.h:541
void factorize(const MatrixType &a)
Definition SimplicialCholesky.h:595
const MatrixU matrixU() const
Definition SimplicialCholesky.h:569
SimplicialNonHermitianLLT()
Definition SimplicialCholesky.h:557
SimplicialNonHermitianLLT & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:575
Scalar determinant() const
Definition SimplicialCholesky.h:598
const MatrixL matrixL() const
Definition SimplicialCholesky.h:563
void analyzePattern(const MatrixType &a)
Definition SimplicialCholesky.h:586
SimplicialNonHermitianLLT(const MatrixType &matrix)
Definition SimplicialCholesky.h:560
A versatible sparse matrix representation.
Definition SparseMatrix.h:121
Index cols() const
Definition SparseMatrix.h:161
Index rows() const
Definition SparseMatrix.h:159
Index nonZeros() const
Definition SparseCompressedBase.h:64
SparseSolverBase()
Definition SparseSolverBase.h:70
ComputationInfo
Definition Constants.h:438
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Definition SimplicialCholesky.h:232