32#ifndef EIGEN_PARDISOSUPPORT_H
33#define EIGEN_PARDISOSUPPORT_H
36#include "./InternalHeaderCheck.h"
40template <
typename MatrixType_>
42template <
typename MatrixType_,
int Options = Upper>
44template <
typename MatrixType_,
int Options = Upper>
48template <
typename IndexType>
49struct pardiso_run_selector {
50 static IndexType run(_MKL_DSS_HANDLE_t pt, IndexType maxfct, IndexType mnum, IndexType type, IndexType phase,
51 IndexType n,
void* a, IndexType* ia, IndexType* ja, IndexType* perm, IndexType nrhs,
52 IndexType* iparm, IndexType msglvl,
void* b,
void* x) {
54 ::pardiso(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
59struct pardiso_run_selector<long long int> {
60 typedef long long int IndexType;
61 static IndexType run(_MKL_DSS_HANDLE_t pt, IndexType maxfct, IndexType mnum, IndexType type, IndexType phase,
62 IndexType n,
void* a, IndexType* ia, IndexType* ja, IndexType* perm, IndexType nrhs,
63 IndexType* iparm, IndexType msglvl,
void* b,
void* x) {
65 ::pardiso_64(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
70template <
class Pardiso>
73template <
typename MatrixType_>
74struct pardiso_traits<PardisoLU<MatrixType_> > {
75 typedef MatrixType_ MatrixType;
76 typedef typename MatrixType_::Scalar Scalar;
77 typedef typename MatrixType_::RealScalar RealScalar;
78 typedef typename MatrixType_::StorageIndex StorageIndex;
81template <
typename MatrixType_,
int Options>
82struct pardiso_traits<PardisoLLT<MatrixType_, Options> > {
83 typedef MatrixType_ MatrixType;
84 typedef typename MatrixType_::Scalar Scalar;
85 typedef typename MatrixType_::RealScalar RealScalar;
86 typedef typename MatrixType_::StorageIndex StorageIndex;
89template <
typename MatrixType_,
int Options>
90struct pardiso_traits<PardisoLDLT<MatrixType_, Options> > {
91 typedef MatrixType_ MatrixType;
92 typedef typename MatrixType_::Scalar Scalar;
93 typedef typename MatrixType_::RealScalar RealScalar;
94 typedef typename MatrixType_::StorageIndex StorageIndex;
99template <
class Derived>
100class PardisoImpl :
public SparseSolverBase<Derived> {
104 using Base::m_isInitialized;
106 typedef internal::pardiso_traits<Derived> Traits;
109 using Base::_solve_impl;
111 typedef typename Traits::MatrixType MatrixType;
112 typedef typename Traits::Scalar Scalar;
113 typedef typename Traits::RealScalar RealScalar;
114 typedef typename Traits::StorageIndex StorageIndex;
115 typedef SparseMatrix<Scalar, RowMajor, StorageIndex> SparseMatrixType;
116 typedef Matrix<Scalar, Dynamic, 1> VectorType;
117 typedef Matrix<StorageIndex, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
118 typedef Matrix<StorageIndex, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
119 typedef Array<StorageIndex, 64, 1, DontAlign> ParameterType;
120 enum { ScalarIsComplex = NumTraits<Scalar>::IsComplex, ColsAtCompileTime =
Dynamic, MaxColsAtCompileTime =
Dynamic };
122 PardisoImpl() : m_analysisIsOk(false), m_factorizationIsOk(false) {
123 eigen_assert((
sizeof(StorageIndex) >=
sizeof(_INTEGER_t) &&
sizeof(StorageIndex) <= 8) &&
124 "Non-supported index type");
127 m_isInitialized =
false;
130 ~PardisoImpl() { pardisoRelease(); }
132 inline Index cols()
const {
return m_size; }
133 inline Index rows()
const {
return m_size; }
140 ComputationInfo info()
const {
141 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
148 ParameterType& pardisoParameterArray() {
return m_iparm; }
156 Derived& analyzePattern(
const MatrixType& matrix);
164 Derived& factorize(
const MatrixType& matrix);
166 Derived& compute(
const MatrixType& matrix);
168 template <
typename Rhs,
typename Dest>
169 void _solve_impl(
const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest)
const;
172 void pardisoRelease() {
175 internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, -1,
176 internal::convert_index<StorageIndex>(m_size), 0, 0, 0,
177 m_perm.
data(), 0, m_iparm.
data(), m_msglvl, NULL, NULL);
178 m_isInitialized =
false;
182 void pardisoInit(
int type) {
184 bool symmetric = std::abs(m_type) < 10;
195 m_iparm[10] = symmetric ? 0 : 1;
197 m_iparm[12] = symmetric ? 0 : 1;
209 m_iparm[27] = (
sizeof(RealScalar) == 4) ? 1 : 0;
214 memset(m_pt, 0,
sizeof(m_pt));
220 void manageErrorCode(Index error)
const {
234 mutable SparseMatrixType m_matrix;
235 mutable ComputationInfo m_info;
236 bool m_analysisIsOk, m_factorizationIsOk;
237 StorageIndex m_type, m_msglvl;
238 mutable void* m_pt[64];
239 mutable ParameterType m_iparm;
240 mutable IntColVectorType m_perm;
244template <
class Derived>
245Derived& PardisoImpl<Derived>::compute(
const MatrixType& a) {
247 eigen_assert(a.rows() == a.cols());
250 m_perm.setZero(m_size);
251 derived().getMatrix(a);
254 error = internal::pardiso_run_selector<StorageIndex>::run(
255 m_pt, 1, 1, m_type, 12, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
256 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
257 manageErrorCode(error);
260 m_isInitialized =
true;
264template <
class Derived>
265Derived& PardisoImpl<Derived>::analyzePattern(
const MatrixType& a) {
267 eigen_assert(m_size == a.cols());
270 m_perm.setZero(m_size);
271 derived().getMatrix(a);
274 error = internal::pardiso_run_selector<StorageIndex>::run(
275 m_pt, 1, 1, m_type, 11, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
276 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
278 manageErrorCode(error);
280 m_factorizationIsOk =
false;
281 m_isInitialized =
true;
285template <
class Derived>
286Derived& PardisoImpl<Derived>::factorize(
const MatrixType& a) {
287 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
288 eigen_assert(m_size == a.rows() && m_size == a.cols());
290 derived().getMatrix(a);
293 error = internal::pardiso_run_selector<StorageIndex>::run(
294 m_pt, 1, 1, m_type, 22, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
295 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
297 manageErrorCode(error);
302template <
class Derived>
303template <
typename BDerived,
typename XDerived>
304void PardisoImpl<Derived>::_solve_impl(
const MatrixBase<BDerived>& b, MatrixBase<XDerived>& x)
const {
312 Index nrhs = Index(b.cols());
313 eigen_assert(m_size == b.rows());
314 eigen_assert(((MatrixBase<BDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) &&
315 "Row-major right hand sides are not supported");
316 eigen_assert(((MatrixBase<XDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) &&
317 "Row-major matrices of unknowns are not supported");
318 eigen_assert(((nrhs == 1) || b.outerStride() == b.rows()));
329 Scalar* rhs_ptr =
const_cast<Scalar*
>(b.derived().data());
330 Matrix<Scalar, Dynamic, Dynamic, ColMajor> tmp;
333 if (rhs_ptr == x.derived().data()) {
335 rhs_ptr = tmp.data();
339 error = internal::pardiso_run_selector<StorageIndex>::run(
340 m_pt, 1, 1, m_type, 33, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
341 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), internal::convert_index<StorageIndex>(nrhs),
342 m_iparm.data(), m_msglvl, rhs_ptr, x.derived().data());
344 manageErrorCode(error);
364template <
typename MatrixType>
365class PardisoLU :
public PardisoImpl<PardisoLU<MatrixType> > {
367 typedef PardisoImpl<PardisoLU> Base;
368 using Base::m_matrix;
369 using Base::pardisoInit;
370 friend class PardisoImpl<
PardisoLU<MatrixType> >;
373 typedef typename Base::Scalar Scalar;
374 typedef typename Base::RealScalar RealScalar;
379 PardisoLU() : Base() { pardisoInit(Base::ScalarIsComplex ? 13 : 11); }
381 explicit PardisoLU(
const MatrixType& matrix) : Base() {
382 pardisoInit(Base::ScalarIsComplex ? 13 : 11);
387 void getMatrix(
const MatrixType& matrix) {
412template <
typename MatrixType,
int UpLo_>
413class PardisoLLT :
public PardisoImpl<PardisoLLT<MatrixType, UpLo_> > {
415 typedef PardisoImpl<PardisoLLT<MatrixType, UpLo_> > Base;
416 using Base::m_matrix;
417 using Base::pardisoInit;
418 friend class PardisoImpl<
PardisoLLT<MatrixType, UpLo_> >;
421 typedef typename Base::Scalar Scalar;
422 typedef typename Base::RealScalar RealScalar;
423 typedef typename Base::StorageIndex StorageIndex;
424 enum { UpLo = UpLo_ };
427 PardisoLLT() : Base() { pardisoInit(Base::ScalarIsComplex ? 4 : 2); }
429 explicit PardisoLLT(
const MatrixType& matrix) : Base() {
430 pardisoInit(Base::ScalarIsComplex ? 4 : 2);
435 void getMatrix(
const MatrixType& matrix) {
438 m_matrix.
resize(matrix.rows(), matrix.cols());
439 m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().
twistedBy(p_null);
466template <
typename MatrixType,
int Options>
467class PardisoLDLT :
public PardisoImpl<PardisoLDLT<MatrixType, Options> > {
469 typedef PardisoImpl<PardisoLDLT<MatrixType, Options> > Base;
470 using Base::m_matrix;
471 using Base::pardisoInit;
472 friend class PardisoImpl<
PardisoLDLT<MatrixType, Options> >;
475 typedef typename Base::Scalar Scalar;
476 typedef typename Base::RealScalar RealScalar;
477 typedef typename Base::StorageIndex StorageIndex;
481 PardisoLDLT() : Base() { pardisoInit(Base::ScalarIsComplex ? (
bool(Options &
Symmetric) ? 6 : -4) : -2); }
483 explicit PardisoLDLT(
const MatrixType& matrix) : Base() {
484 pardisoInit(Base::ScalarIsComplex ? (
bool(Options &
Symmetric) ? 6 : -4) : -2);
488 void getMatrix(
const MatrixType& matrix) {
491 m_matrix.
resize(matrix.rows(), matrix.cols());
492 m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().
twistedBy(p_null);
A sparse direct Cholesky (LDLT) factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:467
A sparse direct Cholesky (LLT) factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:413
A sparse direct LU factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:365
Permutation matrix.
Definition PermutationMatrix.h:280
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:563
const Scalar * data() const
Definition PlainObjectBase.h:273
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition SparseMatrixBase.h:325
void makeCompressed()
Definition SparseMatrix.h:588
void resize(Index rows, Index cols)
Definition SparseMatrix.h:730
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:84
SparseSolverBase()
Definition SparseSolverBase.h:70
@ Symmetric
Definition Constants.h:229
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
@ NumericalIssue
Definition Constants.h:442
@ InvalidInput
Definition Constants.h:447
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition Core:137
const int Dynamic
Definition Constants.h:25