10#ifndef EIGEN_SIMPLICIAL_CHOLESKY_H
11#define EIGEN_SIMPLICIAL_CHOLESKY_H
14#include "./InternalHeaderCheck.h"
18enum SimplicialCholeskyMode { SimplicialCholeskyLLT, SimplicialCholeskyLDLT };
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) {
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& ) { pmat = &input; }
50template <
typename Derived>
53 using Base::m_isInitialized;
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;
68 enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
75 : m_info(
Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {}
78 : m_info(
Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {
79 derived().compute(matrix);
82 ~SimplicialCholeskyBase() {}
84 Derived& derived() {
return *
static_cast<Derived*
>(
this); }
85 const Derived& derived()
const {
return *
static_cast<const Derived*
>(
this); }
87 inline Index cols()
const {
return m_matrix.cols(); }
88 inline Index rows()
const {
return m_matrix.rows(); }
96 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
118 Derived&
setShift(
const DiagonalScalar& offset,
const DiagonalScalar& scale = 1) {
119 m_shiftOffset = offset;
120 m_shiftScale = scale;
124#ifndef EIGEN_PARSED_BY_DOXYGEN
126 template <
typename Stream>
127 void dumpMemory(Stream& s) {
130 << ((total += (m_matrix.
cols() + 1) *
sizeof(
int) + m_matrix.
nonZeros() * (
sizeof(
int) +
sizeof(
Scalar))) >> 20)
133 s <<
" diag: " << ((total += m_diag.size() *
sizeof(
Scalar)) >> 20) <<
"Mb"
135 s <<
" tree: " << ((total += m_parent.size() *
sizeof(
int)) >> 20) <<
"Mb"
137 s <<
" nonzeros: " << ((total += m_workSpace.size() *
sizeof(
int)) >> 20) <<
"Mb"
139 s <<
" perm: " << ((total += m_P.
size() *
sizeof(
int)) >> 20) <<
"Mb"
141 s <<
" perm^-1: " << ((total += m_Pinv.
size() *
sizeof(
int)) >> 20) <<
"Mb"
143 s <<
" TOTAL: " << (total >> 20) <<
"Mb"
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());
163 derived().matrixL().solveInPlace(dest);
165 if (m_diag.size() > 0) dest = m_diag.asDiagonal().inverse() * dest;
168 derived().matrixU().solveInPlace(dest);
170 if (m_P.
size() > 0) dest = m_Pinv * dest;
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);
182 template <
bool DoLDLT,
bool NonHermitian>
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);
193 template <
bool DoLDLT,
bool NonHermitian>
195 eigen_assert(a.rows() == a.cols());
196 Index size = a.cols();
197 CholMatrixType tmp(size, size);
198 ConstCholMatrixPtr pmat;
202 internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, tmp);
204 internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, tmp, m_P.
indices().data());
208 factorize_preordered<DoLDLT, NonHermitian>(*pmat);
211 template <
bool DoLDLT,
bool NonHermitian>
212 void factorize_preordered(
const CholMatrixType& a);
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);
223 void analyzePattern_preordered(
const CholMatrixType& a,
bool doLDLT);
225 template <
bool NonHermitian>
226 void ordering(
const MatrixType& a, ConstCholMatrixPtr& pmat, CholMatrixType& ap);
228 inline DiagonalScalar getDiag(Scalar x) {
return internal::traits<Derived>::getDiag(x); }
229 inline Scalar getSymm(Scalar x) {
return internal::traits<Derived>::getSymm(x); }
233 inline bool operator()(
const Index& row,
const Index& col,
const Scalar&)
const {
return row != col; }
237 bool m_factorizationIsOk;
240 CholMatrixType m_matrix;
247 DiagonalScalar m_shiftOffset;
248 DiagonalScalar m_shiftScale;
251template <
typename MatrixType_,
int UpLo_ =
Lower,
254template <
typename MatrixType_,
int UpLo_ =
Lower,
257template <
typename MatrixType_,
int UpLo_ =
Lower,
260template <
typename MatrixType_,
int UpLo_ =
Lower,
263template <
typename MatrixType_,
int UpLo_ =
Lower,
269template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
270struct traits<
SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
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); }
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); }
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; }
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; }
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); }
370template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
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;
381 typedef internal::traits<SimplicialLLT> Traits;
382 typedef typename Traits::MatrixL MatrixL;
383 typedef typename Traits::MatrixU MatrixU;
393 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
394 return Traits::getL(Base::m_matrix);
399 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
400 return Traits::getU(Base::m_matrix);
428 Scalar detL = Base::m_matrix.diagonal().prod();
429 return numext::abs2(detL);
453template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
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;
464 typedef internal::traits<SimplicialLDLT> Traits;
465 typedef typename Traits::MatrixL MatrixL;
466 typedef typename Traits::MatrixU MatrixU;
477 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
482 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
483 return Traits::getL(Base::m_matrix);
488 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
489 return Traits::getU(Base::m_matrix);
539template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
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;
551 typedef internal::traits<SimplicialNonHermitianLLT> Traits;
552 typedef typename Traits::MatrixL MatrixL;
553 typedef typename Traits::MatrixU MatrixU;
564 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
565 return Traits::getL(Base::m_matrix);
570 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LLT not factorized");
571 return Traits::getU(Base::m_matrix);
599 Scalar detL = Base::m_matrix.diagonal().prod();
624template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
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;
636 typedef internal::traits<SimplicialNonHermitianLDLT> Traits;
637 typedef typename Traits::MatrixL MatrixL;
638 typedef typename Traits::MatrixU MatrixU;
649 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
654 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
655 return Traits::getL(Base::m_matrix);
660 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial LDLT not factorized");
661 return Traits::getU(Base::m_matrix);
697template <
typename MatrixType_,
int UpLo_,
typename Ordering_>
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;
708 typedef internal::traits<SimplicialLDLT<MatrixType, UpLo> > LDLTTraits;
709 typedef internal::traits<SimplicialLLT<MatrixType, UpLo> > LLTTraits;
712 SimplicialCholesky() : Base(), m_LDLT(
true) {}
714 explicit SimplicialCholesky(
const MatrixType& matrix) : Base(), m_LDLT(
true) {
compute(matrix); }
716 SimplicialCholesky& setMode(SimplicialCholeskyMode mode) {
718 case SimplicialCholeskyLLT:
721 case SimplicialCholeskyLDLT:
731 inline const VectorType vectorD()
const {
732 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial Cholesky not factorized");
735 inline const CholMatrixType rawMatrix()
const {
736 eigen_assert(Base::m_factorizationIsOk &&
"Simplicial Cholesky not factorized");
737 return Base::m_matrix;
741 SimplicialCholesky&
compute(
const MatrixType& matrix) {
777 template <
typename Rhs,
typename Dest>
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());
784 if (Base::m_info !=
Success)
return;
786 if (Base::m_P.size() > 0)
787 dest = Base::m_P * b;
791 if (Base::m_matrix.nonZeros() > 0)
794 LDLTTraits::getL(Base::m_matrix).solveInPlace(dest);
796 LLTTraits::getL(Base::m_matrix).solveInPlace(dest);
799 if (Base::m_diag.size() > 0) dest = Base::m_diag.
real().asDiagonal().inverse() * dest;
801 if (Base::m_matrix.nonZeros() > 0)
804 LDLTTraits::getU(Base::m_matrix).solveInPlace(dest);
806 LLTTraits::getU(Base::m_matrix).solveInPlace(dest);
809 if (Base::m_P.size() > 0) dest = Base::m_Pinv * dest;
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);
818 Scalar determinant()
const {
820 return Base::m_diag.prod();
822 Scalar detL = Diagonal<const CholMatrixType>(Base::m_matrix).prod();
823 return numext::abs2(detL);
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();
841 internal::permute_symm_to_fullsymm<UpLo, NonHermitian>(a, C, NULL);
843 OrderingType ordering;
847 if (m_Pinv.size() > 0)
848 m_P = m_Pinv.inverse();
852 ap.resize(size, size);
853 internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, m_P.indices().data());
857 if (
int(UpLo) ==
int(
Lower) || MatrixType::IsRowMajor) {
859 ap.resize(size, size);
860 internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, NULL);
862 internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, ap);
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
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