$darkmode
Eigen  5.0.1-dev
PaStiXSupport.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_PASTIXSUPPORT_H
11 #define EIGEN_PASTIXSUPPORT_H
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 #if defined(DCOMPLEX)
19 #define PASTIX_COMPLEX COMPLEX
20 #define PASTIX_DCOMPLEX DCOMPLEX
21 #else
22 #define PASTIX_COMPLEX std::complex<float>
23 #define PASTIX_DCOMPLEX std::complex<double>
24 #endif
25 
34 template <typename MatrixType_, bool IsStrSym = false>
35 class PastixLU;
36 template <typename MatrixType_, int Options>
37 class PastixLLT;
38 template <typename MatrixType_, int Options>
39 class PastixLDLT;
40 
41 namespace internal {
42 
43 template <class Pastix>
44 struct pastix_traits;
45 
46 template <typename MatrixType_>
47 struct pastix_traits<PastixLU<MatrixType_> > {
48  typedef MatrixType_ MatrixType;
49  typedef typename MatrixType_::Scalar Scalar;
50  typedef typename MatrixType_::RealScalar RealScalar;
51  typedef typename MatrixType_::StorageIndex StorageIndex;
52 };
53 
54 template <typename MatrixType_, int Options>
55 struct pastix_traits<PastixLLT<MatrixType_, Options> > {
56  typedef MatrixType_ MatrixType;
57  typedef typename MatrixType_::Scalar Scalar;
58  typedef typename MatrixType_::RealScalar RealScalar;
59  typedef typename MatrixType_::StorageIndex StorageIndex;
60 };
61 
62 template <typename MatrixType_, int Options>
63 struct pastix_traits<PastixLDLT<MatrixType_, Options> > {
64  typedef MatrixType_ MatrixType;
65  typedef typename MatrixType_::Scalar Scalar;
66  typedef typename MatrixType_::RealScalar RealScalar;
67  typedef typename MatrixType_::StorageIndex StorageIndex;
68 };
69 
70 inline void eigen_pastix(pastix_data_t **pastix_data, int pastix_comm, int n, int *ptr, int *idx, float *vals,
71  int *perm, int *invp, float *x, int nbrhs, int *iparm, double *dparm) {
72  if (n == 0) {
73  ptr = NULL;
74  idx = NULL;
75  vals = NULL;
76  }
77  if (nbrhs == 0) {
78  x = NULL;
79  nbrhs = 1;
80  }
81  s_pastix(pastix_data, pastix_comm, n, ptr, idx, vals, perm, invp, x, nbrhs, iparm, dparm);
82 }
83 
84 inline void eigen_pastix(pastix_data_t **pastix_data, int pastix_comm, int n, int *ptr, int *idx, double *vals,
85  int *perm, int *invp, double *x, int nbrhs, int *iparm, double *dparm) {
86  if (n == 0) {
87  ptr = NULL;
88  idx = NULL;
89  vals = NULL;
90  }
91  if (nbrhs == 0) {
92  x = NULL;
93  nbrhs = 1;
94  }
95  d_pastix(pastix_data, pastix_comm, n, ptr, idx, vals, perm, invp, x, nbrhs, iparm, dparm);
96 }
97 
98 inline void eigen_pastix(pastix_data_t **pastix_data, int pastix_comm, int n, int *ptr, int *idx,
99  std::complex<float> *vals, int *perm, int *invp, std::complex<float> *x, int nbrhs, int *iparm,
100  double *dparm) {
101  if (n == 0) {
102  ptr = NULL;
103  idx = NULL;
104  vals = NULL;
105  }
106  if (nbrhs == 0) {
107  x = NULL;
108  nbrhs = 1;
109  }
110  c_pastix(pastix_data, pastix_comm, n, ptr, idx, reinterpret_cast<PASTIX_COMPLEX *>(vals), perm, invp,
111  reinterpret_cast<PASTIX_COMPLEX *>(x), nbrhs, iparm, dparm);
112 }
113 
114 inline void eigen_pastix(pastix_data_t **pastix_data, int pastix_comm, int n, int *ptr, int *idx,
115  std::complex<double> *vals, int *perm, int *invp, std::complex<double> *x, int nbrhs,
116  int *iparm, double *dparm) {
117  if (n == 0) {
118  ptr = NULL;
119  idx = NULL;
120  vals = NULL;
121  }
122  if (nbrhs == 0) {
123  x = NULL;
124  nbrhs = 1;
125  }
126  z_pastix(pastix_data, pastix_comm, n, ptr, idx, reinterpret_cast<PASTIX_DCOMPLEX *>(vals), perm, invp,
127  reinterpret_cast<PASTIX_DCOMPLEX *>(x), nbrhs, iparm, dparm);
128 }
129 
130 // Convert the matrix to Fortran-style Numbering
131 template <typename MatrixType>
132 void c_to_fortran_numbering(MatrixType &mat) {
133  if (!(mat.outerIndexPtr()[0])) {
134  int i;
135  for (i = 0; i <= mat.rows(); ++i) ++mat.outerIndexPtr()[i];
136  for (i = 0; i < mat.nonZeros(); ++i) ++mat.innerIndexPtr()[i];
137  }
138 }
139 
140 // Convert to C-style Numbering
141 template <typename MatrixType>
142 void fortran_to_c_numbering(MatrixType &mat) {
143  // Check the Numbering
144  if (mat.outerIndexPtr()[0] == 1) { // Convert to C-style numbering
145  int i;
146  for (i = 0; i <= mat.rows(); ++i) --mat.outerIndexPtr()[i];
147  for (i = 0; i < mat.nonZeros(); ++i) --mat.innerIndexPtr()[i];
148  }
149 }
150 } // namespace internal
151 
152 // This is the base class to interface with PaStiX functions.
153 // Users should not used this class directly.
154 template <class Derived>
155 class PastixBase : public SparseSolverBase<Derived> {
156  protected:
157  typedef SparseSolverBase<Derived> Base;
158  using Base::derived;
159  using Base::m_isInitialized;
160 
161  public:
162  using Base::_solve_impl;
163 
164  typedef typename internal::pastix_traits<Derived>::MatrixType MatrixType_;
165  typedef MatrixType_ MatrixType;
166  typedef typename MatrixType::Scalar Scalar;
167  typedef typename MatrixType::RealScalar RealScalar;
168  typedef typename MatrixType::StorageIndex StorageIndex;
169  typedef Matrix<Scalar, Dynamic, 1> Vector;
170  typedef SparseMatrix<Scalar, ColMajor> ColSpMatrix;
171  enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
172 
173  public:
174  PastixBase() : m_initisOk(false), m_analysisIsOk(false), m_factorizationIsOk(false), m_pastixdata(0), m_size(0) {
175  init();
176  }
177 
178  ~PastixBase() { clean(); }
179 
180  template <typename Rhs, typename Dest>
181  bool _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &x) const;
182 
188  Array<StorageIndex, IPARM_SIZE, 1> &iparm() { return m_iparm; }
189 
194  int &iparm(int idxparam) { return m_iparm(idxparam); }
195 
200  Array<double, DPARM_SIZE, 1> &dparm() { return m_dparm; }
201 
205  double &dparm(int idxparam) { return m_dparm(idxparam); }
206 
207  inline Index cols() const { return m_size; }
208  inline Index rows() const { return m_size; }
209 
218  ComputationInfo info() const {
219  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
220  return m_info;
221  }
222 
223  protected:
224  // Initialize the Pastix data structure, check the matrix
225  void init();
226 
227  // Compute the ordering and the symbolic factorization
228  void analyzePattern(ColSpMatrix &mat);
229 
230  // Compute the numerical factorization
231  void factorize(ColSpMatrix &mat);
232 
233  // Free all the data allocated by Pastix
234  void clean() {
235  eigen_assert(m_initisOk && "The Pastix structure should be allocated first");
236  m_iparm(IPARM_START_TASK) = API_TASK_CLEAN;
237  m_iparm(IPARM_END_TASK) = API_TASK_CLEAN;
238  internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, (Scalar *)0, m_perm.data(), m_invp.data(), 0, 0,
239  m_iparm.data(), m_dparm.data());
240  }
241 
242  void compute(ColSpMatrix &mat);
243 
244  int m_initisOk;
245  int m_analysisIsOk;
246  int m_factorizationIsOk;
247  mutable ComputationInfo m_info;
248  mutable pastix_data_t *m_pastixdata; // Data structure for pastix
249  mutable int m_comm; // The MPI communicator identifier
250  mutable Array<int, IPARM_SIZE, 1> m_iparm; // integer vector for the input parameters
251  mutable Array<double, DPARM_SIZE, 1> m_dparm; // Scalar vector for the input parameters
252  mutable Matrix<StorageIndex, Dynamic, 1> m_perm; // Permutation vector
253  mutable Matrix<StorageIndex, Dynamic, 1> m_invp; // Inverse permutation vector
254  mutable int m_size; // Size of the matrix
255 };
256 
261 template <class Derived>
262 void PastixBase<Derived>::init() {
263  m_size = 0;
264  m_iparm.setZero(IPARM_SIZE);
265  m_dparm.setZero(DPARM_SIZE);
266 
267  m_iparm(IPARM_MODIFY_PARAMETER) = API_NO;
268  pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, 0, 0, 0, 0, 1, m_iparm.data(), m_dparm.data());
269 
270  m_iparm[IPARM_MATRIX_VERIFICATION] = API_NO;
271  m_iparm[IPARM_VERBOSE] = API_VERBOSE_NOT;
272  m_iparm[IPARM_ORDERING] = API_ORDER_SCOTCH;
273  m_iparm[IPARM_INCOMPLETE] = API_NO;
274  m_iparm[IPARM_OOC_LIMIT] = 2000;
275  m_iparm[IPARM_RHS_MAKING] = API_RHS_B;
276  m_iparm(IPARM_MATRIX_VERIFICATION) = API_NO;
277 
278  m_iparm(IPARM_START_TASK) = API_TASK_INIT;
279  m_iparm(IPARM_END_TASK) = API_TASK_INIT;
280  internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, 0, 0, 0, (Scalar *)0, 0, 0, 0, 0, m_iparm.data(),
281  m_dparm.data());
282 
283  // Check the returned error
284  if (m_iparm(IPARM_ERROR_NUMBER)) {
285  m_info = InvalidInput;
286  m_initisOk = false;
287  } else {
288  m_info = Success;
289  m_initisOk = true;
290  }
291 }
292 
293 template <class Derived>
294 void PastixBase<Derived>::compute(ColSpMatrix &mat) {
295  eigen_assert(mat.rows() == mat.cols() && "The input matrix should be squared");
296 
297  analyzePattern(mat);
298  factorize(mat);
299 
300  m_iparm(IPARM_MATRIX_VERIFICATION) = API_NO;
301 }
302 
303 template <class Derived>
304 void PastixBase<Derived>::analyzePattern(ColSpMatrix &mat) {
305  eigen_assert(m_initisOk && "The initialization of PaSTiX failed");
306 
307  // clean previous calls
308  if (m_size > 0) clean();
309 
310  m_size = internal::convert_index<int>(mat.rows());
311  m_perm.resize(m_size);
312  m_invp.resize(m_size);
313 
314  m_iparm(IPARM_START_TASK) = API_TASK_ORDERING;
315  m_iparm(IPARM_END_TASK) = API_TASK_ANALYSE;
316  internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, m_size, mat.outerIndexPtr(), mat.innerIndexPtr(),
317  mat.valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
318 
319  // Check the returned error
320  if (m_iparm(IPARM_ERROR_NUMBER)) {
321  m_info = NumericalIssue;
322  m_analysisIsOk = false;
323  } else {
324  m_info = Success;
325  m_analysisIsOk = true;
326  }
327 }
328 
329 template <class Derived>
330 void PastixBase<Derived>::factorize(ColSpMatrix &mat) {
331  // if(&m_cpyMat != &mat) m_cpyMat = mat;
332  eigen_assert(m_analysisIsOk && "The analysis phase should be called before the factorization phase");
333  m_iparm(IPARM_START_TASK) = API_TASK_NUMFACT;
334  m_iparm(IPARM_END_TASK) = API_TASK_NUMFACT;
335  m_size = internal::convert_index<int>(mat.rows());
336 
337  internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, m_size, mat.outerIndexPtr(), mat.innerIndexPtr(),
338  mat.valuePtr(), m_perm.data(), m_invp.data(), 0, 0, m_iparm.data(), m_dparm.data());
339 
340  // Check the returned error
341  if (m_iparm(IPARM_ERROR_NUMBER)) {
342  m_info = NumericalIssue;
343  m_factorizationIsOk = false;
344  m_isInitialized = false;
345  } else {
346  m_info = Success;
347  m_factorizationIsOk = true;
348  m_isInitialized = true;
349  }
350 }
351 
352 /* Solve the system */
353 template <typename Base>
354 template <typename Rhs, typename Dest>
355 bool PastixBase<Base>::_solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &x) const {
356  eigen_assert(m_isInitialized && "The matrix should be factorized first");
357  EIGEN_STATIC_ASSERT((Dest::Flags & RowMajorBit) == 0, THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
358  int rhs = 1;
359 
360  x = b; /* on return, x is overwritten by the computed solution */
361 
362  for (int i = 0; i < b.cols(); i++) {
363  m_iparm[IPARM_START_TASK] = API_TASK_SOLVE;
364  m_iparm[IPARM_END_TASK] = API_TASK_REFINE;
365 
366  internal::eigen_pastix(&m_pastixdata, MPI_COMM_WORLD, internal::convert_index<int>(x.rows()), 0, 0, 0,
367  m_perm.data(), m_invp.data(), &x(0, i), rhs, m_iparm.data(), m_dparm.data());
368  }
369 
370  // Check the returned error
371  m_info = m_iparm(IPARM_ERROR_NUMBER) == 0 ? Success : NumericalIssue;
372 
373  return m_iparm(IPARM_ERROR_NUMBER) == 0;
374 }
375 
397 template <typename MatrixType_, bool IsStrSym>
398 class PastixLU : public PastixBase<PastixLU<MatrixType_> > {
399  public:
400  typedef MatrixType_ MatrixType;
401  typedef PastixBase<PastixLU<MatrixType> > Base;
402  typedef typename Base::ColSpMatrix ColSpMatrix;
403  typedef typename MatrixType::StorageIndex StorageIndex;
404 
405  public:
406  PastixLU() : Base() { init(); }
407 
408  explicit PastixLU(const MatrixType &matrix) : Base() {
409  init();
410  compute(matrix);
411  }
417  void compute(const MatrixType &matrix) {
418  m_structureIsUptodate = false;
419  ColSpMatrix temp;
420  grabMatrix(matrix, temp);
421  Base::compute(temp);
422  }
428  void analyzePattern(const MatrixType &matrix) {
429  m_structureIsUptodate = false;
430  ColSpMatrix temp;
431  grabMatrix(matrix, temp);
432  Base::analyzePattern(temp);
433  }
434 
440  void factorize(const MatrixType &matrix) {
441  ColSpMatrix temp;
442  grabMatrix(matrix, temp);
443  Base::factorize(temp);
444  }
445 
446  protected:
447  void init() {
448  m_structureIsUptodate = false;
449  m_iparm(IPARM_SYM) = API_SYM_NO;
450  m_iparm(IPARM_FACTORIZATION) = API_FACT_LU;
451  }
452 
453  void grabMatrix(const MatrixType &matrix, ColSpMatrix &out) {
454  if (IsStrSym)
455  out = matrix;
456  else {
457  if (!m_structureIsUptodate) {
458  // update the transposed structure
459  m_transposedStructure = matrix.transpose();
460 
461  // Set the elements of the matrix to zero
462  for (Index j = 0; j < m_transposedStructure.outerSize(); ++j)
463  for (typename ColSpMatrix::InnerIterator it(m_transposedStructure, j); it; ++it) it.valueRef() = 0.0;
464 
465  m_structureIsUptodate = true;
466  }
467 
468  out = m_transposedStructure + matrix;
469  }
470  internal::c_to_fortran_numbering(out);
471  }
472 
473  using Base::m_dparm;
474  using Base::m_iparm;
475 
476  ColSpMatrix m_transposedStructure;
477  bool m_structureIsUptodate;
478 };
479 
496 template <typename MatrixType_, int UpLo_>
497 class PastixLLT : public PastixBase<PastixLLT<MatrixType_, UpLo_> > {
498  public:
499  typedef MatrixType_ MatrixType;
500  typedef PastixBase<PastixLLT<MatrixType, UpLo_> > Base;
501  typedef typename Base::ColSpMatrix ColSpMatrix;
502 
503  public:
504  enum { UpLo = UpLo_ };
505  PastixLLT() : Base() { init(); }
506 
507  explicit PastixLLT(const MatrixType &matrix) : Base() {
508  init();
509  compute(matrix);
510  }
511 
515  void compute(const MatrixType &matrix) {
516  ColSpMatrix temp;
517  grabMatrix(matrix, temp);
518  Base::compute(temp);
519  }
520 
525  void analyzePattern(const MatrixType &matrix) {
526  ColSpMatrix temp;
527  grabMatrix(matrix, temp);
528  Base::analyzePattern(temp);
529  }
533  void factorize(const MatrixType &matrix) {
534  ColSpMatrix temp;
535  grabMatrix(matrix, temp);
536  Base::factorize(temp);
537  }
538 
539  protected:
540  using Base::m_iparm;
541 
542  void init() {
543  m_iparm(IPARM_SYM) = API_SYM_YES;
544  m_iparm(IPARM_FACTORIZATION) = API_FACT_LLT;
545  }
546 
547  void grabMatrix(const MatrixType &matrix, ColSpMatrix &out) {
548  out.resize(matrix.rows(), matrix.cols());
549  // Pastix supports only lower, column-major matrices
550  out.template selfadjointView<Lower>() = matrix.template selfadjointView<UpLo>();
551  internal::c_to_fortran_numbering(out);
552  }
553 };
554 
571 template <typename MatrixType_, int UpLo_>
572 class PastixLDLT : public PastixBase<PastixLDLT<MatrixType_, UpLo_> > {
573  public:
574  typedef MatrixType_ MatrixType;
575  typedef PastixBase<PastixLDLT<MatrixType, UpLo_> > Base;
576  typedef typename Base::ColSpMatrix ColSpMatrix;
577 
578  public:
579  enum { UpLo = UpLo_ };
580  PastixLDLT() : Base() { init(); }
581 
582  explicit PastixLDLT(const MatrixType &matrix) : Base() {
583  init();
584  compute(matrix);
585  }
586 
590  void compute(const MatrixType &matrix) {
591  ColSpMatrix temp;
592  grabMatrix(matrix, temp);
593  Base::compute(temp);
594  }
595 
600  void analyzePattern(const MatrixType &matrix) {
601  ColSpMatrix temp;
602  grabMatrix(matrix, temp);
603  Base::analyzePattern(temp);
604  }
608  void factorize(const MatrixType &matrix) {
609  ColSpMatrix temp;
610  grabMatrix(matrix, temp);
611  Base::factorize(temp);
612  }
613 
614  protected:
615  using Base::m_iparm;
616 
617  void init() {
618  m_iparm(IPARM_SYM) = API_SYM_YES;
619  m_iparm(IPARM_FACTORIZATION) = API_FACT_LDLT;
620  }
621 
622  void grabMatrix(const MatrixType &matrix, ColSpMatrix &out) {
623  // Pastix supports only lower, column-major matrices
624  out.resize(matrix.rows(), matrix.cols());
625  out.template selfadjointView<Lower>() = matrix.template selfadjointView<UpLo>();
626  internal::c_to_fortran_numbering(out);
627  }
628 };
629 
630 } // end namespace Eigen
631 
632 #endif
void analyzePattern(const MatrixType &matrix)
Definition: PaStiXSupport.h:600
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
void factorize(const MatrixType &matrix)
Definition: PaStiXSupport.h:440
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
void analyzePattern(const MatrixType &matrix)
Definition: PaStiXSupport.h:428
const unsigned int RowMajorBit
Definition: Constants.h:70
void compute(const MatrixType &matrix)
Definition: PaStiXSupport.h:515
void compute(const MatrixType &matrix)
Definition: PaStiXSupport.h:417
Matrix< Type, Size, 1 > Vector
[c++11] Size×1 vector of type Type.
Definition: Matrix.h:522
void compute(const MatrixType &matrix)
Definition: PaStiXSupport.h:590
A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library...
Definition: PaStiXSupport.h:39
Definition: Constants.h:442
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Index outerSize() const
Definition: SparseMatrix.h:166
void analyzePattern(const MatrixType &matrix)
Definition: PaStiXSupport.h:525
Definition: Constants.h:447
Definition: Constants.h:440
Interface to the PaStix solver.
Definition: PaStiXSupport.h:35
A sparse direct supernodal Cholesky (LLT) factorization and solver based on the PaStiX library...
Definition: PaStiXSupport.h:37
void factorize(const MatrixType &matrix)
Definition: PaStiXSupport.h:608
ComputationInfo
Definition: Constants.h:438
void factorize(const MatrixType &matrix)
Definition: PaStiXSupport.h:533