WPILibC++ 2023.4.3-108-ge5452e3
Eigen::SparseLU< _MatrixType, _OrderingType > Class Template Reference

Sparse supernodal LU factorization for general matrices. More...

#include </home/runner/work/allwpilib/allwpilib/wpimath/src/main/native/thirdparty/eigen/include/Eigen/src/SparseLU/SparseLU.h>

Inheritance diagram for Eigen::SparseLU< _MatrixType, _OrderingType >:
Eigen::SparseSolverBase< SparseLU< _MatrixType, _OrderingType > > Eigen::internal::SparseLUImpl< _MatrixType::Scalar, _MatrixType::StorageIndex > Eigen::internal::noncopyable

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef _OrderingType OrderingType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexNCMatrix
 
typedef internal::MappedSuperNodalMatrix< Scalar, StorageIndexSCMatrix
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef PermutationMatrix< Dynamic, Dynamic, StorageIndexPermutationType
 
typedef internal::SparseLUImpl< Scalar, StorageIndexBase
 
- Public Types inherited from Eigen::internal::SparseLUImpl< _MatrixType::Scalar, _MatrixType::StorageIndex >
typedef Matrix< _MatrixType::Scalar, Dynamic, 1 > ScalarVector
 
typedef Matrix< _MatrixType::StorageIndex, Dynamic, 1 > IndexVector
 
typedef Matrix< _MatrixType::Scalar, Dynamic, Dynamic, ColMajor > ScalarMatrix
 
typedef Map< ScalarMatrix, 0, OuterStride<> > MappedMatrixBlock
 
typedef ScalarVector::RealScalar RealScalar
 
typedef Ref< Matrix< _MatrixType::Scalar, Dynamic, 1 > > BlockScalarVector
 
typedef Ref< Matrix< _MatrixType::StorageIndex, Dynamic, 1 > > BlockIndexVector
 
typedef LU_GlobalLU_t< IndexVector, ScalarVectorGlobalLU_t
 
typedef SparseMatrix< _MatrixType::Scalar, ColMajor, _MatrixType::StorageIndex > MatrixType
 

Public Member Functions

 SparseLU ()
 
 SparseLU (const MatrixType &matrix)
 
 ~SparseLU ()
 
void analyzePattern (const MatrixType &matrix)
 Compute the column permutation to minimize the fill-in. More...
 
void factorize (const MatrixType &matrix)
 
void simplicialfactorize (const MatrixType &matrix)
 
void compute (const MatrixType &matrix)
 Compute the symbolic and numeric factorization of the input sparse matrix. More...
 
const SparseLUTransposeView< false, SparseLU< _MatrixType, _OrderingType > > transpose ()
 
const SparseLUTransposeView< true, SparseLU< _MatrixType, _OrderingType > > adjoint ()
 
Index rows () const
 
Index cols () const
 
void isSymmetric (bool sym)
 Indicate that the pattern of the input matrix is symmetric. More...
 
SparseLUMatrixLReturnType< SCMatrixmatrixL () const
 
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU () const
 
const PermutationTyperowsPermutation () const
 
const PermutationTypecolsPermutation () const
 
void setPivotThreshold (const RealScalar &thresh)
 Set the threshold used for a diagonal entry to be an acceptable pivot. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
std::string lastErrorMessage () const
 
template<typename Rhs , typename Dest >
bool _solve_impl (const MatrixBase< Rhs > &B, MatrixBase< Dest > &X_base) const
 
Scalar absDeterminant ()
 
Scalar logAbsDeterminant () const
 
Scalar signDeterminant ()
 
Scalar determinant ()
 
Index nnzL () const
 
Index nnzU () const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SparseLU< _MatrixType, _OrderingType > >
 SparseSolverBase ()
 Default constructor. More...
 
 ~SparseSolverBase ()
 
SparseLU< _MatrixType, _OrderingType > & derived ()
 
const SparseLU< _MatrixType, _OrderingType > & derived () const
 
const Solve< SparseLU< _MatrixType, _OrderingType >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SparseLU< _MatrixType, _OrderingType >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< SparseLU< _MatrixType, _OrderingType > > APIBase
 

Protected Member Functions

void initperfvalues ()
 
- Protected Member Functions inherited from Eigen::internal::SparseLUImpl< _MatrixType::Scalar, _MatrixType::StorageIndex >
Index expand (VectorType &vec, Index &length, Index nbElts, Index keep_prev, Index &num_expansions)
 Expand the existing storage to accommodate more fill-ins. More...
 
Index memInit (Index m, Index n, Index annz, Index lwork, Index fillratio, Index panel_size, GlobalLU_t &glu)
 Allocate various working space for the numerical factorization phase. More...
 
Index memXpand (VectorType &vec, Index &maxlen, Index nbElts, MemType memtype, Index &num_expansions)
 Expand the existing storage. More...
 
void heap_relax_snode (const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
 Identify the initial relaxed supernodes. More...
 
void relax_snode (const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
 Identify the initial relaxed supernodes. More...
 
Index snode_dfs (const Index jcol, const Index kcol, const MatrixType &mat, IndexVector &xprune, IndexVector &marker, GlobalLU_t &glu)
 
Index snode_bmod (const Index jcol, const Index fsupc, ScalarVector &dense, GlobalLU_t &glu)
 
Index pivotL (const Index jcol, const RealScalar &diagpivotthresh, IndexVector &perm_r, IndexVector &iperm_c, Index &pivrow, GlobalLU_t &glu)
 Performs the numerical pivotin on the current column of L, and the CDIV operation. More...
 
void dfs_kernel (const _MatrixType::StorageIndex jj, IndexVector &perm_r, Index &nseg, IndexVector &panel_lsub, IndexVector &segrep, Ref< IndexVector > repfnz_col, IndexVector &xprune, Ref< IndexVector > marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu, Index &nextl_col, Index krow, Traits &traits)
 
void panel_dfs (const Index m, const Index w, const Index jcol, MatrixType &A, IndexVector &perm_r, Index &nseg, ScalarVector &dense, IndexVector &panel_lsub, IndexVector &segrep, IndexVector &repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
 Performs a symbolic factorization on a panel of columns [jcol, jcol+w) More...
 
void panel_bmod (const Index m, const Index w, const Index jcol, const Index nseg, ScalarVector &dense, ScalarVector &tempv, IndexVector &segrep, IndexVector &repfnz, GlobalLU_t &glu)
 Performs numeric block updates (sup-panel) in topological order. More...
 
Index column_dfs (const Index m, const Index jcol, IndexVector &perm_r, Index maxsuper, Index &nseg, BlockIndexVector lsub_col, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
 Performs a symbolic factorization on column jcol and decide the supernode boundary. More...
 
Index column_bmod (const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector &tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t &glu)
 Performs numeric block updates (sup-col) in topological order. More...
 
Index copy_to_ucol (const Index jcol, const Index nseg, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &perm_r, BlockScalarVector dense, GlobalLU_t &glu)
 Performs numeric block updates (sup-col) in topological order. More...
 
void pruneL (const Index jcol, const IndexVector &perm_r, const Index pivrow, const Index nseg, const IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, GlobalLU_t &glu)
 Prunes the L-structure. More...
 
void countnz (const Index n, Index &nnzL, Index &nnzU, GlobalLU_t &glu)
 Count Nonzero elements in the factors. More...
 
void fixupL (const Index n, const IndexVector &perm_r, GlobalLU_t &glu)
 Fix up the data storage lsub for L-subscripts. More...
 

Protected Attributes

ComputationInfo m_info
 
bool m_factorizationIsOk
 
bool m_analysisIsOk
 
std::string m_lastError
 
NCMatrix m_mat
 
SCMatrix m_Lstore
 
MappedSparseMatrix< Scalar, ColMajor, StorageIndexm_Ustore
 
PermutationType m_perm_c
 
PermutationType m_perm_r
 
IndexVector m_etree
 
Base::GlobalLU_t m_glu
 
bool m_symmetricmode
 
internal::perfvalues m_perfv
 
RealScalar m_diagpivotthresh
 
Index m_nnzL
 
Index m_nnzU
 
Index m_detPermR
 
Index m_detPermC
 
- Protected Attributes inherited from Eigen::SparseSolverBase< SparseLU< _MatrixType, _OrderingType > >
bool m_isInitialized
 

Detailed Description

template<typename _MatrixType, typename _OrderingType>
class Eigen::SparseLU< _MatrixType, _OrderingType >

Sparse supernodal LU factorization for general matrices.

This class implements the supernodal LU factorization for general matrices. It uses the main techniques from the sequential SuperLU package (http://crd-legacy.lbl.gov/~xiaoye/SuperLU/). It handles transparently real and complex arithmetic with single and double precision, depending on the scalar type of your input matrix. The code has been optimized to provide BLAS-3 operations during supernode-panel updates. It benefits directly from the built-in high-performant Eigen BLAS routines. Moreover, when the size of a supernode is very small, the BLAS calls are avoided to enable a better optimization from the compiler. For best performance, you should compile it with NDEBUG flag to avoid the numerous bounds checking on vectors.

An important parameter of this class is the ordering method. It is used to reorder the columns (and eventually the rows) of the matrix to reduce the number of new elements that are created during numerical factorization. The cheapest method available is COLAMD. See the OrderingMethods module for the list of built-in and external ordering methods.

Simple example with key steps

VectorXd x(n), b(n);
// fill A and b;
// Compute the ordering permutation vector from the structural pattern of A
solver.analyzePattern(A);
// Compute the numerical factorization
solver.factorize(A);
//Use the factors to solve the linear system
x = solver.solve(b);
Definition: Ordering.h:115
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:132
void factorize(const MatrixType &matrix)
Definition: SparseLU.h:595
void analyzePattern(const MatrixType &matrix)
Compute the column permutation to minimize the fill-in.
Definition: SparseLU.h:510
A versatible sparse matrix representation.
Definition: SparseMatrix.h:98
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:88
b
Definition: data.h:44
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.
Note
Unlike the initial SuperLU implementation, there is no step to equilibrate the matrix. For badly scaled matrices, this step can be useful to reduce the pivoting during factorization. If this is the case for your matrices, you can try the basic scaling method at "unsupported/Eigen/src/IterativeSolvers/Scaling.h"
Template Parameters
_MatrixTypeThe type of the sparse matrix. It must be a column-major SparseMatrix<>
_OrderingTypeThe ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD

\implsparsesolverconcept

See also
TutorialSparseSolverConcept
OrderingMethods module

Member Typedef Documentation

◆ APIBase

template<typename _MatrixType , typename _OrderingType >
typedef SparseSolverBase<SparseLU<_MatrixType,_OrderingType> > Eigen::SparseLU< _MatrixType, _OrderingType >::APIBase
protected

◆ Base

template<typename _MatrixType , typename _OrderingType >
typedef internal::SparseLUImpl<Scalar, StorageIndex> Eigen::SparseLU< _MatrixType, _OrderingType >::Base

◆ IndexVector

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<StorageIndex,Dynamic,1> Eigen::SparseLU< _MatrixType, _OrderingType >::IndexVector

◆ MatrixType

template<typename _MatrixType , typename _OrderingType >
typedef _MatrixType Eigen::SparseLU< _MatrixType, _OrderingType >::MatrixType

◆ NCMatrix

template<typename _MatrixType , typename _OrderingType >
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SparseLU< _MatrixType, _OrderingType >::NCMatrix

◆ OrderingType

template<typename _MatrixType , typename _OrderingType >
typedef _OrderingType Eigen::SparseLU< _MatrixType, _OrderingType >::OrderingType

◆ PermutationType

template<typename _MatrixType , typename _OrderingType >
typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> Eigen::SparseLU< _MatrixType, _OrderingType >::PermutationType

◆ RealScalar

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::RealScalar Eigen::SparseLU< _MatrixType, _OrderingType >::RealScalar

◆ Scalar

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::Scalar

◆ ScalarVector

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<Scalar,Dynamic,1> Eigen::SparseLU< _MatrixType, _OrderingType >::ScalarVector

◆ SCMatrix

template<typename _MatrixType , typename _OrderingType >
typedef internal::MappedSuperNodalMatrix<Scalar, StorageIndex> Eigen::SparseLU< _MatrixType, _OrderingType >::SCMatrix

◆ StorageIndex

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::StorageIndex Eigen::SparseLU< _MatrixType, _OrderingType >::StorageIndex

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType , typename _OrderingType >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Constructor & Destructor Documentation

◆ SparseLU() [1/2]

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::SparseLU ( )
inline

◆ SparseLU() [2/2]

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::SparseLU ( const MatrixType matrix)
inlineexplicit

◆ ~SparseLU()

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::~SparseLU ( )
inline

Member Function Documentation

◆ _solve_impl()

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs , typename Dest >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::_solve_impl ( const MatrixBase< Rhs > &  B,
MatrixBase< Dest > &  X_base 
) const
inline

◆ absDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::absDeterminant ( )
inline
Returns
the absolute value of the determinant of the matrix of which *this is the QR decomposition.
Warning
a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow. One way to work around that is to use logAbsDeterminant() instead.
See also
logAbsDeterminant(), signDeterminant()

◆ adjoint()

template<typename _MatrixType , typename _OrderingType >
const SparseLUTransposeView< true, SparseLU< _MatrixType, _OrderingType > > Eigen::SparseLU< _MatrixType, _OrderingType >::adjoint ( )
inline
Returns
an expression of the adjoint of the factored matrix

A typical usage is to solve for the adjoint problem A' x = b:

solver.compute(A);
x = solver.adjoint().solve(b);

For real scalar types, this function is equivalent to transpose().

See also
transpose(), solve()

◆ analyzePattern()

template<typename MatrixType , typename OrderingType >
void Eigen::SparseLU< MatrixType, OrderingType >::analyzePattern ( const MatrixType mat)

Compute the column permutation to minimize the fill-in.

  • Apply this permutation to the input matrix -
  • Compute the column elimination tree on the permuted matrix
  • Postorder the elimination tree and the column permutation

◆ cols()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::cols ( ) const
inline

◆ colsPermutation()

template<typename _MatrixType , typename _OrderingType >
const PermutationType & Eigen::SparseLU< _MatrixType, _OrderingType >::colsPermutation ( ) const
inline
Returns
a reference to the column matrix permutation \( P_c^T \) such that \(P_r A P_c^T = L U\)
See also
rowsPermutation()

◆ compute()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::compute ( const MatrixType matrix)
inline

Compute the symbolic and numeric factorization of the input sparse matrix.

The input matrix should be in column-major storage.

◆ determinant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::determinant ( )
inline
Returns
The determinant of the matrix.
See also
absDeterminant(), logAbsDeterminant()

◆ factorize()

template<typename MatrixType , typename OrderingType >
void Eigen::SparseLU< MatrixType, OrderingType >::factorize ( const MatrixType matrix)
  • Numerical factorization
  • Interleaved with the symbolic factorization On exit, info is

    = 0: successful factorization

‍0: if info = i, and i is

<= A->ncol: U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

> A->ncol: number of bytes allocated when memory allocation failure occurred, plus A->ncol. If lwork = -1, it is the estimated amount of space needed, plus A->ncol.

◆ info()

template<typename _MatrixType , typename _OrderingType >
ComputationInfo Eigen::SparseLU< _MatrixType, _OrderingType >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the LU factorization reports a problem, zero diagonal for instance InvalidInput if the input matrix is invalid
See also
iparm()

◆ initperfvalues()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::initperfvalues ( )
inlineprotected

◆ isSymmetric()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::isSymmetric ( bool  sym)
inline

Indicate that the pattern of the input matrix is symmetric.

◆ lastErrorMessage()

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseLU< _MatrixType, _OrderingType >::lastErrorMessage ( ) const
inline
Returns
A string describing the type of error

◆ logAbsDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::logAbsDeterminant ( ) const
inline
Returns
the natural log of the absolute value of the determinant of the matrix of which **this is the QR decomposition
Note
This method is useful to work around the risk of overflow/underflow that's inherent to the determinant computation.
See also
absDeterminant(), signDeterminant()

◆ matrixL()

template<typename _MatrixType , typename _OrderingType >
SparseLUMatrixLReturnType< SCMatrix > Eigen::SparseLU< _MatrixType, _OrderingType >::matrixL ( ) const
inline
Returns
an expression of the matrix L, internally stored as supernodes The only operation available with this expression is the triangular solve
y = b; matrixL().solveInPlace(y);
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:243
const Scalar & y
Definition: MathFunctions.h:821

◆ matrixU()

template<typename _MatrixType , typename _OrderingType >
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > Eigen::SparseLU< _MatrixType, _OrderingType >::matrixU ( ) const
inline
Returns
an expression of the matrix U, The only operation available with this expression is the triangular solve
y = b; matrixU().solveInPlace(y);
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition: SparseLU.h:253

◆ nnzL()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::nnzL ( ) const
inline

◆ nnzU()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::nnzU ( ) const
inline

◆ rows()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::rows ( ) const
inline

◆ rowsPermutation()

template<typename _MatrixType , typename _OrderingType >
const PermutationType & Eigen::SparseLU< _MatrixType, _OrderingType >::rowsPermutation ( ) const
inline
Returns
a reference to the row matrix permutation \( P_r \) such that \(P_r A P_c^T = L U\)
See also
colsPermutation()

◆ setPivotThreshold()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::setPivotThreshold ( const RealScalar thresh)
inline

Set the threshold used for a diagonal entry to be an acceptable pivot.

◆ signDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::signDeterminant ( )
inline
Returns
A number representing the sign of the determinant
See also
absDeterminant(), logAbsDeterminant()

◆ simplicialfactorize()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::simplicialfactorize ( const MatrixType matrix)

◆ transpose()

template<typename _MatrixType , typename _OrderingType >
const SparseLUTransposeView< false, SparseLU< _MatrixType, _OrderingType > > Eigen::SparseLU< _MatrixType, _OrderingType >::transpose ( )
inline
Returns
an expression of the transposed of the factored matrix.

A typical usage is to solve for the transposed problem A^T x = b:

solver.compute(A);
x = solver.transpose().solve(b);
See also
adjoint(), solve()

Member Data Documentation

◆ m_analysisIsOk

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_analysisIsOk
protected

◆ m_detPermC

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_detPermC
protected

◆ m_detPermR

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_detPermR
protected

◆ m_diagpivotthresh

template<typename _MatrixType , typename _OrderingType >
RealScalar Eigen::SparseLU< _MatrixType, _OrderingType >::m_diagpivotthresh
protected

◆ m_etree

template<typename _MatrixType , typename _OrderingType >
IndexVector Eigen::SparseLU< _MatrixType, _OrderingType >::m_etree
protected

◆ m_factorizationIsOk

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_factorizationIsOk
protected

◆ m_glu

template<typename _MatrixType , typename _OrderingType >
Base::GlobalLU_t Eigen::SparseLU< _MatrixType, _OrderingType >::m_glu
protected

◆ m_info

template<typename _MatrixType , typename _OrderingType >
ComputationInfo Eigen::SparseLU< _MatrixType, _OrderingType >::m_info
mutableprotected

◆ m_lastError

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseLU< _MatrixType, _OrderingType >::m_lastError
protected

◆ m_Lstore

template<typename _MatrixType , typename _OrderingType >
SCMatrix Eigen::SparseLU< _MatrixType, _OrderingType >::m_Lstore
protected

◆ m_mat

template<typename _MatrixType , typename _OrderingType >
NCMatrix Eigen::SparseLU< _MatrixType, _OrderingType >::m_mat
protected

◆ m_nnzL

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_nnzL
protected

◆ m_nnzU

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_nnzU
protected

◆ m_perfv

template<typename _MatrixType , typename _OrderingType >
internal::perfvalues Eigen::SparseLU< _MatrixType, _OrderingType >::m_perfv
protected

◆ m_perm_c

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseLU< _MatrixType, _OrderingType >::m_perm_c
protected

◆ m_perm_r

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseLU< _MatrixType, _OrderingType >::m_perm_r
protected

◆ m_symmetricmode

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_symmetricmode
protected

◆ m_Ustore

template<typename _MatrixType , typename _OrderingType >
MappedSparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SparseLU< _MatrixType, _OrderingType >::m_Ustore
protected

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