PardisoSupport.h 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545
  1. /*
  2. Copyright (c) 2011, Intel Corporation. All rights reserved.
  3. Redistribution and use in source and binary forms, with or without modification,
  4. are permitted provided that the following conditions are met:
  5. * Redistributions of source code must retain the above copyright notice, this
  6. list of conditions and the following disclaimer.
  7. * Redistributions in binary form must reproduce the above copyright notice,
  8. this list of conditions and the following disclaimer in the documentation
  9. and/or other materials provided with the distribution.
  10. * Neither the name of Intel Corporation nor the names of its contributors may
  11. be used to endorse or promote products derived from this software without
  12. specific prior written permission.
  13. THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
  14. ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
  15. WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
  16. DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
  17. ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
  18. (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
  19. LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
  20. ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
  21. (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
  22. SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
  23. ********************************************************************************
  24. * Content : Eigen bindings to Intel(R) MKL PARDISO
  25. ********************************************************************************
  26. */
  27. #ifndef EIGEN_PARDISOSUPPORT_H
  28. #define EIGEN_PARDISOSUPPORT_H
  29. namespace Eigen {
  30. template<typename _MatrixType> class PardisoLU;
  31. template<typename _MatrixType, int Options=Upper> class PardisoLLT;
  32. template<typename _MatrixType, int Options=Upper> class PardisoLDLT;
  33. namespace internal
  34. {
  35. template<typename IndexType>
  36. struct pardiso_run_selector
  37. {
  38. static IndexType run( _MKL_DSS_HANDLE_t pt, IndexType maxfct, IndexType mnum, IndexType type, IndexType phase, IndexType n, void *a,
  39. IndexType *ia, IndexType *ja, IndexType *perm, IndexType nrhs, IndexType *iparm, IndexType msglvl, void *b, void *x)
  40. {
  41. IndexType error = 0;
  42. ::pardiso(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
  43. return error;
  44. }
  45. };
  46. template<>
  47. struct pardiso_run_selector<long long int>
  48. {
  49. typedef long long int IndexType;
  50. static IndexType run( _MKL_DSS_HANDLE_t pt, IndexType maxfct, IndexType mnum, IndexType type, IndexType phase, IndexType n, void *a,
  51. IndexType *ia, IndexType *ja, IndexType *perm, IndexType nrhs, IndexType *iparm, IndexType msglvl, void *b, void *x)
  52. {
  53. IndexType error = 0;
  54. ::pardiso_64(pt, &maxfct, &mnum, &type, &phase, &n, a, ia, ja, perm, &nrhs, iparm, &msglvl, b, x, &error);
  55. return error;
  56. }
  57. };
  58. template<class Pardiso> struct pardiso_traits;
  59. template<typename _MatrixType>
  60. struct pardiso_traits< PardisoLU<_MatrixType> >
  61. {
  62. typedef _MatrixType MatrixType;
  63. typedef typename _MatrixType::Scalar Scalar;
  64. typedef typename _MatrixType::RealScalar RealScalar;
  65. typedef typename _MatrixType::StorageIndex StorageIndex;
  66. };
  67. template<typename _MatrixType, int Options>
  68. struct pardiso_traits< PardisoLLT<_MatrixType, Options> >
  69. {
  70. typedef _MatrixType MatrixType;
  71. typedef typename _MatrixType::Scalar Scalar;
  72. typedef typename _MatrixType::RealScalar RealScalar;
  73. typedef typename _MatrixType::StorageIndex StorageIndex;
  74. };
  75. template<typename _MatrixType, int Options>
  76. struct pardiso_traits< PardisoLDLT<_MatrixType, Options> >
  77. {
  78. typedef _MatrixType MatrixType;
  79. typedef typename _MatrixType::Scalar Scalar;
  80. typedef typename _MatrixType::RealScalar RealScalar;
  81. typedef typename _MatrixType::StorageIndex StorageIndex;
  82. };
  83. } // end namespace internal
  84. template<class Derived>
  85. class PardisoImpl : public SparseSolverBase<Derived>
  86. {
  87. protected:
  88. typedef SparseSolverBase<Derived> Base;
  89. using Base::derived;
  90. using Base::m_isInitialized;
  91. typedef internal::pardiso_traits<Derived> Traits;
  92. public:
  93. using Base::_solve_impl;
  94. typedef typename Traits::MatrixType MatrixType;
  95. typedef typename Traits::Scalar Scalar;
  96. typedef typename Traits::RealScalar RealScalar;
  97. typedef typename Traits::StorageIndex StorageIndex;
  98. typedef SparseMatrix<Scalar,RowMajor,StorageIndex> SparseMatrixType;
  99. typedef Matrix<Scalar,Dynamic,1> VectorType;
  100. typedef Matrix<StorageIndex, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
  101. typedef Matrix<StorageIndex, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
  102. typedef Array<StorageIndex,64,1,DontAlign> ParameterType;
  103. enum {
  104. ScalarIsComplex = NumTraits<Scalar>::IsComplex,
  105. ColsAtCompileTime = Dynamic,
  106. MaxColsAtCompileTime = Dynamic
  107. };
  108. PardisoImpl()
  109. : m_analysisIsOk(false), m_factorizationIsOk(false)
  110. {
  111. eigen_assert((sizeof(StorageIndex) >= sizeof(_INTEGER_t) && sizeof(StorageIndex) <= 8) && "Non-supported index type");
  112. m_iparm.setZero();
  113. m_msglvl = 0; // No output
  114. m_isInitialized = false;
  115. }
  116. ~PardisoImpl()
  117. {
  118. pardisoRelease();
  119. }
  120. inline Index cols() const { return m_size; }
  121. inline Index rows() const { return m_size; }
  122. /** \brief Reports whether previous computation was successful.
  123. *
  124. * \returns \c Success if computation was successful,
  125. * \c NumericalIssue if the matrix appears to be negative.
  126. */
  127. ComputationInfo info() const
  128. {
  129. eigen_assert(m_isInitialized && "Decomposition is not initialized.");
  130. return m_info;
  131. }
  132. /** \warning for advanced usage only.
  133. * \returns a reference to the parameter array controlling PARDISO.
  134. * See the PARDISO manual to know how to use it. */
  135. ParameterType& pardisoParameterArray()
  136. {
  137. return m_iparm;
  138. }
  139. /** Performs a symbolic decomposition on the sparcity of \a matrix.
  140. *
  141. * This function is particularly useful when solving for several problems having the same structure.
  142. *
  143. * \sa factorize()
  144. */
  145. Derived& analyzePattern(const MatrixType& matrix);
  146. /** Performs a numeric decomposition of \a matrix
  147. *
  148. * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
  149. *
  150. * \sa analyzePattern()
  151. */
  152. Derived& factorize(const MatrixType& matrix);
  153. Derived& compute(const MatrixType& matrix);
  154. template<typename Rhs,typename Dest>
  155. void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const;
  156. protected:
  157. void pardisoRelease()
  158. {
  159. if(m_isInitialized) // Factorization ran at least once
  160. {
  161. internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, -1, internal::convert_index<StorageIndex>(m_size),0, 0, 0, m_perm.data(), 0,
  162. m_iparm.data(), m_msglvl, NULL, NULL);
  163. m_isInitialized = false;
  164. }
  165. }
  166. void pardisoInit(int type)
  167. {
  168. m_type = type;
  169. bool symmetric = std::abs(m_type) < 10;
  170. m_iparm[0] = 1; // No solver default
  171. m_iparm[1] = 2; // use Metis for the ordering
  172. m_iparm[2] = 0; // Reserved. Set to zero. (??Numbers of processors, value of OMP_NUM_THREADS??)
  173. m_iparm[3] = 0; // No iterative-direct algorithm
  174. m_iparm[4] = 0; // No user fill-in reducing permutation
  175. m_iparm[5] = 0; // Write solution into x, b is left unchanged
  176. m_iparm[6] = 0; // Not in use
  177. m_iparm[7] = 2; // Max numbers of iterative refinement steps
  178. m_iparm[8] = 0; // Not in use
  179. m_iparm[9] = 13; // Perturb the pivot elements with 1E-13
  180. m_iparm[10] = symmetric ? 0 : 1; // Use nonsymmetric permutation and scaling MPS
  181. m_iparm[11] = 0; // Not in use
  182. m_iparm[12] = symmetric ? 0 : 1; // Maximum weighted matching algorithm is switched-off (default for symmetric).
  183. // Try m_iparm[12] = 1 in case of inappropriate accuracy
  184. m_iparm[13] = 0; // Output: Number of perturbed pivots
  185. m_iparm[14] = 0; // Not in use
  186. m_iparm[15] = 0; // Not in use
  187. m_iparm[16] = 0; // Not in use
  188. m_iparm[17] = -1; // Output: Number of nonzeros in the factor LU
  189. m_iparm[18] = -1; // Output: Mflops for LU factorization
  190. m_iparm[19] = 0; // Output: Numbers of CG Iterations
  191. m_iparm[20] = 0; // 1x1 pivoting
  192. m_iparm[26] = 0; // No matrix checker
  193. m_iparm[27] = (sizeof(RealScalar) == 4) ? 1 : 0;
  194. m_iparm[34] = 1; // C indexing
  195. m_iparm[36] = 0; // CSR
  196. m_iparm[59] = 0; // 0 - In-Core ; 1 - Automatic switch between In-Core and Out-of-Core modes ; 2 - Out-of-Core
  197. memset(m_pt, 0, sizeof(m_pt));
  198. }
  199. protected:
  200. // cached data to reduce reallocation, etc.
  201. void manageErrorCode(Index error) const
  202. {
  203. switch(error)
  204. {
  205. case 0:
  206. m_info = Success;
  207. break;
  208. case -4:
  209. case -7:
  210. m_info = NumericalIssue;
  211. break;
  212. default:
  213. m_info = InvalidInput;
  214. }
  215. }
  216. mutable SparseMatrixType m_matrix;
  217. mutable ComputationInfo m_info;
  218. bool m_analysisIsOk, m_factorizationIsOk;
  219. StorageIndex m_type, m_msglvl;
  220. mutable void *m_pt[64];
  221. mutable ParameterType m_iparm;
  222. mutable IntColVectorType m_perm;
  223. Index m_size;
  224. };
  225. template<class Derived>
  226. Derived& PardisoImpl<Derived>::compute(const MatrixType& a)
  227. {
  228. m_size = a.rows();
  229. eigen_assert(a.rows() == a.cols());
  230. pardisoRelease();
  231. m_perm.setZero(m_size);
  232. derived().getMatrix(a);
  233. Index error;
  234. error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 12, internal::convert_index<StorageIndex>(m_size),
  235. m_matrix.valuePtr(), m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(),
  236. m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
  237. manageErrorCode(error);
  238. m_analysisIsOk = true;
  239. m_factorizationIsOk = true;
  240. m_isInitialized = true;
  241. return derived();
  242. }
  243. template<class Derived>
  244. Derived& PardisoImpl<Derived>::analyzePattern(const MatrixType& a)
  245. {
  246. m_size = a.rows();
  247. eigen_assert(m_size == a.cols());
  248. pardisoRelease();
  249. m_perm.setZero(m_size);
  250. derived().getMatrix(a);
  251. Index error;
  252. error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 11, internal::convert_index<StorageIndex>(m_size),
  253. m_matrix.valuePtr(), m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(),
  254. m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
  255. manageErrorCode(error);
  256. m_analysisIsOk = true;
  257. m_factorizationIsOk = false;
  258. m_isInitialized = true;
  259. return derived();
  260. }
  261. template<class Derived>
  262. Derived& PardisoImpl<Derived>::factorize(const MatrixType& a)
  263. {
  264. eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
  265. eigen_assert(m_size == a.rows() && m_size == a.cols());
  266. derived().getMatrix(a);
  267. Index error;
  268. error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 22, internal::convert_index<StorageIndex>(m_size),
  269. m_matrix.valuePtr(), m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(),
  270. m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
  271. manageErrorCode(error);
  272. m_factorizationIsOk = true;
  273. return derived();
  274. }
  275. template<class Derived>
  276. template<typename BDerived,typename XDerived>
  277. void PardisoImpl<Derived>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived>& x) const
  278. {
  279. if(m_iparm[0] == 0) // Factorization was not computed
  280. {
  281. m_info = InvalidInput;
  282. return;
  283. }
  284. //Index n = m_matrix.rows();
  285. Index nrhs = Index(b.cols());
  286. eigen_assert(m_size==b.rows());
  287. eigen_assert(((MatrixBase<BDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) && "Row-major right hand sides are not supported");
  288. eigen_assert(((MatrixBase<XDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) && "Row-major matrices of unknowns are not supported");
  289. eigen_assert(((nrhs == 1) || b.outerStride() == b.rows()));
  290. // switch (transposed) {
  291. // case SvNoTrans : m_iparm[11] = 0 ; break;
  292. // case SvTranspose : m_iparm[11] = 2 ; break;
  293. // case SvAdjoint : m_iparm[11] = 1 ; break;
  294. // default:
  295. // //std::cerr << "Eigen: transposition option \"" << transposed << "\" not supported by the PARDISO backend\n";
  296. // m_iparm[11] = 0;
  297. // }
  298. Scalar* rhs_ptr = const_cast<Scalar*>(b.derived().data());
  299. Matrix<Scalar,Dynamic,Dynamic,ColMajor> tmp;
  300. // Pardiso cannot solve in-place
  301. if(rhs_ptr == x.derived().data())
  302. {
  303. tmp = b;
  304. rhs_ptr = tmp.data();
  305. }
  306. Index error;
  307. error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 33, internal::convert_index<StorageIndex>(m_size),
  308. m_matrix.valuePtr(), m_matrix.outerIndexPtr(), m_matrix.innerIndexPtr(),
  309. m_perm.data(), internal::convert_index<StorageIndex>(nrhs), m_iparm.data(), m_msglvl,
  310. rhs_ptr, x.derived().data());
  311. manageErrorCode(error);
  312. }
  313. /** \ingroup PardisoSupport_Module
  314. * \class PardisoLU
  315. * \brief A sparse direct LU factorization and solver based on the PARDISO library
  316. *
  317. * This class allows to solve for A.X = B sparse linear problems via a direct LU factorization
  318. * using the Intel MKL PARDISO library. The sparse matrix A must be squared and invertible.
  319. * The vectors or matrices X and B can be either dense or sparse.
  320. *
  321. * By default, it runs in in-core mode. To enable PARDISO's out-of-core feature, set:
  322. * \code solver.pardisoParameterArray()[59] = 1; \endcode
  323. *
  324. * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
  325. *
  326. * \implsparsesolverconcept
  327. *
  328. * \sa \ref TutorialSparseSolverConcept, class SparseLU
  329. */
  330. template<typename MatrixType>
  331. class PardisoLU : public PardisoImpl< PardisoLU<MatrixType> >
  332. {
  333. protected:
  334. typedef PardisoImpl<PardisoLU> Base;
  335. using Base::pardisoInit;
  336. using Base::m_matrix;
  337. friend class PardisoImpl< PardisoLU<MatrixType> >;
  338. public:
  339. typedef typename Base::Scalar Scalar;
  340. typedef typename Base::RealScalar RealScalar;
  341. using Base::compute;
  342. using Base::solve;
  343. PardisoLU()
  344. : Base()
  345. {
  346. pardisoInit(Base::ScalarIsComplex ? 13 : 11);
  347. }
  348. explicit PardisoLU(const MatrixType& matrix)
  349. : Base()
  350. {
  351. pardisoInit(Base::ScalarIsComplex ? 13 : 11);
  352. compute(matrix);
  353. }
  354. protected:
  355. void getMatrix(const MatrixType& matrix)
  356. {
  357. m_matrix = matrix;
  358. m_matrix.makeCompressed();
  359. }
  360. };
  361. /** \ingroup PardisoSupport_Module
  362. * \class PardisoLLT
  363. * \brief A sparse direct Cholesky (LLT) factorization and solver based on the PARDISO library
  364. *
  365. * This class allows to solve for A.X = B sparse linear problems via a LL^T Cholesky factorization
  366. * using the Intel MKL PARDISO library. The sparse matrix A must be selfajoint and positive definite.
  367. * The vectors or matrices X and B can be either dense or sparse.
  368. *
  369. * By default, it runs in in-core mode. To enable PARDISO's out-of-core feature, set:
  370. * \code solver.pardisoParameterArray()[59] = 1; \endcode
  371. *
  372. * \tparam MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
  373. * \tparam UpLo can be any bitwise combination of Upper, Lower. The default is Upper, meaning only the upper triangular part has to be used.
  374. * Upper|Lower can be used to tell both triangular parts can be used as input.
  375. *
  376. * \implsparsesolverconcept
  377. *
  378. * \sa \ref TutorialSparseSolverConcept, class SimplicialLLT
  379. */
  380. template<typename MatrixType, int _UpLo>
  381. class PardisoLLT : public PardisoImpl< PardisoLLT<MatrixType,_UpLo> >
  382. {
  383. protected:
  384. typedef PardisoImpl< PardisoLLT<MatrixType,_UpLo> > Base;
  385. using Base::pardisoInit;
  386. using Base::m_matrix;
  387. friend class PardisoImpl< PardisoLLT<MatrixType,_UpLo> >;
  388. public:
  389. typedef typename Base::Scalar Scalar;
  390. typedef typename Base::RealScalar RealScalar;
  391. typedef typename Base::StorageIndex StorageIndex;
  392. enum { UpLo = _UpLo };
  393. using Base::compute;
  394. PardisoLLT()
  395. : Base()
  396. {
  397. pardisoInit(Base::ScalarIsComplex ? 4 : 2);
  398. }
  399. explicit PardisoLLT(const MatrixType& matrix)
  400. : Base()
  401. {
  402. pardisoInit(Base::ScalarIsComplex ? 4 : 2);
  403. compute(matrix);
  404. }
  405. protected:
  406. void getMatrix(const MatrixType& matrix)
  407. {
  408. // PARDISO supports only upper, row-major matrices
  409. PermutationMatrix<Dynamic,Dynamic,StorageIndex> p_null;
  410. m_matrix.resize(matrix.rows(), matrix.cols());
  411. m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().twistedBy(p_null);
  412. m_matrix.makeCompressed();
  413. }
  414. };
  415. /** \ingroup PardisoSupport_Module
  416. * \class PardisoLDLT
  417. * \brief A sparse direct Cholesky (LDLT) factorization and solver based on the PARDISO library
  418. *
  419. * This class allows to solve for A.X = B sparse linear problems via a LDL^T Cholesky factorization
  420. * using the Intel MKL PARDISO library. The sparse matrix A is assumed to be selfajoint and positive definite.
  421. * For complex matrices, A can also be symmetric only, see the \a Options template parameter.
  422. * The vectors or matrices X and B can be either dense or sparse.
  423. *
  424. * By default, it runs in in-core mode. To enable PARDISO's out-of-core feature, set:
  425. * \code solver.pardisoParameterArray()[59] = 1; \endcode
  426. *
  427. * \tparam MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
  428. * \tparam Options can be any bitwise combination of Upper, Lower, and Symmetric. The default is Upper, meaning only the upper triangular part has to be used.
  429. * Symmetric can be used for symmetric, non-selfadjoint complex matrices, the default being to assume a selfadjoint matrix.
  430. * Upper|Lower can be used to tell both triangular parts can be used as input.
  431. *
  432. * \implsparsesolverconcept
  433. *
  434. * \sa \ref TutorialSparseSolverConcept, class SimplicialLDLT
  435. */
  436. template<typename MatrixType, int Options>
  437. class PardisoLDLT : public PardisoImpl< PardisoLDLT<MatrixType,Options> >
  438. {
  439. protected:
  440. typedef PardisoImpl< PardisoLDLT<MatrixType,Options> > Base;
  441. using Base::pardisoInit;
  442. using Base::m_matrix;
  443. friend class PardisoImpl< PardisoLDLT<MatrixType,Options> >;
  444. public:
  445. typedef typename Base::Scalar Scalar;
  446. typedef typename Base::RealScalar RealScalar;
  447. typedef typename Base::StorageIndex StorageIndex;
  448. using Base::compute;
  449. enum { UpLo = Options&(Upper|Lower) };
  450. PardisoLDLT()
  451. : Base()
  452. {
  453. pardisoInit(Base::ScalarIsComplex ? ( bool(Options&Symmetric) ? 6 : -4 ) : -2);
  454. }
  455. explicit PardisoLDLT(const MatrixType& matrix)
  456. : Base()
  457. {
  458. pardisoInit(Base::ScalarIsComplex ? ( bool(Options&Symmetric) ? 6 : -4 ) : -2);
  459. compute(matrix);
  460. }
  461. void getMatrix(const MatrixType& matrix)
  462. {
  463. // PARDISO supports only upper, row-major matrices
  464. PermutationMatrix<Dynamic,Dynamic,StorageIndex> p_null;
  465. m_matrix.resize(matrix.rows(), matrix.cols());
  466. m_matrix.template selfadjointView<Upper>() = matrix.template selfadjointView<UpLo>().twistedBy(p_null);
  467. m_matrix.makeCompressed();
  468. }
  469. };
  470. } // end namespace Eigen
  471. #endif // EIGEN_PARDISOSUPPORT_H