Eigen::PastixLU< MatrixType_, IsStrSym > Class Template Reference

Interface to the PaStix solver. More...

#include <PaStiXSupport.h>

+ Inheritance diagram for Eigen::PastixLU< MatrixType_, IsStrSym >:

Public Types

typedef MatrixType_ MatrixType
 
typedef PastixBase< PastixLU< MatrixType > > Base
 
typedef Base::ColSpMatrix ColSpMatrix
 
typedef MatrixType::StorageIndex StorageIndex
 
- Public Types inherited from Eigen::PastixBase< PastixLU< MatrixType_ > >
enum  
 
typedef internal::pastix_traits< PastixLU< MatrixType_ > >::MatrixType MatrixType_
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 
typedef SparseMatrix< Scalar, ColMajor > ColSpMatrix
 

Public Member Functions

 PastixLU ()
 
 PastixLU (const MatrixType &matrix)
 
void compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
- Public Member Functions inherited from Eigen::PastixBase< PastixLU< MatrixType_ > >
 PastixBase ()
 
 ~PastixBase ()
 
bool _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &x) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
Array< StorageIndex, IPARM_SIZE, 1 > & iparm ()
 
intiparm (int idxparam)
 
Array< double, DPARM_SIZE, 1 > & dparm ()
 
doubledparm (int idxparam)
 
Index cols () const
 
Index rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
- 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 Member Functions

void init ()
 
void grabMatrix (const MatrixType &matrix, ColSpMatrix &out)
 
- Protected Member Functions inherited from Eigen::PastixBase< PastixLU< MatrixType_ > >
void init ()
 
void analyzePattern (ColSpMatrix &mat)
 
void factorize (ColSpMatrix &mat)
 
void clean ()
 
void compute (ColSpMatrix &mat)
 
PastixLU< MatrixType_ > & derived ()
 
const PastixLU< MatrixType_ > & derived () const
 

Protected Attributes

ColSpMatrix m_transposedStructure
 
bool m_structureIsUptodate
 
Array< double, DPARM_SIZE, 1 > m_dparm
 
Array< int, IPARM_SIZE, 1 > m_iparm
 
- Protected Attributes inherited from Eigen::PastixBase< PastixLU< MatrixType_ > >
int m_initisOk
 
int m_analysisIsOk
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
pastix_data_t * m_pastixdata
 
int m_comm
 
Array< int, IPARM_SIZE, 1 > m_iparm
 
Array< double, DPARM_SIZE, 1 > m_dparm
 
Matrix< StorageIndex, Dynamic, 1 > m_perm
 
Matrix< StorageIndex, Dynamic, 1 > m_invp
 
int m_size
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Additional Inherited Members

- Protected Types inherited from Eigen::PastixBase< PastixLU< MatrixType_ > >
typedef SparseSolverBase< PastixLU< MatrixType_ > > Base
 

Detailed Description

template<typename MatrixType_, bool IsStrSym>
class Eigen::PastixLU< MatrixType_, IsStrSym >

Interface to the PaStix solver.

Sparse direct LU solver based on PaStiX library.

This class is used to solve the linear systems A.X = B via the PaStix library. The matrix can be either real or complex, symmetric or not.

See also
TutorialSparseDirectSolvers

This class is used to solve the linear systems A.X = B with a supernodal LU factorization in the PaStiX library. The matrix A should be squared and nonsingular PaStiX requires that the matrix A has a symmetric structural pattern. This interface can symmetrize the input matrix otherwise. The vectors or matrices X and B can be either dense or sparse.

Template Parameters
MatrixType_the type of the sparse matrix A, it must be a SparseMatrix<>
IsStrSymIndicates if the input matrix has a symmetric pattern, default is false NOTE : Note that if the analysis and factorization phase are called separately, the input matrix will be symmetrized at each call, hence it is advised to symmetrize the matrix in a end-user program and set IsStrSym to true

\implsparsesolverconcept

See also
Sparse solver concept, class SparseLU

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , bool IsStrSym>
typedef PastixBase<PastixLU<MatrixType> > Eigen::PastixLU< MatrixType_, IsStrSym >::Base

◆ ColSpMatrix

template<typename MatrixType_ , bool IsStrSym>
typedef Base::ColSpMatrix Eigen::PastixLU< MatrixType_, IsStrSym >::ColSpMatrix

◆ MatrixType

template<typename MatrixType_ , bool IsStrSym>
typedef MatrixType_ Eigen::PastixLU< MatrixType_, IsStrSym >::MatrixType

◆ StorageIndex

template<typename MatrixType_ , bool IsStrSym>
typedef MatrixType::StorageIndex Eigen::PastixLU< MatrixType_, IsStrSym >::StorageIndex

Constructor & Destructor Documentation

◆ PastixLU() [1/2]

template<typename MatrixType_ , bool IsStrSym>
Eigen::PastixLU< MatrixType_, IsStrSym >::PastixLU ( )
inline
406 : Base() { init(); }
void init()
Definition: PaStiXSupport.h:447
PastixBase< PastixLU< MatrixType > > Base
Definition: PaStiXSupport.h:401

References Eigen::PastixLU< MatrixType_, IsStrSym >::init().

◆ PastixLU() [2/2]

template<typename MatrixType_ , bool IsStrSym>
Eigen::PastixLU< MatrixType_, IsStrSym >::PastixLU ( const MatrixType matrix)
inlineexplicit
408  : Base() {
409  init();
410  compute(matrix);
411  }
void compute(const MatrixType &matrix)
Definition: PaStiXSupport.h:417
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::PastixLU< MatrixType_, IsStrSym >::compute(), Eigen::PastixLU< MatrixType_, IsStrSym >::init(), and matrix().

Member Function Documentation

◆ analyzePattern()

template<typename MatrixType_ , bool IsStrSym>
void Eigen::PastixLU< MatrixType_, IsStrSym >::analyzePattern ( const MatrixType matrix)
inline

Compute the LU symbolic factorization of matrix using its sparsity pattern. Several ordering methods can be used at this step. See the PaStiX user's manual. The result of this operation can be used with successive matrices having the same pattern as matrix

See also
factorize()
428  {
429  m_structureIsUptodate = false;
430  ColSpMatrix temp;
431  grabMatrix(matrix, temp);
432  Base::analyzePattern(temp);
433  }
void analyzePattern(ColSpMatrix &mat)
Definition: PaStiXSupport.h:304
bool m_structureIsUptodate
Definition: PaStiXSupport.h:477
void grabMatrix(const MatrixType &matrix, ColSpMatrix &out)
Definition: PaStiXSupport.h:453
Base::ColSpMatrix ColSpMatrix
Definition: PaStiXSupport.h:402

References Eigen::PastixBase< Derived >::analyzePattern(), Eigen::PastixLU< MatrixType_, IsStrSym >::grabMatrix(), Eigen::PastixLU< MatrixType_, IsStrSym >::m_structureIsUptodate, and matrix().

◆ compute()

template<typename MatrixType_ , bool IsStrSym>
void Eigen::PastixLU< MatrixType_, IsStrSym >::compute ( const MatrixType matrix)
inline

Compute the LU supernodal factorization of matrix. iparm and dparm can be used to tune the PaStiX parameters. see the PaStiX user's manual

See also
analyzePattern() factorize()
417  {
418  m_structureIsUptodate = false;
419  ColSpMatrix temp;
420  grabMatrix(matrix, temp);
421  Base::compute(temp);
422  }
void compute(ColSpMatrix &mat)
Definition: PaStiXSupport.h:294

References Eigen::PastixBase< Derived >::compute(), Eigen::PastixLU< MatrixType_, IsStrSym >::grabMatrix(), Eigen::PastixLU< MatrixType_, IsStrSym >::m_structureIsUptodate, and matrix().

Referenced by Eigen::PastixLU< MatrixType_, IsStrSym >::PastixLU().

◆ factorize()

template<typename MatrixType_ , bool IsStrSym>
void Eigen::PastixLU< MatrixType_, IsStrSym >::factorize ( const MatrixType matrix)
inline

Compute the LU supernodal factorization of matrix WARNING The matrix matrix should have the same structural pattern as the same used in the analysis phase.

See also
analyzePattern()
440  {
441  ColSpMatrix temp;
442  grabMatrix(matrix, temp);
443  Base::factorize(temp);
444  }
void factorize(ColSpMatrix &mat)
Definition: PaStiXSupport.h:330

References Eigen::PastixBase< Derived >::factorize(), Eigen::PastixLU< MatrixType_, IsStrSym >::grabMatrix(), and matrix().

◆ grabMatrix()

template<typename MatrixType_ , bool IsStrSym>
void Eigen::PastixLU< MatrixType_, IsStrSym >::grabMatrix ( const MatrixType matrix,
ColSpMatrix out 
)
inlineprotected
453  {
454  if (IsStrSym)
455  out = matrix;
456  else {
457  if (!m_structureIsUptodate) {
458  // update the transposed structure
459  m_transposedStructure = matrix.transpose();
460 
461  // Set the elements of the matrix to zero
462  for (Index j = 0; j < m_transposedStructure.outerSize(); ++j)
463  for (typename ColSpMatrix::InnerIterator it(m_transposedStructure, j); it; ++it) it.valueRef() = 0.0;
464 
465  m_structureIsUptodate = true;
466  }
467 
469  }
471  }
ColSpMatrix m_transposedStructure
Definition: PaStiXSupport.h:476
Index outerSize() const
Definition: SparseMatrix.h:166
Base::InnerIterator InnerIterator
Definition: SparseMatrix.h:138
void c_to_fortran_numbering(MatrixType &mat)
Definition: PaStiXSupport.h:132
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
std::ofstream out("Result.txt")
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Eigen::internal::c_to_fortran_numbering(), j, Eigen::PastixLU< MatrixType_, IsStrSym >::m_structureIsUptodate, Eigen::PastixLU< MatrixType_, IsStrSym >::m_transposedStructure, matrix(), out(), and Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::outerSize().

Referenced by Eigen::PastixLU< MatrixType_, IsStrSym >::analyzePattern(), Eigen::PastixLU< MatrixType_, IsStrSym >::compute(), and Eigen::PastixLU< MatrixType_, IsStrSym >::factorize().

◆ init()

template<typename MatrixType_ , bool IsStrSym>
void Eigen::PastixLU< MatrixType_, IsStrSym >::init ( )
inlineprotected
447  {
448  m_structureIsUptodate = false;
449  m_iparm(IPARM_SYM) = API_SYM_NO;
450  m_iparm(IPARM_FACTORIZATION) = API_FACT_LU;
451  }
Array< int, IPARM_SIZE, 1 > m_iparm
Definition: PaStiXSupport.h:250

References Eigen::PastixLU< MatrixType_, IsStrSym >::m_iparm, and Eigen::PastixLU< MatrixType_, IsStrSym >::m_structureIsUptodate.

Referenced by Eigen::PastixLU< MatrixType_, IsStrSym >::PastixLU().

Member Data Documentation

◆ m_dparm

template<typename MatrixType_ , bool IsStrSym>
Array<double, DPARM_SIZE, 1> Eigen::PastixBase< Derived >::m_dparm
mutableprotected

◆ m_iparm

template<typename MatrixType_ , bool IsStrSym>
Array<int, IPARM_SIZE, 1> Eigen::PastixBase< Derived >::m_iparm
mutableprotected

◆ m_structureIsUptodate

◆ m_transposedStructure

template<typename MatrixType_ , bool IsStrSym>
ColSpMatrix Eigen::PastixLU< MatrixType_, IsStrSym >::m_transposedStructure
protected

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