Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
PardisoSupport.h
1/*
2 Copyright (c) 2011, Intel Corporation. All rights reserved.
3
4 Redistribution and use in source and binary forms, with or without modification,
5 are permitted provided that the following conditions are met:
6
7 * Redistributions of source code must retain the above copyright notice, this
8 list of conditions and the following disclaimer.
9 * Redistributions in binary form must reproduce the above copyright notice,
10 this list of conditions and the following disclaimer in the documentation
11 and/or other materials provided with the distribution.
12 * Neither the name of Intel Corporation nor the names of its contributors may
13 be used to endorse or promote products derived from this software without
14 specific prior written permission.
15
16 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
17 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
18 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
19 DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
20 ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
21 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
22 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
23 ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
24 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
25 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26
27 ********************************************************************************
28 * Content : Eigen bindings to Intel(R) MKL PARDISO
29 ********************************************************************************
30*/
31
32#ifndef EIGEN_PARDISOSUPPORT_H
33#define EIGEN_PARDISOSUPPORT_H
34
35// IWYU pragma: private
36#include "./InternalHeaderCheck.h"
37
38namespace Eigen {
39
40template <typename MatrixType_>
41class PardisoLU;
42template <typename MatrixType_, int Options = Upper>
43class PardisoLLT;
44template <typename MatrixType_, int Options = Upper>
45class PardisoLDLT;
46
47namespace internal {
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) {
53 IndexType error = 0;
54 ::pardiso(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
55 return error;
56 }
57};
58template <>
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) {
64 IndexType error = 0;
65 ::pardiso_64(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
66 return error;
67 }
68};
69
70template <class Pardiso>
71struct pardiso_traits;
72
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;
79};
80
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;
87};
88
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;
95};
96
97} // end namespace internal
98
99template <class Derived>
100class PardisoImpl : public SparseSolverBase<Derived> {
101 protected:
102 typedef SparseSolverBase<Derived> Base;
103 using Base::derived;
104 using Base::m_isInitialized;
105
106 typedef internal::pardiso_traits<Derived> Traits;
107
108 public:
109 using Base::_solve_impl;
110
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 };
121
122 PardisoImpl() : m_analysisIsOk(false), m_factorizationIsOk(false) {
123 eigen_assert((sizeof(StorageIndex) >= sizeof(_INTEGER_t) && sizeof(StorageIndex) <= 8) &&
124 "Non-supported index type");
125 m_iparm.setZero();
126 m_msglvl = 0; // No output
127 m_isInitialized = false;
128 }
129
130 ~PardisoImpl() { pardisoRelease(); }
131
132 inline Index cols() const { return m_size; }
133 inline Index rows() const { return m_size; }
134
140 ComputationInfo info() const {
141 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
142 return m_info;
143 }
144
148 ParameterType& pardisoParameterArray() { return m_iparm; }
149
156 Derived& analyzePattern(const MatrixType& matrix);
157
165 Derived& factorize(const MatrixType& matrix);
166
167 Derived& compute(const MatrixType& matrix);
168
169 template <typename Rhs, typename Dest>
170 void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const;
171
172 protected:
173 void pardisoRelease() {
174 if (m_isInitialized) // Factorization ran at least once
175 {
176 internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, -1,
177 internal::convert_index<StorageIndex>(m_size), 0, 0, 0,
178 m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
179 m_isInitialized = false;
180 }
181 }
182
183 void pardisoInit(int type) {
184 m_type = type;
185 bool symmetric = std::abs(m_type) < 10;
186 m_iparm[0] = 1; // No solver default
187 m_iparm[1] = 2; // use Metis for the ordering
188 m_iparm[2] = 0; // Reserved. Set to zero. (??Numbers of processors, value of OMP_NUM_THREADS??)
189 m_iparm[3] = 0; // No iterative-direct algorithm
190 m_iparm[4] = 0; // No user fill-in reducing permutation
191 m_iparm[5] = 0; // Write solution into x, b is left unchanged
192 m_iparm[6] = 0; // Not in use
193 m_iparm[7] = 2; // Max numbers of iterative refinement steps
194 m_iparm[8] = 0; // Not in use
195 m_iparm[9] = 13; // Perturb the pivot elements with 1E-13
196 m_iparm[10] = symmetric ? 0 : 1; // Use nonsymmetric permutation and scaling MPS
197 m_iparm[11] = 0; // Not in use
198 m_iparm[12] = symmetric ? 0 : 1; // Maximum weighted matching algorithm is switched-off (default for symmetric).
199 // Try m_iparm[12] = 1 in case of inappropriate accuracy
200 m_iparm[13] = 0; // Output: Number of perturbed pivots
201 m_iparm[14] = 0; // Not in use
202 m_iparm[15] = 0; // Not in use
203 m_iparm[16] = 0; // Not in use
204 m_iparm[17] = -1; // Output: Number of nonzeros in the factor LU
205 m_iparm[18] = -1; // Output: Mflops for LU factorization
206 m_iparm[19] = 0; // Output: Numbers of CG Iterations
207
208 m_iparm[20] = 0; // 1x1 pivoting
209 m_iparm[26] = 0; // No matrix checker
210 m_iparm[27] = (sizeof(RealScalar) == 4) ? 1 : 0;
211 m_iparm[34] = 1; // C indexing
212 m_iparm[36] = 0; // CSR
213 m_iparm[59] = 0; // 0 - In-Core ; 1 - Automatic switch between In-Core and Out-of-Core modes ; 2 - Out-of-Core
214
215 memset(m_pt, 0, sizeof(m_pt));
216 }
217
218 protected:
219 // cached data to reduce reallocation, etc.
220
221 void manageErrorCode(Index error) const {
222 switch (error) {
223 case 0:
224 m_info = Success;
225 break;
226 case -4:
227 case -7:
228 m_info = NumericalIssue;
229 break;
230 default:
231 m_info = InvalidInput;
232 }
233 }
234
235 mutable SparseMatrixType m_matrix;
236 mutable ComputationInfo m_info;
237 bool m_analysisIsOk, m_factorizationIsOk;
238 StorageIndex m_type, m_msglvl;
239 mutable void* m_pt[64];
240 mutable ParameterType m_iparm;
241 mutable IntColVectorType m_perm;
242 Index m_size;
243};
244
245template <class Derived>
246Derived& PardisoImpl<Derived>::compute(const MatrixType& a) {
247 m_size = a.rows();
248 eigen_assert(a.rows() == a.cols());
249
250 pardisoRelease();
251 m_perm.setZero(m_size);
252 derived().getMatrix(a);
253
254 Index error;
255 error = internal::pardiso_run_selector<StorageIndex>::run(
256 m_pt, 1, 1, m_type, 12, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
257 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
258 manageErrorCode(error);
259 m_analysisIsOk = m_info == Eigen::Success;
260 m_factorizationIsOk = m_info == Eigen::Success;
261 m_isInitialized = true;
262 return derived();
263}
264
265template <class Derived>
266Derived& PardisoImpl<Derived>::analyzePattern(const MatrixType& a) {
267 m_size = a.rows();
268 eigen_assert(m_size == a.cols());
269
270 pardisoRelease();
271 m_perm.setZero(m_size);
272 derived().getMatrix(a);
273
274 Index error;
275 error = internal::pardiso_run_selector<StorageIndex>::run(
276 m_pt, 1, 1, m_type, 11, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
277 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
278
279 manageErrorCode(error);
280 m_analysisIsOk = m_info == Eigen::Success;
281 m_factorizationIsOk = false;
282 m_isInitialized = true;
283 return derived();
284}
285
286template <class Derived>
287Derived& PardisoImpl<Derived>::factorize(const MatrixType& a) {
288 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
289 eigen_assert(m_size == a.rows() && m_size == a.cols());
290
291 derived().getMatrix(a);
292
293 Index error;
294 error = internal::pardiso_run_selector<StorageIndex>::run(
295 m_pt, 1, 1, m_type, 22, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
296 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
297
298 manageErrorCode(error);
299 m_factorizationIsOk = m_info == Eigen::Success;
300 return derived();
301}
302
303template <class Derived>
304template <typename BDerived, typename XDerived>
305void PardisoImpl<Derived>::_solve_impl(const MatrixBase<BDerived>& b, MatrixBase<XDerived>& x) const {
306 if (m_iparm[0] == 0) // Factorization was not computed
307 {
308 m_info = InvalidInput;
309 return;
310 }
311
312 // Index n = m_matrix.rows();
313 Index nrhs = Index(b.cols());
314 eigen_assert(m_size == b.rows());
315 eigen_assert(((MatrixBase<BDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) &&
316 "Row-major right hand sides are not supported");
317 eigen_assert(((MatrixBase<XDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) &&
318 "Row-major matrices of unknowns are not supported");
319 eigen_assert(((nrhs == 1) || b.outerStride() == b.rows()));
320
321 // switch (transposed) {
322 // case SvNoTrans : m_iparm[11] = 0 ; break;
323 // case SvTranspose : m_iparm[11] = 2 ; break;
324 // case SvAdjoint : m_iparm[11] = 1 ; break;
325 // default:
326 // //std::cerr << "Eigen: transposition option \"" << transposed << "\" not supported by the PARDISO backend\n";
327 // m_iparm[11] = 0;
328 // }
329
330 Scalar* rhs_ptr = const_cast<Scalar*>(b.derived().data());
332
333 // Pardiso cannot solve in-place
334 if (rhs_ptr == x.derived().data()) {
335 tmp = b;
336 rhs_ptr = tmp.data();
337 }
338
339 Index error;
340 error = internal::pardiso_run_selector<StorageIndex>::run(
341 m_pt, 1, 1, m_type, 33, internal::convert_index<StorageIndex>(m_size), m_matrix.valuePtr(),
342 m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(), m_perm.data(), internal::convert_index<StorageIndex>(nrhs),
343 m_iparm.data(), m_msglvl, rhs_ptr, x.derived().data());
344
345 manageErrorCode(error);
346}
347
365template <typename MatrixType>
366class PardisoLU : public PardisoImpl<PardisoLU<MatrixType> > {
367 protected:
368 typedef PardisoImpl<PardisoLU> Base;
369 using Base::m_matrix;
370 using Base::pardisoInit;
371 friend class PardisoImpl<PardisoLU<MatrixType> >;
372
373 public:
374 typedef typename Base::Scalar Scalar;
375 typedef typename Base::RealScalar RealScalar;
376
377 using Base::compute;
378 using Base::solve;
379
380 PardisoLU() : Base() { pardisoInit(Base::ScalarIsComplex ? 13 : 11); }
381
382 explicit PardisoLU(const MatrixType& matrix) : Base() {
383 pardisoInit(Base::ScalarIsComplex ? 13 : 11);
384 compute(matrix);
385 }
386
387 protected:
388 void getMatrix(const MatrixType& matrix) {
389 m_matrix = matrix;
390 m_matrix.makeCompressed();
391 }
392};
393
413template <typename MatrixType, int UpLo_>
414class PardisoLLT : public PardisoImpl<PardisoLLT<MatrixType, UpLo_> > {
415 protected:
416 typedef PardisoImpl<PardisoLLT<MatrixType, UpLo_> > Base;
417 using Base::m_matrix;
418 using Base::pardisoInit;
419 friend class PardisoImpl<PardisoLLT<MatrixType, UpLo_> >;
420
421 public:
422 typedef typename Base::Scalar Scalar;
423 typedef typename Base::RealScalar RealScalar;
424 typedef typename Base::StorageIndex StorageIndex;
425 enum { UpLo = UpLo_ };
426 using Base::compute;
427
428 PardisoLLT() : Base() { pardisoInit(Base::ScalarIsComplex ? 4 : 2); }
429
430 explicit PardisoLLT(const MatrixType& matrix) : Base() {
431 pardisoInit(Base::ScalarIsComplex ? 4 : 2);
432 compute(matrix);
433 }
434
435 protected:
436 void getMatrix(const MatrixType& matrix) {
437 // PARDISO supports only upper, row-major matrices
439 m_matrix.resize(matrix.rows(), matrix.cols());
440 m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().twistedBy(p_null);
441 m_matrix.makeCompressed();
442 }
443};
444
467template <typename MatrixType, int Options>
468class PardisoLDLT : public PardisoImpl<PardisoLDLT<MatrixType, Options> > {
469 protected:
470 typedef PardisoImpl<PardisoLDLT<MatrixType, Options> > Base;
471 using Base::m_matrix;
472 using Base::pardisoInit;
473 friend class PardisoImpl<PardisoLDLT<MatrixType, Options> >;
474
475 public:
476 typedef typename Base::Scalar Scalar;
477 typedef typename Base::RealScalar RealScalar;
478 typedef typename Base::StorageIndex StorageIndex;
479 using Base::compute;
480 enum { UpLo = Options & (Upper | Lower) };
481
482 PardisoLDLT() : Base() { pardisoInit(Base::ScalarIsComplex ? (bool(Options & Symmetric) ? 6 : -4) : -2); }
483
484 explicit PardisoLDLT(const MatrixType& matrix) : Base() {
485 pardisoInit(Base::ScalarIsComplex ? (bool(Options & Symmetric) ? 6 : -4) : -2);
486 compute(matrix);
487 }
488
489 void getMatrix(const MatrixType& matrix) {
490 // PARDISO supports only upper, row-major matrices
492 m_matrix.resize(matrix.rows(), matrix.cols());
493 m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().twistedBy(p_null);
494 m_matrix.makeCompressed();
495 }
496};
497
498} // end namespace Eigen
499
500#endif // EIGEN_PARDISOSUPPORT_H
@ Flags
Definition DenseBase.h:161
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
A sparse direct Cholesky (LDLT) factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:468
A sparse direct Cholesky (LLT) factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:414
A sparse direct LU factorization and solver based on the PARDISO library.
Definition PardisoSupport.h:366
void resize(Index newSize)
Definition PermutationMatrix.h:122
Permutation matrix.
Definition PermutationMatrix.h:283
A base class for sparse solvers.
Definition SparseSolverBase.h:67
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:84
SparseSolverBase()
Definition SparseSolverBase.h:70
ComputationInfo
Definition Constants.h:438
@ 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
const unsigned int RowMajorBit
Definition Constants.h:70
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const int Dynamic
Definition Constants.h:25