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

A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library. More...

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

+ Inheritance diagram for Eigen::PastixLDLT< _MatrixType, _UpLo >:
+ Collaboration diagram for Eigen::PastixLDLT< _MatrixType, _UpLo >:

Public Types

enum  { UpLo = _UpLo }
 
typedef _MatrixType MatrixType
 
typedef PastixBase< PastixLDLT< MatrixType, _UpLo > > Base
 
typedef Base::ColSpMatrix ColSpMatrix
 
enum  
 
typedef internal::pastix_traits< PastixLDLT< _MatrixType, _UpLo > >::MatrixType _MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 

Public Member Functions

 PastixLDLT ()
 
 PastixLDLT (const MatrixType &matrix)
 
void compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
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 ()
 
int & iparm (int idxparam)
 
Array< double, DPARM_SIZE, 1 > & dparm ()
 
double & dparm (int idxparam)
 
Index cols () const
 
Index rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
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
 

Protected Member Functions

void init ()
 
void grabMatrix (const MatrixType &matrix, ColSpMatrix &out)
 
void analyzePattern (ColSpMatrix &mat)
 
void factorize (ColSpMatrix &mat)
 
void clean ()
 
void compute (ColSpMatrix &mat)
 
PastixLDLT< _MatrixType, _UpLo > & derived ()
 
const PastixLDLT< _MatrixType, _UpLo > & derived () const
 

Protected Attributes

Array< int, IPARM_SIZE, 1 > m_iparm
 
int m_initisOk
 
int m_analysisIsOk
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
pastix_data_t * m_pastixdata
 
int m_comm
 
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
 

Detailed Description

template<typename _MatrixType, int _UpLo>
class Eigen::PastixLDLT< _MatrixType, _UpLo >

A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library.

This class is used to solve the linear systems A.X = B via a LDL^T supernodal Cholesky factorization available in the PaStiX library. The matrix A should be symmetric and positive definite WARNING Selfadjoint complex matrices are not supported in the current version of PaStiX The vectors or matrices X and B can be either dense or sparse

Template Parameters
MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
UpLoThe part of the matrix to use : Lower or Upper. The default is Lower as required by PaStiX

\implsparsesolverconcept

See also
TutorialSparseSolverConcept, class SimplicialLDLT

Member Typedef Documentation

◆ _MatrixType

typedef internal::pastix_traits<PastixLDLT< _MatrixType, _UpLo > >::MatrixType Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::_MatrixType
inherited

◆ Base

template<typename _MatrixType , int _UpLo>
typedef PastixBase<PastixLDLT<MatrixType, _UpLo> > Eigen::PastixLDLT< _MatrixType, _UpLo >::Base

◆ ColSpMatrix

template<typename _MatrixType , int _UpLo>
typedef Base::ColSpMatrix Eigen::PastixLDLT< _MatrixType, _UpLo >::ColSpMatrix

◆ MatrixType

template<typename _MatrixType , int _UpLo>
typedef _MatrixType Eigen::PastixLDLT< _MatrixType, _UpLo >::MatrixType

◆ RealScalar

typedef MatrixType::RealScalar Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::RealScalar
inherited

◆ Scalar

typedef MatrixType::Scalar Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::Scalar
inherited

◆ StorageIndex

typedef MatrixType::StorageIndex Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::StorageIndex
inherited

◆ Vector

typedef Matrix<Scalar,Dynamic,1> Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::Vector
inherited

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
inherited
144 {
145 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
146 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
147 };
@ ColsAtCompileTime
Definition PaStiXSupport.h:145
@ MaxColsAtCompileTime
Definition PaStiXSupport.h:146

◆ anonymous enum

template<typename _MatrixType , int _UpLo>
anonymous enum
Enumerator
UpLo 
616{ UpLo = _UpLo };
@ UpLo
Definition PaStiXSupport.h:616

Constructor & Destructor Documentation

◆ PastixLDLT() [1/2]

template<typename _MatrixType , int _UpLo>
Eigen::PastixLDLT< _MatrixType, _UpLo >::PastixLDLT ( )
inline
617 :Base()
618 {
619 init();
620 }
PastixBase< PastixLDLT< MatrixType, _UpLo > > Base
Definition PaStiXSupport.h:612
void init()
Definition PaStiXSupport.h:661

References Eigen::PastixLDLT< _MatrixType, _UpLo >::init().

+ Here is the call graph for this function:

◆ PastixLDLT() [2/2]

template<typename _MatrixType , int _UpLo>
Eigen::PastixLDLT< _MatrixType, _UpLo >::PastixLDLT ( const MatrixType matrix)
inlineexplicit
622 :Base()
623 {
624 init();
625 compute(matrix);
626 }
void compute(const MatrixType &matrix)
Definition PaStiXSupport.h:631

References Eigen::PastixLDLT< _MatrixType, _UpLo >::compute(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::init().

+ Here is the call graph for this function:

Member Function Documentation

◆ _solve_impl() [1/2]

bool Eigen::PastixBase< Base >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  x 
) const
inherited
368{
369 eigen_assert(m_isInitialized && "The matrix should be factorized first");
370 EIGEN_STATIC_ASSERT((Dest::Flags&RowMajorBit)==0,
371 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
372 int rhs = 1;
373
374 x = b; /* on return, x is overwritten by the computed solution */
375
376 for (int i = 0; i < b.cols(); i++){
377 m_iparm[IPARM_START_TASK] = API_TASK_SOLVE;
378 m_iparm[IPARM_END_TASK] = API_TASK_REFINE;
379
380 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, internal::convert_index<int>(x.rows()), 0, 0, 0,
381 m_perm.data(), m_invp.data(), &x(0, i), rhs, m_iparm.data(), m_dparm.data());
382 }
383
384 // Check the returned error
385 m_info = m_iparm(IPARM_ERROR_NUMBER)==0 ? Success : NumericalIssue;
386
387 return m_iparm(IPARM_ERROR_NUMBER)==0;
388}
#define eigen_assert(x)
Definition Macros.h:579
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition StaticAssert.h:124
Array< double, DPARM_SIZE, 1 > m_dparm
Definition PaStiXSupport.h:248
Matrix< StorageIndex, Dynamic, 1 > m_invp
Definition PaStiXSupport.h:250
Matrix< StorageIndex, Dynamic, 1 > m_perm
Definition PaStiXSupport.h:249
Array< int, IPARM_SIZE, 1 > m_iparm
Definition PaStiXSupport.h:247
ComputationInfo m_info
Definition PaStiXSupport.h:244
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Definition PlainObjectBase.h:255
@ NumericalIssue
Definition Constants.h:434
@ Success
Definition Constants.h:432
TCoord< P > x(const P &p)
Definition geometry_traits.hpp:297

◆ _solve_impl() [2/2]

void Eigen::SparseSolverBase< PastixLDLT< _MatrixType, _UpLo > >::_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 }
PastixLDLT< _MatrixType, _UpLo > & derived()
Definition SparseSolverBase.h:79

◆ analyzePattern() [1/2]

void Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::analyzePattern ( ColSpMatrix mat)
protectedinherited
308{
309 eigen_assert(m_initisOk && "The initialization of PaSTiX failed");
310
311 // clean previous calls
312 if(m_size>0)
313 clean();
314
315 m_size = internal::convert_index<int>(mat.rows());
316 m_perm.resize(m_size);
317 m_invp.resize(m_size);
318
319 m_iparm(IPARM_START_TASK) = API_TASK_ORDERING;
320 m_iparm(IPARM_END_TASK) = API_TASK_ANALYSE;
321 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, m_size, mat.outerIndexPtr(), mat.innerIndexPtr(),
322 mat.valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
323
324 // Check the returned error
325 if(m_iparm(IPARM_ERROR_NUMBER))
326 {
328 m_analysisIsOk = false;
329 }
330 else
331 {
332 m_info = Success;
333 m_analysisIsOk = true;
334 }
335}
int m_analysisIsOk
Definition PaStiXSupport.h:242
void clean()
Definition PaStiXSupport.h:230
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279

◆ analyzePattern() [2/2]

template<typename _MatrixType , int _UpLo>
void Eigen::PastixLDLT< _MatrixType, _UpLo >::analyzePattern ( const MatrixType matrix)
inline

Compute the LDL^T symbolic factorization of matrix using its sparsity pattern The result of this operation can be used with successive matrices having the same pattern as matrix

See also
factorize()
643 {
644 ColSpMatrix temp;
645 grabMatrix(matrix, temp);
647 }
void analyzePattern(ColSpMatrix &mat)
Definition PaStiXSupport.h:307
Base::ColSpMatrix ColSpMatrix
Definition PaStiXSupport.h:613
void grabMatrix(const MatrixType &matrix, ColSpMatrix &out)
Definition PaStiXSupport.h:667

References Eigen::PastixBase< Derived >::analyzePattern(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::grabMatrix().

+ Here is the call graph for this function:

◆ clean()

void Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::clean ( )
inlineprotectedinherited
231 {
232 eigen_assert(m_initisOk && "The Pastix structure should be allocated first");
233 m_iparm(IPARM_START_TASK) = API_TASK_CLEAN;
234 m_iparm(IPARM_END_TASK) = API_TASK_CLEAN;
235 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, (Scalar*)0,
236 m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
237 }

◆ cols()

Index Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::cols ( ) const
inlineinherited
201{ return m_size; }

◆ compute() [1/2]

void Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::compute ( ColSpMatrix mat)
protectedinherited
296{
297 eigen_assert(mat.rows() == mat.cols() && "The input matrix should be squared");
298
299 analyzePattern(mat);
300 factorize(mat);
301
302 m_iparm(IPARM_MATRIX_VERIFICATION) = API_NO;
303}
void factorize(ColSpMatrix &mat)
Definition PaStiXSupport.h:338

◆ compute() [2/2]

template<typename _MatrixType , int _UpLo>
void Eigen::PastixLDLT< _MatrixType, _UpLo >::compute ( const MatrixType matrix)
inline

Compute the L and D factors of the LDL^T factorization of matrix

See also
analyzePattern() factorize()
632 {
633 ColSpMatrix temp;
634 grabMatrix(matrix, temp);
635 Base::compute(temp);
636 }
void compute(ColSpMatrix &mat)
Definition PaStiXSupport.h:295

References Eigen::PastixBase< Derived >::compute(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::grabMatrix().

Referenced by Eigen::PastixLDLT< _MatrixType, _UpLo >::PastixLDLT().

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

◆ derived() [1/2]

PastixLDLT< _MatrixType, _UpLo > & Eigen::SparseSolverBase< PastixLDLT< _MatrixType, _UpLo > >::derived ( )
inlineprotectedinherited
79{ return *static_cast<Derived*>(this); }

◆ derived() [2/2]

const PastixLDLT< _MatrixType, _UpLo > & Eigen::SparseSolverBase< PastixLDLT< _MatrixType, _UpLo > >::derived ( ) const
inlineprotectedinherited
80{ return *static_cast<const Derived*>(this); }

◆ dparm() [1/2]

Array< double, DPARM_SIZE, 1 > & Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::dparm ( )
inlineinherited

Returns a reference to the double vector DPARM of PaStiX parameters The statistics related to the different phases of factorization and solve are saved here as well

See also
analyzePattern() factorize()
188 {
189 return m_dparm;
190 }

◆ dparm() [2/2]

double & Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::dparm ( int  idxparam)
inlineinherited

Return a reference to a particular index parameter of the DPARM vector

See also
dparm()
197 {
198 return m_dparm(idxparam);
199 }

◆ factorize() [1/2]

void Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::factorize ( ColSpMatrix mat)
protectedinherited
339{
340// if(&m_cpyMat != &mat) m_cpyMat = mat;
341 eigen_assert(m_analysisIsOk && "The analysis phase should be called before the factorization phase");
342 m_iparm(IPARM_START_TASK) = API_TASK_NUMFACT;
343 m_iparm(IPARM_END_TASK) = API_TASK_NUMFACT;
344 m_size = internal::convert_index<int>(mat.rows());
345
346 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, m_size, mat.outerIndexPtr(), mat.innerIndexPtr(),
347 mat.valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
348
349 // Check the returned error
350 if(m_iparm(IPARM_ERROR_NUMBER))
351 {
353 m_factorizationIsOk = false;
354 m_isInitialized = false;
355 }
356 else
357 {
358 m_info = Success;
359 m_factorizationIsOk = true;
360 m_isInitialized = true;
361 }
362}
int m_factorizationIsOk
Definition PaStiXSupport.h:243
bool m_isInitialized
Definition SparseSolverBase.h:119

◆ factorize() [2/2]

template<typename _MatrixType , int _UpLo>
void Eigen::PastixLDLT< _MatrixType, _UpLo >::factorize ( const MatrixType matrix)
inline

Compute the LDL^T supernodal numerical factorization of matrix

652 {
653 ColSpMatrix temp;
654 grabMatrix(matrix, temp);
655 Base::factorize(temp);
656 }

References Eigen::PastixBase< Derived >::factorize(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::grabMatrix().

+ Here is the call graph for this function:

◆ grabMatrix()

template<typename _MatrixType , int _UpLo>
void Eigen::PastixLDLT< _MatrixType, _UpLo >::grabMatrix ( const MatrixType matrix,
ColSpMatrix out 
)
inlineprotected
668 {
669 // Pastix supports only lower, column-major matrices
670 out.resize(matrix.rows(), matrix.cols());
671 out.template selfadjointView<Lower>() = matrix.template selfadjointView<UpLo>();
673 }
void c_to_fortran_numbering(MatrixType &mat)
Definition PaStiXSupport.h:97

References Eigen::internal::c_to_fortran_numbering(), and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::resize().

Referenced by Eigen::PastixLDLT< _MatrixType, _UpLo >::analyzePattern(), Eigen::PastixLDLT< _MatrixType, _UpLo >::compute(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::factorize().

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

◆ info()

ComputationInfo Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::info ( ) const
inlineinherited

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the PaStiX reports a problem InvalidInput if the input matrix is invalid
See also
iparm()
213 {
214 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
215 return m_info;
216 }

◆ init()

template<typename _MatrixType , int _UpLo>
void Eigen::PastixLDLT< _MatrixType, _UpLo >::init ( )
inlineprotected
662 {
663 m_iparm(IPARM_SYM) = API_SYM_YES;
664 m_iparm(IPARM_FACTORIZATION) = API_FACT_LDLT;
665 }
Array< int, IPARM_SIZE, 1 > m_iparm
Definition PaStiXSupport.h:247

References Eigen::PastixLDLT< _MatrixType, _UpLo >::m_iparm.

Referenced by Eigen::PastixLDLT< _MatrixType, _UpLo >::PastixLDLT(), and Eigen::PastixLDLT< _MatrixType, _UpLo >::PastixLDLT().

+ Here is the caller graph for this function:

◆ iparm() [1/2]

Array< StorageIndex, IPARM_SIZE, 1 > & Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::iparm ( )
inlineinherited

Returns a reference to the integer vector IPARM of PaStiX parameters to modify the default parameters. The statistics related to the different phases of factorization and solve are saved here as well

See also
analyzePattern() factorize()
170 {
171 return m_iparm;
172 }

◆ iparm() [2/2]

int & Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::iparm ( int  idxparam)
inlineinherited

Return a reference to a particular index parameter of the IPARM vector

See also
iparm()
179 {
180 return m_iparm(idxparam);
181 }

◆ rows()

Index Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::rows ( ) const
inlineinherited
202{ return m_size; }

◆ 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
Derived & derived()
Definition SparseSolverBase.h:79
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

int Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_analysisIsOk
protectedinherited

◆ m_comm

int Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_comm
mutableprotectedinherited

◆ m_dparm

Array<double,DPARM_SIZE,1> Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_dparm
mutableprotectedinherited

◆ m_factorizationIsOk

int Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_factorizationIsOk
protectedinherited

◆ m_info

ComputationInfo Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_info
mutableprotectedinherited

◆ m_initisOk

int Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_initisOk
protectedinherited

◆ m_invp

Matrix<StorageIndex,Dynamic,1> Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_invp
mutableprotectedinherited

◆ m_iparm

template<typename _MatrixType , int _UpLo>
Array<int,IPARM_SIZE,1> Eigen::PastixBase< Derived >::m_iparm
mutableprotected

◆ m_isInitialized

bool Eigen::SparseSolverBase< PastixLDLT< _MatrixType, _UpLo > >::m_isInitialized
mutableprotectedinherited

◆ m_pastixdata

pastix_data_t* Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_pastixdata
mutableprotectedinherited

◆ m_perm

Matrix<StorageIndex,Dynamic,1> Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_perm
mutableprotectedinherited

◆ m_size

int Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >::m_size
mutableprotectedinherited

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