Eigen::CholmodBase< MatrixType_, UpLo_, Derived > Class Template Reference

The base class for the direct Cholesky factorization of Cholmod. More...

#include <CholmodSupport.h>

+ Inheritance diagram for Eigen::CholmodBase< MatrixType_, UpLo_, Derived >:

Public Types

enum  { UpLo = UpLo_ }
 
enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType CholMatrixType
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

 CholmodBase ()
 
 CholmodBase (const MatrixType &matrix)
 
 ~CholmodBase ()
 
StorageIndex cols () const
 
StorageIndex rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
Derived & compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
cholmod_common & cholmod ()
 
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
template<typename RhsDerived , typename DestDerived >
void _solve_impl (const SparseMatrixBase< RhsDerived > &b, SparseMatrixBase< DestDerived > &dest) const
 
Derived & setShift (const RealScalar &offset)
 
Scalar determinant () const
 
Scalar logDeterminant () const
 
template<typename Stream >
void dumpMemory (Stream &)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< Derived > Base
 

Protected Member Functions

Derived & derived ()
 
const Derived & derived () const
 

Protected Attributes

cholmod_common m_cholmod
 
cholmod_factor * m_cholmodFactor
 
double m_shiftOffset [2]
 
ComputationInfo m_info
 
int m_factorizationIsOk
 
int m_analysisIsOk
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Detailed Description

template<typename MatrixType_, int UpLo_, typename Derived>
class Eigen::CholmodBase< MatrixType_, UpLo_, Derived >

The base class for the direct Cholesky factorization of Cholmod.

See also
class CholmodSupernodalLLT, class CholmodSimplicialLDLT, class CholmodSimplicialLLT

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef SparseSolverBase<Derived> Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::Base
protected

◆ CholMatrixType

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholMatrixType

◆ MatrixType

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType_ Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::MatrixType

◆ RealScalar

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::RealScalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::RealScalar

◆ Scalar

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::Scalar

◆ StorageIndex

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::StorageIndex Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::StorageIndex

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Derived >
anonymous enum
Enumerator
UpLo 
253 { UpLo = UpLo_ };
@ UpLo
Definition: CholmodSupport.h:253

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Derived >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
258 { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
@ ColsAtCompileTime
Definition: CholmodSupport.h:258
@ MaxColsAtCompileTime
Definition: CholmodSupport.h:258

Constructor & Destructor Documentation

◆ CholmodBase() [1/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholmodBase ( )
inline
262  EIGEN_STATIC_ASSERT((internal::is_same<double, RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
263  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
264  internal::cm_start<StorageIndex>(m_cholmod);
265  }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
int m_factorizationIsOk
Definition: CholmodSupport.h:464
ComputationInfo m_info
Definition: CholmodSupport.h:463
cholmod_factor * m_cholmodFactor
Definition: CholmodSupport.h:461
double m_shiftOffset[2]
Definition: CholmodSupport.h:462
cholmod_common m_cholmod
Definition: CholmodSupport.h:460
int m_analysisIsOk
Definition: CholmodSupport.h:465
@ Success
Definition: Constants.h:440
@ value
Definition: Meta.h:206

References EIGEN_STATIC_ASSERT, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, and Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset.

◆ CholmodBase() [2/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholmodBase ( const MatrixType matrix)
inlineexplicit
269  EIGEN_STATIC_ASSERT((internal::is_same<double, RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
270  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
271  internal::cm_start<StorageIndex>(m_cholmod);
272  compute(matrix);
273  }
Derived & compute(const MatrixType &matrix)
Definition: CholmodSupport.h:294
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

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::compute(), EIGEN_STATIC_ASSERT, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset, and matrix().

◆ ~CholmodBase()

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::~CholmodBase ( )
inline
275  {
276  if (m_cholmodFactor) internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
277  internal::cm_finish<StorageIndex>(m_cholmod);
278  }

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, and Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor.

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
template<typename Rhs , typename Dest >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
inline
346  {
348  "The decomposition is not in a valid state for solving, you must first call either compute() or "
349  "symbolic()/numeric()");
350  const Index size = m_cholmodFactor->n;
352  eigen_assert(size == b.rows());
353 
354  // Cholmod needs column-major storage without inner-stride, which corresponds to the default behavior of Ref.
355  Ref<const Matrix<typename Rhs::Scalar, Dynamic, Dynamic, ColMajor> > b_ref(b.derived());
356 
357  cholmod_dense b_cd = viewAsCholmod(b_ref);
358  cholmod_dense* x_cd = internal::cm_solve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cd, m_cholmod);
359  if (!x_cd) {
360  this->m_info = NumericalIssue;
361  return;
362  }
363  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
364  // NOTE Actually, the copy can be avoided by calling cholmod_solve2 instead of cholmod_solve
366  b.rows(), b.cols());
367  internal::cm_free_dense<StorageIndex>(x_cd, m_cholmod);
368  }
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:966
#define eigen_assert(x)
Definition: Macros.h:910
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:595
@ NumericalIssue
Definition: Constants.h:442
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)
Definition: CholmodSupport.h:64

References b, eigen_assert, EIGEN_UNUSED_VARIABLE, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info, Eigen::PlainObjectBase< Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ > >::Map(), Eigen::NumericalIssue, size, and Eigen::viewAsCholmod().

◆ _solve_impl() [2/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
template<typename RhsDerived , typename DestDerived >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::_solve_impl ( const SparseMatrixBase< RhsDerived > &  b,
SparseMatrixBase< DestDerived > &  dest 
) const
inline
372  {
374  "The decomposition is not in a valid state for solving, you must first call either compute() or "
375  "symbolic()/numeric()");
376  const Index size = m_cholmodFactor->n;
378  eigen_assert(size == b.rows());
379 
380  // note: cs stands for Cholmod Sparse
381  Ref<SparseMatrix<typename RhsDerived::Scalar, ColMajor, typename RhsDerived::StorageIndex> > b_ref(
382  b.const_cast_derived());
383  cholmod_sparse b_cs = viewAsCholmod(b_ref);
384  cholmod_sparse* x_cs = internal::cm_spsolve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cs, m_cholmod);
385  if (!x_cs) {
386  this->m_info = NumericalIssue;
387  return;
388  }
389  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
390  // NOTE cholmod_spsolve in fact just calls the dense solver for blocks of 4 columns at a time (similar to Eigen's
391  // sparse solver)
392  dest.derived() = viewAsEigen<typename DestDerived::Scalar, typename DestDerived::StorageIndex>(*x_cs);
393  internal::cm_free_sparse<StorageIndex>(x_cs, m_cholmod);
394  }

References b, Eigen::SparseMatrixBase< Derived >::derived(), eigen_assert, EIGEN_UNUSED_VARIABLE, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info, Eigen::NumericalIssue, size, and Eigen::viewAsCholmod().

◆ analyzePattern()

template<typename MatrixType_ , int UpLo_, typename Derived >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::analyzePattern ( const MatrixType matrix)
inline

Performs a symbolic decomposition on the sparsity pattern of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
306  {
307  if (m_cholmodFactor) {
308  internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
309  m_cholmodFactor = 0;
310  }
311  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
312  m_cholmodFactor = internal::cm_analyze<StorageIndex>(A, m_cholmod);
313 
314  this->m_isInitialized = true;
315  this->m_info = Success;
316  m_analysisIsOk = true;
317  m_factorizationIsOk = false;
318  }
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
bool m_isInitialized
Definition: SparseSolverBase.h:110

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_analysisIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_isInitialized, matrix(), Eigen::Success, and Eigen::viewAsCholmod().

Referenced by Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::compute().

◆ cholmod()

template<typename MatrixType_ , int UpLo_, typename Derived >
cholmod_common& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::cholmod ( )
inline

Returns a reference to the Cholmod's configuration structure to get a full control over the performed operations. See the Cholmod user guide for details.

341 { return m_cholmod; }

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod.

◆ cols()

◆ compute()

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::compute ( const MatrixType matrix)
inline

Computes the sparse Cholesky decomposition of matrix

294  {
296  factorize(matrix);
297  return derived();
298  }
void factorize(const MatrixType &matrix)
Definition: CholmodSupport.h:327
Derived & derived()
Definition: SparseSolverBase.h:76
void analyzePattern(const MatrixType &matrix)
Definition: CholmodSupport.h:306

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::analyzePattern(), Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::derived(), Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::factorize(), and matrix().

Referenced by Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholmodBase().

◆ derived() [1/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

◆ derived() [2/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
const Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected
77 { return *static_cast<const Derived*>(this); }

◆ determinant()

template<typename MatrixType_ , int UpLo_, typename Derived >
Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization
412  {
413  using std::exp;
414  return exp(logDeterminant());
415  }
Scalar logDeterminant() const
Definition: CholmodSupport.h:418
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 exp(const bfloat16 &a)
Definition: BFloat16.h:615

References Eigen::bfloat16_impl::exp(), and Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::logDeterminant().

◆ dumpMemory()

template<typename MatrixType_ , int UpLo_, typename Derived >
template<typename Stream >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::dumpMemory ( Stream )
inline
457 {}

◆ factorize()

template<typename MatrixType_ , int UpLo_, typename Derived >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::factorize ( const MatrixType matrix)
inline

Performs a numeric decomposition of matrix

The given matrix must have the same sparsity pattern as the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
327  {
328  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
329  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
330  internal::cm_factorize_p<StorageIndex>(&A, m_shiftOffset, 0, 0, m_cholmodFactor, m_cholmod);
331 
332  // If the factorization failed, either the input matrix was zero (so m_cholmodFactor == nullptr), or minor is the
333  // column at which it failed. On success minor == n.
334  this->m_info =
335  (m_cholmodFactor != nullptr && m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
336  m_factorizationIsOk = true;
337  }

References eigen_assert, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_analysisIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset, matrix(), Eigen::NumericalIssue, Eigen::Success, and Eigen::viewAsCholmod().

Referenced by Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::compute().

◆ info()

template<typename MatrixType_ , int UpLo_, typename Derived >
ComputationInfo Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the matrix.appears to be negative.
288  {
289  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
290  return m_info;
291  }

References eigen_assert, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info, and Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_isInitialized.

◆ logDeterminant()

template<typename MatrixType_ , int UpLo_, typename Derived >
Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::logDeterminant ( ) const
inline
Returns
the log determinant of the underlying matrix from the current factorization
418  {
419  using numext::real;
420  using std::log;
422  "The decomposition is not in a valid state for solving, you must first call either compute() or "
423  "symbolic()/numeric()");
424 
425  RealScalar logDet = 0;
426  Scalar* x = static_cast<Scalar*>(m_cholmodFactor->x);
427  if (m_cholmodFactor->is_super) {
428  // Supernodal factorization stored as a packed list of dense column-major blocks,
429  // as described by the following structure:
430 
431  // super[k] == index of the first column of the j-th super node
432  StorageIndex* super = static_cast<StorageIndex*>(m_cholmodFactor->super);
433  // pi[k] == offset to the description of row indices
434  StorageIndex* pi = static_cast<StorageIndex*>(m_cholmodFactor->pi);
435  // px[k] == offset to the respective dense block
436  StorageIndex* px = static_cast<StorageIndex*>(m_cholmodFactor->px);
437 
438  Index nb_super_nodes = m_cholmodFactor->nsuper;
439  for (Index k = 0; k < nb_super_nodes; ++k) {
440  StorageIndex ncols = super[k + 1] - super[k];
441  StorageIndex nrows = pi[k + 1] - pi[k];
442 
443  Map<const Array<Scalar, 1, Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows + 1));
444  logDet += sk.real().log().sum();
445  }
446  } else {
447  // Simplicial factorization stored as standard CSC matrix.
448  StorageIndex* p = static_cast<StorageIndex*>(m_cholmodFactor->p);
450  for (Index k = 0; k < size; ++k) logDet += log(real(x[p[k]]));
451  }
452  if (m_cholmodFactor->is_ll) logDet *= 2.0;
453  return logDet;
454  }
float * p
Definition: Tutorial_Map_using.cpp:9
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixType::StorageIndex StorageIndex
Definition: CholmodSupport.h:257
float real
Definition: datatypes.h:10
RealScalar RealScalar * px
Definition: level1_cplx_impl.h:27
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 log(const bfloat16 &a)
Definition: BFloat16.h:618
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:486
const Mdouble pi
Definition: ExtendedMath.h:23
list x
Definition: plotDoE.py:28

References eigen_assert, k, Eigen::bfloat16_impl::log(), Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor, Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk, p, constants::pi, px, Eigen::real(), size, and plotDoE::x.

Referenced by Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::determinant().

◆ rows()

◆ setShift()

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::setShift ( const RealScalar offset)
inline

Sets the shift parameter that will be used to adjust the diagonal coefficients during the numerical factorization.

During the numerical factorization, an offset term is added to the diagonal coefficients:
d_ii = offset + d_ii

The default is offset=0.

Returns
a reference to *this.
406  {
407  m_shiftOffset[0] = double(offset);
408  return derived();
409  }

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::derived(), and Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset.

Member Data Documentation

◆ m_analysisIsOk

template<typename MatrixType_ , int UpLo_, typename Derived >
int Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_analysisIsOk
protected

◆ m_cholmod

◆ m_cholmodFactor

◆ m_factorizationIsOk

◆ m_info

◆ m_isInitialized

template<typename MatrixType_ , int UpLo_, typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_shiftOffset

template<typename MatrixType_ , int UpLo_, typename Derived >
double Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset[2]
protected

The documentation for this class was generated from the following file: