10 #ifndef EIGEN_PASTIXSUPPORT_H 11 #define EIGEN_PASTIXSUPPORT_H 23 template<
typename _MatrixType,
bool IsStrSym = false>
class PastixLU;
24 template<
typename _MatrixType,
int Options>
class PastixLLT;
25 template<
typename _MatrixType,
int Options>
class PastixLDLT;
32 template<
typename _MatrixType>
35 typedef _MatrixType MatrixType;
37 typedef typename _MatrixType::RealScalar RealScalar;
38 typedef typename _MatrixType::Index Index;
41 template<
typename _MatrixType,
int Options>
44 typedef _MatrixType MatrixType;
46 typedef typename _MatrixType::RealScalar RealScalar;
47 typedef typename _MatrixType::Index Index;
50 template<
typename _MatrixType,
int Options>
53 typedef _MatrixType MatrixType;
55 typedef typename _MatrixType::RealScalar RealScalar;
56 typedef typename _MatrixType::Index Index;
59 void eigen_pastix(pastix_data_t **pastix_data,
int pastix_comm,
int n,
int *ptr,
int *idx,
float *vals,
int *perm,
int * invp,
float *x,
int nbrhs,
int *iparm,
double *dparm)
61 if (n == 0) { ptr = NULL; idx = NULL; vals = NULL; }
62 if (nbrhs == 0) {x = NULL; nbrhs=1;}
63 s_pastix(pastix_data, pastix_comm, n, ptr, idx, vals, perm, invp, x, nbrhs, iparm, dparm);
66 void eigen_pastix(pastix_data_t **pastix_data,
int pastix_comm,
int n,
int *ptr,
int *idx,
double *vals,
int *perm,
int * invp,
double *x,
int nbrhs,
int *iparm,
double *dparm)
68 if (n == 0) { ptr = NULL; idx = NULL; vals = NULL; }
69 if (nbrhs == 0) {x = NULL; nbrhs=1;}
70 d_pastix(pastix_data, pastix_comm, n, ptr, idx, vals, perm, invp, x, nbrhs, iparm, dparm);
73 void eigen_pastix(pastix_data_t **pastix_data,
int pastix_comm,
int n,
int *ptr,
int *idx, std::complex<float> *vals,
int *perm,
int * invp, std::complex<float> *x,
int nbrhs,
int *iparm,
double *dparm)
75 if (n == 0) { ptr = NULL; idx = NULL; vals = NULL; }
76 if (nbrhs == 0) {x = NULL; nbrhs=1;}
77 c_pastix(pastix_data, pastix_comm, n, ptr, idx, reinterpret_cast<COMPLEX*>(vals), perm, invp, reinterpret_cast<COMPLEX*>(x), nbrhs, iparm, dparm);
80 void eigen_pastix(pastix_data_t **pastix_data,
int pastix_comm,
int n,
int *ptr,
int *idx, std::complex<double> *vals,
int *perm,
int * invp, std::complex<double> *x,
int nbrhs,
int *iparm,
double *dparm)
82 if (n == 0) { ptr = NULL; idx = NULL; vals = NULL; }
83 if (nbrhs == 0) {x = NULL; nbrhs=1;}
84 z_pastix(pastix_data, pastix_comm, n, ptr, idx, reinterpret_cast<DCOMPLEX*>(vals), perm, invp, reinterpret_cast<DCOMPLEX*>(x), nbrhs, iparm, dparm);
88 template <
typename MatrixType>
89 void c_to_fortran_numbering (MatrixType& mat)
91 if ( !(mat.outerIndexPtr()[0]) )
94 for(i = 0; i <= mat.rows(); ++i)
95 ++mat.outerIndexPtr()[i];
96 for(i = 0; i < mat.nonZeros(); ++i)
97 ++mat.innerIndexPtr()[i];
102 template <
typename MatrixType>
103 void fortran_to_c_numbering (MatrixType& mat)
106 if ( mat.outerIndexPtr()[0] == 1 )
109 for(i = 0; i <= mat.rows(); ++i)
110 --mat.outerIndexPtr()[i];
111 for(i = 0; i < mat.nonZeros(); ++i)
112 --mat.innerIndexPtr()[i];
119 template <
class Derived>
124 typedef _MatrixType MatrixType;
126 typedef typename MatrixType::RealScalar RealScalar;
127 typedef typename MatrixType::Index Index;
133 PastixBase() : m_initisOk(
false), m_analysisIsOk(
false), m_factorizationIsOk(
false), m_isInitialized(
false), m_pastixdata(0), m_size(0)
147 template<
typename Rhs>
151 eigen_assert(m_isInitialized &&
"Pastix solver is not initialized.");
152 eigen_assert(rows()==b.rows()
153 &&
"PastixBase::solve(): invalid number of rows of the right hand side matrix b");
157 template<
typename Rhs,
typename Dest>
162 return *
static_cast<Derived*
>(
this);
164 const Derived& derived()
const 166 return *
static_cast<const Derived*
>(
this);
185 return m_iparm(idxparam);
203 return m_dparm(idxparam);
206 inline Index cols()
const {
return m_size; }
207 inline Index rows()
const {
return m_size; }
219 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
227 template<
typename Rhs>
231 eigen_assert(m_isInitialized &&
"Pastix LU, LLT or LDLT is not initialized.");
232 eigen_assert(rows()==b.
rows()
233 &&
"PastixBase::solve(): invalid number of rows of the right hand side matrix b");
243 void analyzePattern(ColSpMatrix& mat);
246 void factorize(ColSpMatrix& mat);
251 eigen_assert(m_initisOk &&
"The Pastix structure should be allocated first");
252 m_iparm(IPARM_START_TASK) = API_TASK_CLEAN;
253 m_iparm(IPARM_END_TASK) = API_TASK_CLEAN;
254 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, (Scalar*)0,
255 m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
258 void compute(ColSpMatrix& mat);
262 int m_factorizationIsOk;
263 bool m_isInitialized;
265 mutable pastix_data_t *m_pastixdata;
278 template <
class Derived>
282 m_iparm.setZero(IPARM_SIZE);
283 m_dparm.setZero(DPARM_SIZE);
285 m_iparm(IPARM_MODIFY_PARAMETER) = API_NO;
286 pastix(&m_pastixdata, MPI_COMM_WORLD,
288 0, 0, 0, 1, m_iparm.data(), m_dparm.data());
290 m_iparm[IPARM_MATRIX_VERIFICATION] = API_NO;
291 m_iparm[IPARM_VERBOSE] = 2;
292 m_iparm[IPARM_ORDERING] = API_ORDER_SCOTCH;
293 m_iparm[IPARM_INCOMPLETE] = API_NO;
294 m_iparm[IPARM_OOC_LIMIT] = 2000;
295 m_iparm[IPARM_RHS_MAKING] = API_RHS_B;
296 m_iparm(IPARM_MATRIX_VERIFICATION) = API_NO;
298 m_iparm(IPARM_START_TASK) = API_TASK_INIT;
299 m_iparm(IPARM_END_TASK) = API_TASK_INIT;
300 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, (Scalar*)0,
301 0, 0, 0, 0, m_iparm.data(), m_dparm.data());
304 if(m_iparm(IPARM_ERROR_NUMBER)) {
314 template <
class Derived>
317 eigen_assert(mat.
rows() == mat.
cols() &&
"The input matrix should be squared");
322 m_iparm(IPARM_MATRIX_VERIFICATION) = API_NO;
323 m_isInitialized = m_factorizationIsOk;
327 template <
class Derived>
330 eigen_assert(m_initisOk &&
"The initialization of PaSTiX failed");
337 m_perm.resize(m_size);
338 m_invp.resize(m_size);
340 m_iparm(IPARM_START_TASK) = API_TASK_ORDERING;
341 m_iparm(IPARM_END_TASK) = API_TASK_ANALYSE;
343 mat.
valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
346 if(m_iparm(IPARM_ERROR_NUMBER))
349 m_analysisIsOk =
false;
354 m_analysisIsOk =
true;
358 template <
class Derived>
362 eigen_assert(m_analysisIsOk &&
"The analysis phase should be called before the factorization phase");
363 m_iparm(IPARM_START_TASK) = API_TASK_NUMFACT;
364 m_iparm(IPARM_END_TASK) = API_TASK_NUMFACT;
368 mat.
valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
371 if(m_iparm(IPARM_ERROR_NUMBER))
374 m_factorizationIsOk =
false;
375 m_isInitialized =
false;
380 m_factorizationIsOk =
true;
381 m_isInitialized =
true;
386 template<
typename Base>
387 template<
typename Rhs,
typename Dest>
390 eigen_assert(m_isInitialized &&
"The matrix should be factorized first");
392 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
397 for (
int i = 0; i < b.cols(); i++){
398 m_iparm[IPARM_START_TASK] = API_TASK_SOLVE;
399 m_iparm[IPARM_END_TASK] = API_TASK_REFINE;
401 internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, x.rows(), 0, 0, 0,
402 m_perm.data(), m_invp.data(), &x(0, i), rhs, m_iparm.data(), m_dparm.data());
408 return m_iparm(IPARM_ERROR_NUMBER)==0;
430 template<
typename _MatrixType,
bool IsStrSym>
434 typedef _MatrixType MatrixType;
437 typedef typename MatrixType::Index Index;
445 PastixLU(
const MatrixType& matrix):Base()
457 m_structureIsUptodate =
false;
459 grabMatrix(matrix, temp);
469 m_structureIsUptodate =
false;
471 grabMatrix(matrix, temp);
472 Base::analyzePattern(temp);
483 grabMatrix(matrix, temp);
484 Base::factorize(temp);
490 m_structureIsUptodate =
false;
491 m_iparm(IPARM_SYM) = API_SYM_NO;
492 m_iparm(IPARM_FACTORIZATION) = API_FACT_LU;
495 void grabMatrix(
const MatrixType& matrix, ColSpMatrix& out)
501 if(!m_structureIsUptodate)
504 m_transposedStructure = matrix.transpose();
507 for (Index j=0; j<m_transposedStructure.outerSize(); ++j)
508 for(
typename ColSpMatrix::InnerIterator it(m_transposedStructure, j); it; ++it)
511 m_structureIsUptodate =
true;
514 out = m_transposedStructure + matrix;
516 internal::c_to_fortran_numbering(out);
522 ColSpMatrix m_transposedStructure;
523 bool m_structureIsUptodate;
540 template<
typename _MatrixType,
int _UpLo>
544 typedef _MatrixType MatrixType;
549 enum { UpLo = _UpLo };
555 PastixLLT(
const MatrixType& matrix):Base()
567 grabMatrix(matrix, temp);
578 grabMatrix(matrix, temp);
579 Base::analyzePattern(temp);
587 grabMatrix(matrix, temp);
588 Base::factorize(temp);
595 m_iparm(IPARM_SYM) = API_SYM_YES;
596 m_iparm(IPARM_FACTORIZATION) = API_FACT_LLT;
599 void grabMatrix(
const MatrixType& matrix, ColSpMatrix& out)
602 out.template selfadjointView<Lower>() = matrix.template selfadjointView<UpLo>();
603 internal::c_to_fortran_numbering(out);
621 template<
typename _MatrixType,
int _UpLo>
625 typedef _MatrixType MatrixType;
630 enum { UpLo = _UpLo };
648 grabMatrix(matrix, temp);
659 grabMatrix(matrix, temp);
660 Base::analyzePattern(temp);
668 grabMatrix(matrix, temp);
669 Base::factorize(temp);
677 m_iparm(IPARM_SYM) = API_SYM_YES;
678 m_iparm(IPARM_FACTORIZATION) = API_FACT_LDLT;
681 void grabMatrix(
const MatrixType& matrix, ColSpMatrix& out)
684 out.template selfadjointView<Lower>() = matrix.template selfadjointView<UpLo>();
685 internal::c_to_fortran_numbering(out);
691 template<
typename _MatrixType,
typename Rhs>
696 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
698 template<
typename Dest>
void evalTo(Dest& dst)
const 700 dec()._solve(rhs(),dst);
704 template<
typename _MatrixType,
typename Rhs>
709 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
711 template<
typename Dest>
void evalTo(Dest& dst)
const 713 this->defaultEvalTo(dst);
Eigen::Matrix< Scalar, n, 1 > Vector
A vector of length n.
Definition: FlexibleKalmanBase.h:98
Definition: gtest_unittest.cc:5031
Definition: ForwardDeclarations.h:124
const Scalar * valuePtr() const
Definition: SparseMatrix.h:131
void compute(const MatrixType &matrix)
Compute the LU supernodal factorization of matrix.
Definition: PaStiXSupport.h:455
void factorize(const MatrixType &matrix)
Compute the LU supernodal factorization of matrix WARNING The matrix matrix should have the same stru...
Definition: PaStiXSupport.h:480
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: TestIMU_Common.h:87
void factorize(const MatrixType &matrix)
Compute the LDL^T supernodal numerical factorization of matrix.
Definition: PaStiXSupport.h:665
Definition: SparseSolve.h:18
Index cols() const
Definition: SparseMatrix.h:121
void compute(const MatrixType &matrix)
Compute the L factor of the LL^T supernodal factorization of matrix.
Definition: PaStiXSupport.h:564
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition: Constants.h:53
void compute(const MatrixType &matrix)
Compute the L and D factors of the LDL^T factorization of matrix.
Definition: PaStiXSupport.h:645
const Index * outerIndexPtr() const
Definition: SparseMatrix.h:149
Array< RealScalar, IPARM_SIZE, 1 > & dparm()
Returns a reference to the double vector DPARM of PaStiX parameters The statistics related to the dif...
Definition: PaStiXSupport.h:192
void init()
Initialize the PaStiX data structure.
Definition: PaStiXSupport.h:279
double & dparm(int idxparam)
Return a reference to a particular index parameter of the DPARM vector.
Definition: PaStiXSupport.h:201
void analyzePattern(const MatrixType &matrix)
Compute the LL^T symbolic factorization of matrix using its sparsity pattern The result of this opera...
Definition: PaStiXSupport.h:575
Index rows() const
Definition: SparseMatrixBase.h:159
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:239
A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library...
Definition: PaStiXSupport.h:25
Definition: PaStiXSupport.h:120
The provided data did not satisfy the prerequisites.
Definition: Constants.h:378
Index rows() const
Definition: SparseMatrix.h:119
int & iparm(int idxparam)
Return a reference to a particular index parameter of the IPARM vector.
Definition: PaStiXSupport.h:183
Definition: PaStiXSupport.h:30
The inputs are invalid, or the algorithm has been improperly called.
Definition: Constants.h:383
Computation was successful.
Definition: Constants.h:376
void analyzePattern(const MatrixType &matrix)
Compute the LU symbolic factorization of matrix using its sparsity pattern.
Definition: PaStiXSupport.h:467
Definition: BandTriangularSolver.h:13
Definition: SparseSolve.h:17
Array< Index, IPARM_SIZE, 1 > & iparm()
Returns a reference to the integer vector IPARM of PaStiX parameters to modify the default parameters...
Definition: PaStiXSupport.h:174
const Index * innerIndexPtr() const
Definition: SparseMatrix.h:140
void factorize(const MatrixType &matrix)
Compute the LL^T supernodal numerical factorization of matrix.
Definition: PaStiXSupport.h:584
const internal::solve_retval< PastixBase, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: PaStiXSupport.h:149
General-purpose arrays with easy API for coefficient-wise operations.
Definition: Array.h:42
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: PaStiXSupport.h:217
Definition: ForwardDeclarations.h:125
Interface to the PaStix solver.
Definition: PaStiXSupport.h:23
A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library...
Definition: PaStiXSupport.h:24
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:374
void analyzePattern(const MatrixType &matrix)
Compute the LDL^T symbolic factorization of matrix using its sparsity pattern The result of this oper...
Definition: PaStiXSupport.h:656
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const internal::sparse_solve_retval< PastixBase, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
Definition: PaStiXSupport.h:229
Definition: osvr_print_tree.cpp:52
double Scalar
Common scalar type.
Definition: FlexibleKalmanBase.h:48