10#ifndef EIGEN_SIMPLICIAL_CHOLESKY_H
11#define EIGEN_SIMPLICIAL_CHOLESKY_H
15enum SimplicialCholeskyMode {
16 SimplicialCholeskyLLT,
17 SimplicialCholeskyLDLT
21 template<
typename CholMatrixType,
typename InputMatrixType>
31 template<
typename MatrixType>
54template<
typename Derived>
58 using Base::m_isInitialized;
64 typedef typename MatrixType::Scalar Scalar;
65 typedef typename MatrixType::RealScalar RealScalar;
66 typedef typename MatrixType::StorageIndex StorageIndex;
73 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
74 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
84 m_factorizationIsOk(
false),
85 m_analysisIsOk(
false),
92 m_factorizationIsOk(
false),
93 m_analysisIsOk(
false),
97 derived().compute(matrix);
100 ~SimplicialCholeskyBase()
104 Derived& derived() {
return *
static_cast<Derived*
>(
this); }
105 const Derived& derived()
const {
return *
static_cast<const Derived*
>(
this); }
107 inline Index cols()
const {
return m_matrix.
cols(); }
108 inline Index rows()
const {
return m_matrix.
rows(); }
117 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
140 Derived&
setShift(
const RealScalar& offset,
const RealScalar& scale = 1)
142 m_shiftOffset = offset;
143 m_shiftScale = scale;
147#ifndef EIGEN_PARSED_BY_DOXYGEN
149 template<
typename Stream>
150 void dumpMemory(Stream& s)
153 s <<
" L: " << ((
total+=(m_matrix.
cols()+1) *
sizeof(
int) + m_matrix.
nonZeros()*(
sizeof(int)+
sizeof(Scalar))) >> 20) <<
"Mb" <<
"\n";
154 s <<
" diag: " << ((
total+=m_diag.size() *
sizeof(Scalar)) >> 20) <<
"Mb" <<
"\n";
155 s <<
" tree: " << ((
total+=m_parent.size() *
sizeof(int)) >> 20) <<
"Mb" <<
"\n";
156 s <<
" nonzeros: " << ((
total+=m_nonZerosPerCol.size() *
sizeof(int)) >> 20) <<
"Mb" <<
"\n";
157 s <<
" perm: " << ((
total+=m_P.
size() *
sizeof(int)) >> 20) <<
"Mb" <<
"\n";
158 s <<
" perm^-1: " << ((
total+=m_Pinv.
size() *
sizeof(int)) >> 20) <<
"Mb" <<
"\n";
159 s <<
" TOTAL: " << (
total>> 20) <<
"Mb" <<
"\n";
163 template<
typename Rhs,
typename Dest>
164 void _solve_impl(
const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest)
const
166 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
167 eigen_assert(m_matrix.
rows()==b.rows());
178 derived().matrixL().solveInPlace(dest);
181 dest = m_diag.asDiagonal().inverse() * dest;
184 derived().matrixU().solveInPlace(dest);
187 dest = m_Pinv * dest;
190 template<
typename Rhs,
typename Dest>
191 void _solve_impl(
const SparseMatrixBase<Rhs> &b, SparseMatrixBase<Dest> &dest)
const
193 internal::solve_sparse_through_dense_panels(derived(), b, dest);
201 template<
bool DoLDLT>
204 eigen_assert(matrix.rows()==matrix.cols());
205 Index size = matrix.cols();
207 ConstCholMatrixPtr
pmat;
208 ordering(matrix,
pmat, tmp);
213 template<
bool DoLDLT>
214 void factorize(
const MatrixType& a)
216 eigen_assert(a.rows()==a.cols());
217 Index size = a.cols();
218 CholMatrixType tmp(size,size);
219 ConstCholMatrixPtr
pmat;
228 tmp.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>().twistedBy(m_P);
232 factorize_preordered<DoLDLT>(*pmat);
235 template<
bool DoLDLT>
236 void factorize_preordered(
const CholMatrixType& a);
238 void analyzePattern(
const MatrixType& a,
bool doLDLT)
240 eigen_assert(a.rows()==a.cols());
241 Index size = a.cols();
242 CholMatrixType tmp(size,size);
243 ConstCholMatrixPtr pmat;
244 ordering(a, pmat, tmp);
245 analyzePattern_preordered(*pmat,doLDLT);
247 void analyzePattern_preordered(
const CholMatrixType& a,
bool doLDLT);
249 void ordering(
const MatrixType& a, ConstCholMatrixPtr &pmat, CholMatrixType& ap);
260 bool m_factorizationIsOk;
270 RealScalar m_shiftOffset;
271 RealScalar m_shiftScale;
274template<
typename _MatrixType,
int _UpLo = Lower,
typename _Ordering = AMDOrdering<
typename _MatrixType::StorageIndex> >
class SimplicialLLT;
275template<
typename _MatrixType,
int _UpLo = Lower,
typename _Ordering = AMDOrdering<
typename _MatrixType::StorageIndex> >
class SimplicialLDLT;
276template<
typename _MatrixType,
int _UpLo = Lower,
typename _Ordering = AMDOrdering<
typename _MatrixType::StorageIndex> >
class SimplicialCholesky;
282 typedef _MatrixType MatrixType;
284 enum { UpLo =
_UpLo };
285 typedef typename MatrixType::Scalar Scalar;
286 typedef typename MatrixType::StorageIndex StorageIndex;
296 typedef _MatrixType MatrixType;
298 enum { UpLo =
_UpLo };
299 typedef typename MatrixType::Scalar Scalar;
300 typedef typename MatrixType::StorageIndex StorageIndex;
310 typedef _MatrixType MatrixType;
312 enum { UpLo =
_UpLo };
337template<
typename _MatrixType,
int _UpLo,
typename _Ordering>
341 typedef _MatrixType MatrixType;
342 enum { UpLo =
_UpLo };
344 typedef typename MatrixType::Scalar Scalar;
345 typedef typename MatrixType::RealScalar RealScalar;
346 typedef typename MatrixType::StorageIndex StorageIndex;
350 typedef typename Traits::MatrixL MatrixL;
351 typedef typename Traits::MatrixU MatrixU;
361 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
362 return Traits::getL(Base::m_matrix);
367 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
368 return Traits::getU(Base::m_matrix);
386 Base::analyzePattern(a,
false);
404 return numext::abs2(
detL);
428template<
typename _MatrixType,
int _UpLo,
typename _Ordering>
432 typedef _MatrixType MatrixType;
433 enum { UpLo =
_UpLo };
435 typedef typename MatrixType::Scalar Scalar;
436 typedef typename MatrixType::RealScalar RealScalar;
437 typedef typename MatrixType::StorageIndex StorageIndex;
441 typedef typename Traits::MatrixL MatrixL;
442 typedef typename Traits::MatrixU MatrixU;
453 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
458 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
459 return Traits::getL(Base::m_matrix);
464 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
465 return Traits::getU(Base::m_matrix);
483 Base::analyzePattern(a,
true);
500 return Base::m_diag.prod();
510template<
typename _MatrixType,
int _UpLo,
typename _Ordering>
514 typedef _MatrixType MatrixType;
515 enum { UpLo =
_UpLo };
517 typedef typename MatrixType::Scalar Scalar;
518 typedef typename MatrixType::RealScalar RealScalar;
519 typedef typename MatrixType::StorageIndex StorageIndex;
529 :
Base(), m_LDLT(
true)
538 case SimplicialCholeskyLLT:
541 case SimplicialCholeskyLDLT:
552 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial Cholesky not factorized");
556 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial Cholesky not factorized");
557 return Base::m_matrix;
578 Base::analyzePattern(a, m_LDLT);
596 template<
typename Rhs,
typename Dest>
599 eigen_assert(Base::m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
600 eigen_assert(Base::m_matrix.rows()==b.rows());
605 if(Base::m_P.size()>0)
606 dest = Base::m_P * b;
610 if(Base::m_matrix.nonZeros()>0)
613 LDLTTraits::getL(Base::m_matrix).solveInPlace(dest);
615 LLTTraits::getL(Base::m_matrix).solveInPlace(dest);
618 if(Base::m_diag.size()>0)
619 dest = Base::m_diag.real().asDiagonal().inverse() * dest;
621 if (Base::m_matrix.nonZeros()>0)
624 LDLTTraits::getU(Base::m_matrix).solveInPlace(dest);
626 LLTTraits::getU(Base::m_matrix).solveInPlace(dest);
629 if(Base::m_P.size()>0)
630 dest = Base::m_Pinv * dest;
634 template<
typename Rhs,
typename Dest>
635 void _solve_impl(
const SparseMatrixBase<Rhs> &b, SparseMatrixBase<Dest> &dest)
const
637 internal::solve_sparse_through_dense_panels(*
this, b, dest);
640 Scalar determinant()
const
644 return Base::m_diag.prod();
648 Scalar detL = Diagonal<const CholMatrixType>(Base::m_matrix).prod();
649 return numext::abs2(detL);
657template<
typename Derived>
658void SimplicialCholeskyBase<Derived>::ordering(
const MatrixType& a, ConstCholMatrixPtr &pmat, CholMatrixType& ap)
660 eigen_assert(a.rows()==a.cols());
661 const Index size = a.rows();
664 if(!internal::is_same<OrderingType,NaturalOrdering<Index> >
::value)
668 C = a.template selfadjointView<UpLo>();
670 OrderingType ordering;
674 if(m_Pinv.size()>0) m_P = m_Pinv.inverse();
677 ap.resize(size,size);
678 ap.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>().twistedBy(m_P);
684 if(
int(UpLo)==
int(
Lower) || MatrixType::IsRowMajor)
687 ap.resize(size,size);
688 ap.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>();
691 internal::simplicial_cholesky_grab_input<CholMatrixType,MatrixType>::run(a, pmat, ap);
EIGEN_DEVICE_FUNC CoeffReturnType value() const
Definition DenseBase.h:526
EIGEN_DEVICE_FUNC Scalar prod() const
Definition Redux.h:493
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
EIGEN_DEVICE_FUNC DiagonalReturnType diagonal()
Definition Diagonal.h:187
EIGEN_DEVICE_FUNC Index size() const
Definition PermutationMatrix.h:97
A base class for direct sparse Cholesky factorizations.
Definition SimplicialCholesky.h:56
SimplicialCholeskyBase()
Default constructor.
Definition SimplicialCholesky.h:82
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SimplicialCholesky.h:115
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP() const
Definition SimplicialCholesky.h:123
Derived & setShift(const RealScalar &offset, const RealScalar &scale=1)
Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical ...
Definition SimplicialCholesky.h:140
void compute(const MatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix.
Definition SimplicialCholesky.h:202
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv() const
Definition SimplicialCholesky.h:128
Definition SimplicialCholesky.h:512
SimplicialCholesky & compute(const MatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix.
Definition SimplicialCholesky.h:561
void analyzePattern(const MatrixType &a)
Performs a symbolic decomposition on the sparcity of matrix.
Definition SimplicialCholesky.h:576
void factorize(const MatrixType &a)
Performs a numeric decomposition of matrix.
Definition SimplicialCholesky.h:587
A direct sparse LDLT Cholesky factorizations without square root.
Definition SimplicialCholesky.h:430
SimplicialLDLT(const MatrixType &matrix)
Constructs and performs the LLT factorization of matrix.
Definition SimplicialCholesky.h:448
SimplicialLDLT()
Default constructor.
Definition SimplicialCholesky.h:445
SimplicialLDLT & compute(const MatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix.
Definition SimplicialCholesky.h:469
void factorize(const MatrixType &a)
Performs a numeric decomposition of matrix.
Definition SimplicialCholesky.h:492
Scalar determinant() const
Definition SimplicialCholesky.h:498
void analyzePattern(const MatrixType &a)
Performs a symbolic decomposition on the sparcity of matrix.
Definition SimplicialCholesky.h:481
const VectorType vectorD() const
Definition SimplicialCholesky.h:452
const MatrixL matrixL() const
Definition SimplicialCholesky.h:457
const MatrixU matrixU() const
Definition SimplicialCholesky.h:463
A direct sparse LLT Cholesky factorizations.
Definition SimplicialCholesky.h:339
const MatrixU matrixU() const
Definition SimplicialCholesky.h:366
SimplicialLLT(const MatrixType &matrix)
Constructs and performs the LLT factorization of matrix.
Definition SimplicialCholesky.h:356
SimplicialLLT & compute(const MatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix.
Definition SimplicialCholesky.h:372
void factorize(const MatrixType &a)
Performs a numeric decomposition of matrix.
Definition SimplicialCholesky.h:395
Scalar determinant() const
Definition SimplicialCholesky.h:401
SimplicialLLT()
Default constructor.
Definition SimplicialCholesky.h:354
void analyzePattern(const MatrixType &a)
Performs a symbolic decomposition on the sparcity of matrix.
Definition SimplicialCholesky.h:384
const MatrixL matrixL() const
Definition SimplicialCholesky.h:360
Index nonZeros() const
Definition SparseCompressedBase.h:56
Index rows() const
Definition SparseMatrix.h:138
Index cols() const
Definition SparseMatrix.h:140
A base class for sparse solvers.
Definition SparseSolverBase.h:68
ComputationInfo
Enum for reporting the status of a computation.
Definition Constants.h:440
@ Lower
View matrix as a lower triangular matrix.
Definition Constants.h:209
@ Upper
View matrix as an upper triangular matrix.
Definition Constants.h:211
@ Success
Computation was successful.
Definition Constants.h:442
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
keeps off-diagonal entries; drops diagonal entries
Definition SimplicialCholesky.h:252
Definition ForwardDeclarations.h:17