Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::SuperLUBase< _MatrixType, Derived > Class Template Reference

The base class for the direct and incomplete LU factorization of SuperLU. More...

#include <src/eigen/Eigen/src/SuperLUSupport/SuperLUSupport.h>

+ Inheritance diagram for Eigen::SuperLUBase< _MatrixType, Derived >:
+ Collaboration diagram for Eigen::SuperLUBase< _MatrixType, Derived >:

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 
typedef Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
 
typedef Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
 
typedef Map< PermutationMatrix< Dynamic, Dynamic, int > > PermutationMap
 
typedef SparseMatrix< ScalarLUMatrixType
 

Public Member Functions

 SuperLUBase ()
 
 ~SuperLUBase ()
 
Index rows () const
 
Index cols () const
 
superlu_options_t & options ()
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
void compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &)
 
template<typename Stream >
void dumpMemory (Stream &)
 
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

void initFactorization (const MatrixType &a)
 
void init ()
 
void extractData () const
 
void clearFactors ()
 
Derived & derived ()
 
const Derived & derived () const
 

Protected Attributes

LUMatrixType m_l
 
LUMatrixType m_u
 
IntColVectorType m_p
 
IntRowVectorType m_q
 
LUMatrixType m_matrix
 
SluMatrix m_sluA
 
SuperMatrix m_sluL
 
SuperMatrix m_sluU
 
SluMatrix m_sluB
 
SluMatrix m_sluX
 
SuperLUStat_t m_sluStat
 
superlu_options_t m_sluOptions
 
std::vector< int > m_sluEtree
 
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
 
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
 
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
 
Matrix< RealScalar, Dynamic, 1 > m_sluBerr
 
char m_sluEqued
 
ComputationInfo m_info
 
int m_factorizationIsOk
 
int m_analysisIsOk
 
bool m_extractedDataAreDirty
 
bool m_isInitialized
 

Private Member Functions

 SuperLUBase (SuperLUBase &)
 

Detailed Description

template<typename _MatrixType, typename Derived>
class Eigen::SuperLUBase< _MatrixType, Derived >

The base class for the direct and incomplete LU factorization of SuperLU.

Member Typedef Documentation

◆ Base

template<typename _MatrixType , typename Derived >
typedef SparseSolverBase<Derived> Eigen::SuperLUBase< _MatrixType, Derived >::Base
protected

◆ IntColVectorType

template<typename _MatrixType , typename Derived >
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> Eigen::SuperLUBase< _MatrixType, Derived >::IntColVectorType

◆ IntRowVectorType

template<typename _MatrixType , typename Derived >
typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> Eigen::SuperLUBase< _MatrixType, Derived >::IntRowVectorType

◆ LUMatrixType

template<typename _MatrixType , typename Derived >
typedef SparseMatrix<Scalar> Eigen::SuperLUBase< _MatrixType, Derived >::LUMatrixType

◆ MatrixType

template<typename _MatrixType , typename Derived >
typedef _MatrixType Eigen::SuperLUBase< _MatrixType, Derived >::MatrixType

◆ PermutationMap

template<typename _MatrixType , typename Derived >
typedef Map<PermutationMatrix<Dynamic,Dynamic,int> > Eigen::SuperLUBase< _MatrixType, Derived >::PermutationMap

◆ RealScalar

template<typename _MatrixType , typename Derived >
typedef MatrixType::RealScalar Eigen::SuperLUBase< _MatrixType, Derived >::RealScalar

◆ Scalar

template<typename _MatrixType , typename Derived >
typedef MatrixType::Scalar Eigen::SuperLUBase< _MatrixType, Derived >::Scalar

◆ StorageIndex

template<typename _MatrixType , typename Derived >
typedef MatrixType::StorageIndex Eigen::SuperLUBase< _MatrixType, Derived >::StorageIndex

◆ Vector

template<typename _MatrixType , typename Derived >
typedef Matrix<Scalar,Dynamic,1> Eigen::SuperLUBase< _MatrixType, Derived >::Vector

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType , typename Derived >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
333 {
334 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
335 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
336 };
@ MaxColsAtCompileTime
Definition SuperLUSupport.h:335
@ ColsAtCompileTime
Definition SuperLUSupport.h:334

Constructor & Destructor Documentation

◆ SuperLUBase() [1/2]

template<typename _MatrixType , typename Derived >
Eigen::SuperLUBase< _MatrixType, Derived >::SuperLUBase ( )
inline
340{}

◆ ~SuperLUBase()

template<typename _MatrixType , typename Derived >
Eigen::SuperLUBase< _MatrixType, Derived >::~SuperLUBase ( )
inline
343 {
344 clearFactors();
345 }
void clearFactors()
Definition SuperLUSupport.h:430

References Eigen::SuperLUBase< _MatrixType, Derived >::clearFactors().

+ Here is the call graph for this function:

◆ SuperLUBase() [2/2]

template<typename _MatrixType , typename Derived >
Eigen::SuperLUBase< _MatrixType, Derived >::SuperLUBase ( SuperLUBase< _MatrixType, Derived > &  )
inlineprivate
467{ }

Member Function Documentation

◆ _solve_impl()

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inlineinherited
112 {
113 internal::solve_sparse_through_dense_panels(derived(), b.derived(), dest.derived());
114 }
Derived & derived()
Definition SparseSolverBase.h:79
enable_if< Rhs::ColsAtCompileTime!=1 &&Dest::ColsAtCompileTime!=1 >::type solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition SparseSolverBase.h:23

References Eigen::SparseSolverBase< Derived >::derived(), Eigen::SparseMatrixBase< Derived >::derived(), and Eigen::internal::solve_sparse_through_dense_panels().

+ Here is the call graph for this function:

◆ analyzePattern()

template<typename _MatrixType , typename Derived >
void Eigen::SuperLUBase< _MatrixType, Derived >::analyzePattern ( const MatrixType )
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()
378 {
379 m_isInitialized = true;
380 m_info = Success;
381 m_analysisIsOk = true;
382 m_factorizationIsOk = false;
383 }
ComputationInfo m_info
Definition SuperLUSupport.h:461
int m_analysisIsOk
Definition SuperLUSupport.h:463
bool m_isInitialized
Definition SparseSolverBase.h:119
int m_factorizationIsOk
Definition SuperLUSupport.h:462
@ Success
Definition Constants.h:432

References Eigen::SuperLUBase< _MatrixType, Derived >::m_analysisIsOk, Eigen::SuperLUBase< _MatrixType, Derived >::m_factorizationIsOk, Eigen::SuperLUBase< _MatrixType, Derived >::m_info, Eigen::SuperLUBase< _MatrixType, Derived >::m_isInitialized, and Eigen::Success.

Referenced by Eigen::SuperLU< _MatrixType >::analyzePattern().

+ Here is the caller graph for this function:

◆ clearFactors()

template<typename _MatrixType , typename Derived >
void Eigen::SuperLUBase< _MatrixType, Derived >::clearFactors ( )
inlineprotected
431 {
432 if(m_sluL.Store)
433 Destroy_SuperNode_Matrix(&m_sluL);
434 if(m_sluU.Store)
435 Destroy_CompCol_Matrix(&m_sluU);
436
437 m_sluL.Store = 0;
438 m_sluU.Store = 0;
439
440 memset(&m_sluL,0,sizeof m_sluL);
441 memset(&m_sluU,0,sizeof m_sluU);
442 }
SuperMatrix m_sluL
Definition SuperLUSupport.h:452
SuperMatrix m_sluU
Definition SuperLUSupport.h:452

References Eigen::SuperLUBase< _MatrixType, Derived >::m_sluL, and Eigen::SuperLUBase< _MatrixType, Derived >::m_sluU.

Referenced by Eigen::SuperLUBase< _MatrixType, Derived >::~SuperLUBase(), and Eigen::SuperLUBase< _MatrixType, Derived >::initFactorization().

+ Here is the caller graph for this function:

◆ cols()

template<typename _MatrixType , typename Derived >
Index Eigen::SuperLUBase< _MatrixType, Derived >::cols ( ) const
inline
348{ return m_matrix.cols(); }
Index cols() const
Definition SparseMatrix.h:138
LUMatrixType m_matrix
Definition SuperLUSupport.h:450

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::cols(), and Eigen::SuperLUBase< _MatrixType, Derived >::m_matrix.

+ Here is the call graph for this function:

◆ compute()

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

Computes the sparse Cholesky decomposition of matrix

366 {
367 derived().analyzePattern(matrix);
368 derived().factorize(matrix);
369 }
Derived & derived()
Definition SparseSolverBase.h:79

References Eigen::SuperLUBase< _MatrixType, Derived >::derived().

Referenced by Eigen::SuperLU< _MatrixType >::SuperLU().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ derived() [1/2]

template<typename _MatrixType , typename Derived >
Derived & Eigen::SparseSolverBase< Derived >::derived ( )
inlineprotected
79{ return *static_cast<Derived*>(this); }

Referenced by Eigen::SuperLUBase< _MatrixType, Derived >::compute().

+ Here is the caller graph for this function:

◆ derived() [2/2]

template<typename _MatrixType , typename Derived >
const Derived & Eigen::SparseSolverBase< Derived >::derived ( ) const
inlineprotected
80{ return *static_cast<const Derived*>(this); }

◆ dumpMemory()

template<typename _MatrixType , typename Derived >
template<typename Stream >
void Eigen::SuperLUBase< _MatrixType, Derived >::dumpMemory ( Stream &  )
inline
387 {}

◆ extractData()

template<typename MatrixType , typename Derived >
void Eigen::SuperLUBase< MatrixType, Derived >::extractData
protected
708{
709 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for extracting factors, you must first call either compute() or analyzePattern()/factorize()");
711 {
712 int upper;
713 int fsupc, istart, nsupr;
714 int lastl = 0, lastu = 0;
715 SCformat *Lstore = static_cast<SCformat*>(m_sluL.Store);
716 NCformat *Ustore = static_cast<NCformat*>(m_sluU.Store);
717 Scalar *SNptr;
718
719 const Index size = m_matrix.rows();
720 m_l.resize(size,size);
721 m_l.resizeNonZeros(Lstore->nnz);
722 m_u.resize(size,size);
723 m_u.resizeNonZeros(Ustore->nnz);
724
725 int* Lcol = m_l.outerIndexPtr();
726 int* Lrow = m_l.innerIndexPtr();
727 Scalar* Lval = m_l.valuePtr();
728
729 int* Ucol = m_u.outerIndexPtr();
730 int* Urow = m_u.innerIndexPtr();
731 Scalar* Uval = m_u.valuePtr();
732
733 Ucol[0] = 0;
734 Ucol[0] = 0;
735
736 /* for each supernode */
737 for (int k = 0; k <= Lstore->nsuper; ++k)
738 {
739 fsupc = L_FST_SUPC(k);
740 istart = L_SUB_START(fsupc);
741 nsupr = L_SUB_START(fsupc+1) - istart;
742 upper = 1;
743
744 /* for each column in the supernode */
745 for (int j = fsupc; j < L_FST_SUPC(k+1); ++j)
746 {
747 SNptr = &((Scalar*)Lstore->nzval)[L_NZ_START(j)];
748
749 /* Extract U */
750 for (int i = U_NZ_START(j); i < U_NZ_START(j+1); ++i)
751 {
752 Uval[lastu] = ((Scalar*)Ustore->nzval)[i];
753 /* Matlab doesn't like explicit zero. */
754 if (Uval[lastu] != 0.0)
755 Urow[lastu++] = U_SUB(i);
756 }
757 for (int i = 0; i < upper; ++i)
758 {
759 /* upper triangle in the supernode */
760 Uval[lastu] = SNptr[i];
761 /* Matlab doesn't like explicit zero. */
762 if (Uval[lastu] != 0.0)
763 Urow[lastu++] = L_SUB(istart+i);
764 }
765 Ucol[j+1] = lastu;
766
767 /* Extract L */
768 Lval[lastl] = 1.0; /* unit diagonal */
769 Lrow[lastl++] = L_SUB(istart + upper - 1);
770 for (int i = upper; i < nsupr; ++i)
771 {
772 Lval[lastl] = SNptr[i];
773 /* Matlab doesn't like explicit zero. */
774 if (Lval[lastl] != 0.0)
775 Lrow[lastl++] = L_SUB(istart+i);
776 }
777 Lcol[j+1] = lastl;
778
779 ++upper;
780 } /* for j ... */
781
782 } /* for k ... */
783
784 // squeeze the matrices :
785 m_l.resizeNonZeros(lastl);
786 m_u.resizeNonZeros(lastu);
787
789 }
790}
#define eigen_assert(x)
Definition Macros.h:579
void resizeNonZeros(Index size)
Definition SparseMatrix.h:644
Index rows() const
Definition SparseMatrix.h:136
const StorageIndex * innerIndexPtr() const
Definition SparseMatrix.h:157
const StorageIndex * outerIndexPtr() const
Definition SparseMatrix.h:166
const Scalar * valuePtr() const
Definition SparseMatrix.h:148
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
bool m_extractedDataAreDirty
Definition SuperLUSupport.h:464
MatrixType::Scalar Scalar
Definition SuperLUSupport.h:325
LUMatrixType m_l
Definition SuperLUSupport.h:445
LUMatrixType m_u
Definition SuperLUSupport.h:446
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

References eigen_assert.

◆ info()

template<typename _MatrixType , typename Derived >
ComputationInfo Eigen::SuperLUBase< _MatrixType, Derived >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.
359 {
360 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
361 return m_info;
362 }

References eigen_assert, Eigen::SuperLUBase< _MatrixType, Derived >::m_info, and Eigen::SuperLUBase< _MatrixType, Derived >::m_isInitialized.

◆ init()

template<typename _MatrixType , typename Derived >
void Eigen::SuperLUBase< _MatrixType, Derived >::init ( )
inlineprotected
421 {
423 m_isInitialized = false;
424 m_sluL.Store = 0;
425 m_sluU.Store = 0;
426 }
@ InvalidInput
Definition Constants.h:439

References Eigen::InvalidInput, Eigen::SuperLUBase< _MatrixType, Derived >::m_info, Eigen::SuperLUBase< _MatrixType, Derived >::m_isInitialized, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluL, and Eigen::SuperLUBase< _MatrixType, Derived >::m_sluU.

Referenced by Eigen::SuperLU< _MatrixType >::init().

+ Here is the caller graph for this function:

◆ initFactorization()

template<typename _MatrixType , typename Derived >
void Eigen::SuperLUBase< _MatrixType, Derived >::initFactorization ( const MatrixType a)
inlineprotected
392 {
393 set_default_options(&this->m_sluOptions);
394
395 const Index size = a.rows();
396 m_matrix = a;
397
399 clearFactors();
400
401 m_p.resize(size);
402 m_q.resize(size);
403 m_sluRscale.resize(size);
404 m_sluCscale.resize(size);
405 m_sluEtree.resize(size);
406
407 // set empty B and X
408 m_sluB.setStorageType(SLU_DN);
410 m_sluB.Mtype = SLU_GE;
411 m_sluB.storage.values = 0;
412 m_sluB.nrow = 0;
413 m_sluB.ncol = 0;
414 m_sluB.storage.lda = internal::convert_index<int>(size);
415 m_sluX = m_sluB;
416
418 }
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
SluMatrix m_sluX
Definition SuperLUSupport.h:453
superlu_options_t m_sluOptions
Definition SuperLUSupport.h:455
SluMatrix m_sluB
Definition SuperLUSupport.h:453
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
Definition SuperLUSupport.h:457
std::vector< int > m_sluEtree
Definition SuperLUSupport.h:456
IntRowVectorType m_q
Definition SuperLUSupport.h:448
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
Definition SuperLUSupport.h:457
SluMatrix m_sluA
Definition SuperLUSupport.h:451
IntColVectorType m_p
Definition SuperLUSupport.h:447
SluMatrix asSluMatrix(MatrixType &mat)
Definition SuperLUSupport.h:291
struct Eigen::SluMatrix::@636 storage
void setScalarType()
Definition SuperLUSupport.h:158
void setStorageType(Stype_t t)
Definition SuperLUSupport.h:145

References Eigen::internal::asSluMatrix(), Eigen::SuperLUBase< _MatrixType, Derived >::clearFactors(), Eigen::SuperLUBase< _MatrixType, Derived >::m_extractedDataAreDirty, Eigen::SuperLUBase< _MatrixType, Derived >::m_matrix, Eigen::SuperLUBase< _MatrixType, Derived >::m_p, Eigen::SuperLUBase< _MatrixType, Derived >::m_q, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluA, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluB, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluCscale, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluEtree, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluOptions, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluRscale, Eigen::SuperLUBase< _MatrixType, Derived >::m_sluX, Eigen::PlainObjectBase< Derived >::resize(), Eigen::SluMatrix::setScalarType(), Eigen::SluMatrix::setStorageType(), and Eigen::SluMatrix::storage.

+ Here is the call graph for this function:

◆ options()

template<typename _MatrixType , typename Derived >
superlu_options_t & Eigen::SuperLUBase< _MatrixType, Derived >::options ( )
inline
Returns
a reference to the Super LU option object to configure the Super LU algorithms.
351{ return m_sluOptions; }

References Eigen::SuperLUBase< _MatrixType, Derived >::m_sluOptions.

◆ rows()

template<typename _MatrixType , typename Derived >
Index Eigen::SuperLUBase< _MatrixType, Derived >::rows ( ) const
inline
347{ return m_matrix.rows(); }

References Eigen::SuperLUBase< _MatrixType, Derived >::m_matrix, and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows().

+ Here is the call graph for this function:

◆ solve() [1/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const MatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
89 {
90 eigen_assert(m_isInitialized && "Solver is not initialized.");
91 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
92 return Solve<Derived, Rhs>(derived(), b.derived());
93 }
bool m_isInitialized
Definition SparseSolverBase.h:119
size_t rows(const T &raster)
Definition MarchingSquares.hpp:55

References Eigen::SparseSolverBase< Derived >::derived(), eigen_assert, and Eigen::SparseSolverBase< Derived >::m_isInitialized.

Referenced by igl::embree::bone_heat(), igl::Frame_field_deformer::compute_optimal_positions(), igl::eigs(), igl::copyleft::comiso::FrameInterpolator::interpolateSymmetric(), igl::slim::solve_weighted_arap(), and igl::copyleft::comiso::NRosyField::solveNoRoundings().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [2/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const SparseMatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
102 {
103 eigen_assert(m_isInitialized && "Solver is not initialized.");
104 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
105 return Solve<Derived, Rhs>(derived(), b.derived());
106 }

Member Data Documentation

◆ m_analysisIsOk

template<typename _MatrixType , typename Derived >
int Eigen::SuperLUBase< _MatrixType, Derived >::m_analysisIsOk
protected

◆ m_extractedDataAreDirty

template<typename _MatrixType , typename Derived >
bool Eigen::SuperLUBase< _MatrixType, Derived >::m_extractedDataAreDirty
mutableprotected

◆ m_factorizationIsOk

template<typename _MatrixType , typename Derived >
int Eigen::SuperLUBase< _MatrixType, Derived >::m_factorizationIsOk
protected

◆ m_info

◆ m_isInitialized

template<typename _MatrixType , typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_l

template<typename _MatrixType , typename Derived >
LUMatrixType Eigen::SuperLUBase< _MatrixType, Derived >::m_l
mutableprotected

◆ m_matrix

template<typename _MatrixType , typename Derived >
LUMatrixType Eigen::SuperLUBase< _MatrixType, Derived >::m_matrix
mutableprotected

◆ m_p

template<typename _MatrixType , typename Derived >
IntColVectorType Eigen::SuperLUBase< _MatrixType, Derived >::m_p
mutableprotected

◆ m_q

template<typename _MatrixType , typename Derived >
IntRowVectorType Eigen::SuperLUBase< _MatrixType, Derived >::m_q
mutableprotected

◆ m_sluA

template<typename _MatrixType , typename Derived >
SluMatrix Eigen::SuperLUBase< _MatrixType, Derived >::m_sluA
mutableprotected

◆ m_sluB

template<typename _MatrixType , typename Derived >
SluMatrix Eigen::SuperLUBase< _MatrixType, Derived >::m_sluB
mutableprotected

◆ m_sluBerr

template<typename _MatrixType , typename Derived >
Matrix<RealScalar,Dynamic,1> Eigen::SuperLUBase< _MatrixType, Derived >::m_sluBerr
protected

◆ m_sluCscale

template<typename _MatrixType , typename Derived >
Matrix<RealScalar,Dynamic,1> Eigen::SuperLUBase< _MatrixType, Derived >::m_sluCscale
protected

◆ m_sluEqued

template<typename _MatrixType , typename Derived >
char Eigen::SuperLUBase< _MatrixType, Derived >::m_sluEqued
mutableprotected

◆ m_sluEtree

template<typename _MatrixType , typename Derived >
std::vector<int> Eigen::SuperLUBase< _MatrixType, Derived >::m_sluEtree
mutableprotected

◆ m_sluFerr

template<typename _MatrixType , typename Derived >
Matrix<RealScalar,Dynamic,1> Eigen::SuperLUBase< _MatrixType, Derived >::m_sluFerr
mutableprotected

◆ m_sluL

template<typename _MatrixType , typename Derived >
SuperMatrix Eigen::SuperLUBase< _MatrixType, Derived >::m_sluL
mutableprotected

◆ m_sluOptions

template<typename _MatrixType , typename Derived >
superlu_options_t Eigen::SuperLUBase< _MatrixType, Derived >::m_sluOptions
mutableprotected

◆ m_sluRscale

template<typename _MatrixType , typename Derived >
Matrix<RealScalar,Dynamic,1> Eigen::SuperLUBase< _MatrixType, Derived >::m_sluRscale
mutableprotected

◆ m_sluStat

template<typename _MatrixType , typename Derived >
SuperLUStat_t Eigen::SuperLUBase< _MatrixType, Derived >::m_sluStat
mutableprotected

◆ m_sluU

template<typename _MatrixType , typename Derived >
SuperMatrix Eigen::SuperLUBase< _MatrixType, Derived >::m_sluU
protected

◆ m_sluX

template<typename _MatrixType , typename Derived >
SluMatrix Eigen::SuperLUBase< _MatrixType, Derived >::m_sluX
protected

◆ m_u

template<typename _MatrixType , typename Derived >
LUMatrixType Eigen::SuperLUBase< _MatrixType, Derived >::m_u
mutableprotected

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