Eigen::UmfPackLU< MatrixType_ > Class Template Reference

A sparse LU factorization and solver based on UmfPack. More...

+ Inheritance diagram for Eigen::UmfPackLU< MatrixType_ >:

Public Types

enum  {
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
 
typedef Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
 
typedef SparseMatrix< ScalarLUMatrixType
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Array< double, UMFPACK_CONTROL, 1 > UmfpackControl
 
typedef Array< double, UMFPACK_INFO, 1 > UmfpackInfo
 
typedef Ref< const UmfpackMatrixType, StandardCompressedFormatUmfpackMatrixRef
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexUmfpackMatrixType
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 

Public Member Functions

template<typename BDerived , typename XDerived >
bool _solve_impl (const MatrixBase< BDerived > &b, MatrixBase< XDerived > &x) const
 
template<typename InputMatrixType >
void analyzePattern (const InputMatrixType &matrix)
 
Index cols () const
 
template<typename InputMatrixType >
void compute (const InputMatrixType &matrix)
 
Scalar determinant () const
 
void extractData () const
 
template<typename InputMatrixType >
void factorize (const InputMatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const LUMatrixTypematrixL () const
 
const LUMatrixTypematrixU () const
 
const IntColVectorTypepermutationP () const
 
const IntRowVectorTypepermutationQ () const
 
void printUmfpackControl ()
 
void printUmfpackInfo ()
 
void printUmfpackStatus ()
 
Index rows () const
 
UmfpackControlumfpackControl ()
 
const UmfpackControlumfpackControl () const
 
int umfpackFactorizeReturncode () const
 
 UmfPackLU ()
 
template<typename InputMatrixType >
 UmfPackLU (const InputMatrixType &matrix)
 
 ~UmfPackLU ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< UmfPackLU< MatrixType_ > >
UmfPackLU< MatrixType_ > & derived ()
 
const UmfPackLU< MatrixType_ > & derived () const
 
const Solve< UmfPackLU< MatrixType_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< UmfPackLU< MatrixType_ >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Protected Types

typedef SparseSolverBase< UmfPackLU< MatrixType_ > > Base
 

Protected Member Functions

void analyzePattern_impl ()
 
void factorize_impl ()
 
template<typename MatrixDerived >
void grab (const EigenBase< MatrixDerived > &A)
 
void grab (const UmfpackMatrixRef &A)
 
void init ()
 

Protected Attributes

int m_analysisIsOk
 
UmfpackControl m_control
 
UmfpackMatrixType m_dummy
 
bool m_extractedDataAreDirty
 
StorageIndex m_fact_errorCode
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
LUMatrixType m_l
 
void * m_numeric
 
IntColVectorType m_p
 
IntRowVectorType m_q
 
void * m_symbolic
 
LUMatrixType m_u
 
UmfpackInfo m_umfpackInfo
 
UmfpackMatrixRef mp_matrix
 
- Protected Attributes inherited from Eigen::SparseSolverBase< UmfPackLU< MatrixType_ > >
bool m_isInitialized
 

Private Member Functions

 UmfPackLU (const UmfPackLU &)
 

Detailed Description

template<typename MatrixType_>
class Eigen::UmfPackLU< MatrixType_ >

A sparse LU factorization and solver based on UmfPack.

This class allows to solve for A.X = B sparse linear problems via a LU factorization using the UmfPack library. The sparse matrix A must be squared and full rank. The vectors or matrices X and B can be either dense or sparse.

Warning
The input matrix A should be in a compressed and column-major form. Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
Template Parameters
MatrixType_the type of the sparse matrix A, it must be a SparseMatrix<>

This class follows the sparse solver concept .

See also
Sparse solver concept, class SparseLU

Definition at line 290 of file UmfPackSupport.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ >
typedef SparseSolverBase<UmfPackLU<MatrixType_> > Eigen::UmfPackLU< MatrixType_ >::Base
protected

Definition at line 293 of file UmfPackSupport.h.

◆ IntColVectorType

template<typename MatrixType_ >
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> Eigen::UmfPackLU< MatrixType_ >::IntColVectorType

Definition at line 303 of file UmfPackSupport.h.

◆ IntRowVectorType

template<typename MatrixType_ >
typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> Eigen::UmfPackLU< MatrixType_ >::IntRowVectorType

Definition at line 302 of file UmfPackSupport.h.

◆ LUMatrixType

template<typename MatrixType_ >
typedef SparseMatrix<Scalar> Eigen::UmfPackLU< MatrixType_ >::LUMatrixType

Definition at line 304 of file UmfPackSupport.h.

◆ MatrixType

template<typename MatrixType_ >
typedef MatrixType_ Eigen::UmfPackLU< MatrixType_ >::MatrixType

Definition at line 297 of file UmfPackSupport.h.

◆ RealScalar

template<typename MatrixType_ >
typedef MatrixType::RealScalar Eigen::UmfPackLU< MatrixType_ >::RealScalar

Definition at line 299 of file UmfPackSupport.h.

◆ Scalar

template<typename MatrixType_ >
typedef MatrixType::Scalar Eigen::UmfPackLU< MatrixType_ >::Scalar

Definition at line 298 of file UmfPackSupport.h.

◆ StorageIndex

template<typename MatrixType_ >
typedef MatrixType::StorageIndex Eigen::UmfPackLU< MatrixType_ >::StorageIndex

Definition at line 300 of file UmfPackSupport.h.

◆ UmfpackControl

template<typename MatrixType_ >
typedef Array<double, UMFPACK_CONTROL, 1> Eigen::UmfPackLU< MatrixType_ >::UmfpackControl

Definition at line 314 of file UmfPackSupport.h.

◆ UmfpackInfo

template<typename MatrixType_ >
typedef Array<double, UMFPACK_INFO, 1> Eigen::UmfPackLU< MatrixType_ >::UmfpackInfo

Definition at line 315 of file UmfPackSupport.h.

◆ UmfpackMatrixRef

template<typename MatrixType_ >
typedef Ref<const UmfpackMatrixType, StandardCompressedFormat> Eigen::UmfPackLU< MatrixType_ >::UmfpackMatrixRef

Definition at line 306 of file UmfPackSupport.h.

◆ UmfpackMatrixType

template<typename MatrixType_ >
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::UmfPackLU< MatrixType_ >::UmfpackMatrixType

Definition at line 305 of file UmfPackSupport.h.

◆ Vector

template<typename MatrixType_ >
typedef Matrix<Scalar,Dynamic,1> Eigen::UmfPackLU< MatrixType_ >::Vector

Definition at line 301 of file UmfPackSupport.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 307 of file UmfPackSupport.h.

307  {
308  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
309  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
310  };

Constructor & Destructor Documentation

◆ UmfPackLU() [1/3]

template<typename MatrixType_ >
Eigen::UmfPackLU< MatrixType_ >::UmfPackLU ( )
inline

Definition at line 317 of file UmfPackSupport.h.

318  : m_dummy(0,0), mp_matrix(m_dummy)
319  {
320  init();
321  }
UmfpackMatrixRef mp_matrix
UmfpackMatrixType m_dummy

◆ UmfPackLU() [2/3]

template<typename MatrixType_ >
template<typename InputMatrixType >
Eigen::UmfPackLU< MatrixType_ >::UmfPackLU ( const InputMatrixType &  matrix)
inlineexplicit

Definition at line 324 of file UmfPackSupport.h.

325  : mp_matrix(matrix)
326  {
327  init();
328  compute(matrix);
329  }
void compute(const InputMatrixType &matrix)

◆ ~UmfPackLU()

template<typename MatrixType_ >
Eigen::UmfPackLU< MatrixType_ >::~UmfPackLU ( )
inline

Definition at line 331 of file UmfPackSupport.h.

332  {
335  }
MatrixType::Scalar Scalar
MatrixType::StorageIndex StorageIndex
void umfpack_free_numeric(void **Numeric, double, int)
void umfpack_free_symbolic(void **Symbolic, double, int)

◆ UmfPackLU() [3/3]

template<typename MatrixType_ >
Eigen::UmfPackLU< MatrixType_ >::UmfPackLU ( const UmfPackLU< MatrixType_ > &  )
inlineprivate

Definition at line 569 of file UmfPackSupport.h.

569 { }

Member Function Documentation

◆ _solve_impl()

template<typename MatrixType >
template<typename BDerived , typename XDerived >
bool Eigen::UmfPackLU< MatrixType >::_solve_impl ( const MatrixBase< BDerived > &  b,
MatrixBase< XDerived > &  x 
) const

Definition at line 611 of file UmfPackSupport.h.

612 {
613  Index rhsCols = b.cols();
614  eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
615  eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
616  eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
617 
618  Scalar* x_ptr = 0;
619  Matrix<Scalar,Dynamic,1> x_tmp;
620  if(x.innerStride()!=1)
621  {
622  x_tmp.resize(x.rows());
623  x_ptr = x_tmp.data();
624  }
625  for (int j=0; j<rhsCols; ++j)
626  {
627  if(x.innerStride()==1)
628  x_ptr = &x.col(j).coeffRef(0);
629  StorageIndex errorCode = umfpack_solve(UMFPACK_A,
630  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
631  x_ptr, &b.const_cast_derived().col(j).coeffRef(0),
633  if(x.innerStride()!=1)
634  x.col(j) = x_tmp;
635  if (errorCode!=0)
636  return false;
637  }
638 
639  return true;
640 }
Array< int, 3, 1 > b
#define eigen_assert(x)
Definition: Macros.h:902
const Scalar * data() const
UmfpackInfo m_umfpackInfo
UmfpackControl m_control
const unsigned int RowMajorBit
Definition: Constants.h:68
int umfpack_solve(int sys, const int Ap[], const int Ai[], const double Ax[], double X[], const double B[], void *Numeric, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
std::ptrdiff_t j

◆ analyzePattern()

template<typename MatrixType_ >
template<typename InputMatrixType >
void Eigen::UmfPackLU< MatrixType_ >::analyzePattern ( const InputMatrixType &  matrix)
inline

Performs a symbolic decomposition on the sparcity of matrix.

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

See also
factorize(), compute()

Definition at line 396 of file UmfPackSupport.h.

397  {
400 
401  grab(matrix.derived());
402 
404  }
void analyzePattern_impl()
void grab(const EigenBase< MatrixDerived > &A)

◆ analyzePattern_impl()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::analyzePattern_impl ( )
inlineprotected

Definition at line 506 of file UmfPackSupport.h.

507  {
508  m_fact_errorCode = umfpack_symbolic(internal::convert_index<StorageIndex>(mp_matrix.rows()),
509  internal::convert_index<StorageIndex>(mp_matrix.cols()),
510  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
512 
513  m_isInitialized = true;
515  m_analysisIsOk = true;
516  m_factorizationIsOk = false;
518  }
ComputationInfo m_info
StorageIndex m_fact_errorCode
@ InvalidInput
Definition: Constants.h:453
@ Success
Definition: Constants.h:446
int umfpack_symbolic(int n_row, int n_col, const int Ap[], const int Ai[], const double Ax[], void **Symbolic, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])

◆ cols()

template<typename MatrixType_ >
Index Eigen::UmfPackLU< MatrixType_ >::cols ( void  ) const
inline

Definition at line 338 of file UmfPackSupport.h.

338 { return mp_matrix.cols(); }

◆ compute()

template<typename MatrixType_ >
template<typename InputMatrixType >
void Eigen::UmfPackLU< MatrixType_ >::compute ( const InputMatrixType &  matrix)
inline

Computes the sparse Cholesky decomposition of matrix Note that the matrix should be column-major, and in compressed format for best performance.

See also
SparseMatrix::makeCompressed().

Definition at line 380 of file UmfPackSupport.h.

381  {
384  grab(matrix.derived());
386  factorize_impl();
387  }

◆ determinant()

template<typename MatrixType >
UmfPackLU< MatrixType >::Scalar Eigen::UmfPackLU< MatrixType >::determinant

Definition at line 602 of file UmfPackSupport.h.

603 {
604  Scalar det;
606  return det;
607 }
int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info[UMFPACK_INFO], int)

◆ extractData()

template<typename MatrixType >
void Eigen::UmfPackLU< MatrixType >::extractData

Definition at line 574 of file UmfPackSupport.h.

575 {
577  {
578  // get size of the data
579  StorageIndex lnz, unz, rows, cols, nz_udiag;
580  umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
581 
582  // allocate data
584  m_l.resizeNonZeros(lnz);
585 
587  m_u.resizeNonZeros(unz);
588 
589  m_p.resize(rows);
590  m_q.resize(cols);
591 
592  // extract
595  m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
596 
597  m_extractedDataAreDirty = false;
598  }
599 }
constexpr void resize(Index rows, Index cols)
const Scalar * valuePtr() const
Definition: SparseMatrix.h:177
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:739
void resizeNonZeros(Index size)
Definition: SparseMatrix.h:758
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:195
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:186
Index cols() const
Index rows() const
IntColVectorType m_p
LUMatrixType m_l
IntRowVectorType m_q
LUMatrixType m_u
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[], int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)

◆ factorize()

template<typename MatrixType_ >
template<typename InputMatrixType >
void Eigen::UmfPackLU< MatrixType_ >::factorize ( const InputMatrixType &  matrix)
inline

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.

See also
analyzePattern(), compute()

Definition at line 446 of file UmfPackSupport.h.

447  {
448  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
449  if(m_numeric)
451 
452  grab(matrix.derived());
453 
454  factorize_impl();
455  }

◆ factorize_impl()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::factorize_impl ( )
inlineprotected

Definition at line 520 of file UmfPackSupport.h.

521  {
522 
523  m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
525 
526  m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
527  m_factorizationIsOk = true;
529  }
@ NumericalIssue
Definition: Constants.h:448
int umfpack_numeric(const int Ap[], const int Ai[], const double Ax[], void *Symbolic, void **Numeric, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])

◆ grab() [1/2]

template<typename MatrixType_ >
template<typename MatrixDerived >
void Eigen::UmfPackLU< MatrixType_ >::grab ( const EigenBase< MatrixDerived > &  A)
inlineprotected

Definition at line 532 of file UmfPackSupport.h.

533  {
535  internal::construct_at(&mp_matrix, A.derived());
536  }
MatrixXcf A
void destroy_at(T *p)
Definition: Memory.h:1264
T * construct_at(T *p, Args &&... args)
Definition: Memory.h:1248

◆ grab() [2/2]

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::grab ( const UmfpackMatrixRef A)
inlineprotected

Definition at line 538 of file UmfPackSupport.h.

539  {
540  if(&(A.derived()) != &mp_matrix)
541  {
544  }
545  }

◆ info()

template<typename MatrixType_ >
ComputationInfo Eigen::UmfPackLU< MatrixType_ >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the matrix.appears to be negative.

Definition at line 345 of file UmfPackSupport.h.

346  {
347  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
348  return m_info;
349  }

◆ init()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::init ( )
inlineprotected

Definition at line 495 of file UmfPackSupport.h.

496  {
498  m_isInitialized = false;
499  m_numeric = 0;
500  m_symbolic = 0;
502 
504  }
void umfpack_defaults(double control[UMFPACK_CONTROL], double, int)

◆ matrixL()

template<typename MatrixType_ >
const LUMatrixType& Eigen::UmfPackLU< MatrixType_ >::matrixL ( ) const
inline

Definition at line 351 of file UmfPackSupport.h.

352  {
354  return m_l;
355  }
void extractData() const

◆ matrixU()

template<typename MatrixType_ >
const LUMatrixType& Eigen::UmfPackLU< MatrixType_ >::matrixU ( ) const
inline

Definition at line 357 of file UmfPackSupport.h.

358  {
360  return m_u;
361  }

◆ permutationP()

template<typename MatrixType_ >
const IntColVectorType& Eigen::UmfPackLU< MatrixType_ >::permutationP ( ) const
inline

Definition at line 363 of file UmfPackSupport.h.

364  {
366  return m_p;
367  }

◆ permutationQ()

template<typename MatrixType_ >
const IntRowVectorType& Eigen::UmfPackLU< MatrixType_ >::permutationQ ( ) const
inline

Definition at line 369 of file UmfPackSupport.h.

370  {
372  return m_q;
373  }

◆ printUmfpackControl()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::printUmfpackControl ( )
inline

Prints the current UmfPack control settings.

See also
umfpackControl()

Definition at line 461 of file UmfPackSupport.h.

462  {
464  }
void umfpack_report_control(double control[UMFPACK_CONTROL], double, int)

◆ printUmfpackInfo()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::printUmfpackInfo ( )
inline

Prints statistics collected by UmfPack.

See also
analyzePattern(), compute()

Definition at line 470 of file UmfPackSupport.h.

471  {
472  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
474  }
void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], double, int)

◆ printUmfpackStatus()

template<typename MatrixType_ >
void Eigen::UmfPackLU< MatrixType_ >::printUmfpackStatus ( )
inline

Prints the status of the previous factorization operation performed by UmfPack (symbolic or numerical factorization).

See also
analyzePattern(), compute()

Definition at line 480 of file UmfPackSupport.h.

480  {
481  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
483  }
void umfpack_report_status(double control[UMFPACK_CONTROL], int status, double, int)

◆ rows()

template<typename MatrixType_ >
Index Eigen::UmfPackLU< MatrixType_ >::rows ( void  ) const
inline

Definition at line 337 of file UmfPackSupport.h.

337 { return mp_matrix.rows(); }

◆ umfpackControl() [1/2]

template<typename MatrixType_ >
UmfpackControl& Eigen::UmfPackLU< MatrixType_ >::umfpackControl ( )
inline

Provides access to the control settings array used by UmfPack.

If this array contains NaN's, the default values are used.

See UMFPACK documentation for details.

Definition at line 434 of file UmfPackSupport.h.

435  {
436  return m_control;
437  }

◆ umfpackControl() [2/2]

template<typename MatrixType_ >
const UmfpackControl& Eigen::UmfPackLU< MatrixType_ >::umfpackControl ( ) const
inline

Provides access to the control settings array used by UmfPack.

If this array contains NaN's, the default values are used.

See UMFPACK documentation for details.

Definition at line 423 of file UmfPackSupport.h.

424  {
425  return m_control;
426  }

◆ umfpackFactorizeReturncode()

template<typename MatrixType_ >
int Eigen::UmfPackLU< MatrixType_ >::umfpackFactorizeReturncode ( ) const
inline

Provides the return status code returned by UmfPack during the numeric factorization.

See also
factorize(), compute()

Definition at line 411 of file UmfPackSupport.h.

412  {
413  eigen_assert(m_numeric && "UmfPackLU: you must first call factorize()");
414  return m_fact_errorCode;
415  }

Member Data Documentation

◆ m_analysisIsOk

template<typename MatrixType_ >
int Eigen::UmfPackLU< MatrixType_ >::m_analysisIsOk
protected

Definition at line 565 of file UmfPackSupport.h.

◆ m_control

template<typename MatrixType_ >
UmfpackControl Eigen::UmfPackLU< MatrixType_ >::m_control
protected

Definition at line 550 of file UmfPackSupport.h.

◆ m_dummy

template<typename MatrixType_ >
UmfpackMatrixType Eigen::UmfPackLU< MatrixType_ >::m_dummy
protected

Definition at line 557 of file UmfPackSupport.h.

◆ m_extractedDataAreDirty

template<typename MatrixType_ >
bool Eigen::UmfPackLU< MatrixType_ >::m_extractedDataAreDirty
mutableprotected

Definition at line 566 of file UmfPackSupport.h.

◆ m_fact_errorCode

template<typename MatrixType_ >
StorageIndex Eigen::UmfPackLU< MatrixType_ >::m_fact_errorCode
protected

Definition at line 549 of file UmfPackSupport.h.

◆ m_factorizationIsOk

template<typename MatrixType_ >
int Eigen::UmfPackLU< MatrixType_ >::m_factorizationIsOk
protected

Definition at line 564 of file UmfPackSupport.h.

◆ m_info

template<typename MatrixType_ >
ComputationInfo Eigen::UmfPackLU< MatrixType_ >::m_info
mutableprotected

Definition at line 563 of file UmfPackSupport.h.

◆ m_l

template<typename MatrixType_ >
LUMatrixType Eigen::UmfPackLU< MatrixType_ >::m_l
mutableprotected

Definition at line 548 of file UmfPackSupport.h.

◆ m_numeric

template<typename MatrixType_ >
void* Eigen::UmfPackLU< MatrixType_ >::m_numeric
protected

Definition at line 560 of file UmfPackSupport.h.

◆ m_p

template<typename MatrixType_ >
IntColVectorType Eigen::UmfPackLU< MatrixType_ >::m_p
mutableprotected

Definition at line 554 of file UmfPackSupport.h.

◆ m_q

template<typename MatrixType_ >
IntRowVectorType Eigen::UmfPackLU< MatrixType_ >::m_q
mutableprotected

Definition at line 555 of file UmfPackSupport.h.

◆ m_symbolic

template<typename MatrixType_ >
void* Eigen::UmfPackLU< MatrixType_ >::m_symbolic
protected

Definition at line 561 of file UmfPackSupport.h.

◆ m_u

template<typename MatrixType_ >
LUMatrixType Eigen::UmfPackLU< MatrixType_ >::m_u
mutableprotected

Definition at line 553 of file UmfPackSupport.h.

◆ m_umfpackInfo

template<typename MatrixType_ >
UmfpackInfo Eigen::UmfPackLU< MatrixType_ >::m_umfpackInfo
mutableprotected

Definition at line 551 of file UmfPackSupport.h.

◆ mp_matrix

template<typename MatrixType_ >
UmfpackMatrixRef Eigen::UmfPackLU< MatrixType_ >::mp_matrix
protected

Definition at line 558 of file UmfPackSupport.h.


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