Eigen::PardisoImpl< Derived > Class Template Reference
+ Inheritance diagram for Eigen::PardisoImpl< Derived >:

Public Types

enum  {
  ScalarIsComplex ,
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef Matrix< StorageIndex, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
 
typedef Matrix< StorageIndex, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
 
typedef Traits::MatrixType MatrixType
 
typedef Array< StorageIndex, 64, 1, DontAlignParameterType
 
typedef Traits::RealScalar RealScalar
 
typedef Traits::Scalar Scalar
 
typedef SparseMatrix< Scalar, RowMajor, StorageIndexSparseMatrixType
 
typedef Traits::StorageIndex StorageIndex
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 

Public Member Functions

template<typename BDerived , typename XDerived >
void _solve_impl (const MatrixBase< BDerived > &b, MatrixBase< XDerived > &x) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
Derived & analyzePattern (const MatrixType &matrix)
 
Index cols () const
 
Derived & compute (const MatrixType &matrix)
 
Derived & factorize (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
 PardisoImpl ()
 
ParameterTypepardisoParameterArray ()
 
Index rows () const
 
 ~PardisoImpl ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Protected Types

typedef SparseSolverBase< Derived > Base
 
typedef internal::pardiso_traits< Derived > Traits
 

Protected Member Functions

Derived & derived ()
 
const Derived & derived () const
 
void manageErrorCode (Index error) const
 
void pardisoInit (int type)
 
void pardisoRelease ()
 

Protected Attributes

bool m_analysisIsOk
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
ParameterType m_iparm
 
bool m_isInitialized
 
SparseMatrixType m_matrix
 
StorageIndex m_msglvl
 
IntColVectorType m_perm
 
void * m_pt [64]
 
Index m_size
 
StorageIndex m_type
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Detailed Description

template<class Derived>
class Eigen::PardisoImpl< Derived >

Definition at line 101 of file PardisoSupport.h.

Member Typedef Documentation

◆ Base

template<class Derived >
typedef SparseSolverBase<Derived> Eigen::PardisoImpl< Derived >::Base
protected

Definition at line 104 of file PardisoSupport.h.

◆ IntColVectorType

template<class Derived >
typedef Matrix<StorageIndex, MatrixType::RowsAtCompileTime, 1> Eigen::PardisoImpl< Derived >::IntColVectorType

Definition at line 119 of file PardisoSupport.h.

◆ IntRowVectorType

template<class Derived >
typedef Matrix<StorageIndex, 1, MatrixType::ColsAtCompileTime> Eigen::PardisoImpl< Derived >::IntRowVectorType

Definition at line 118 of file PardisoSupport.h.

◆ MatrixType

template<class Derived >
typedef Traits::MatrixType Eigen::PardisoImpl< Derived >::MatrixType

Definition at line 112 of file PardisoSupport.h.

◆ ParameterType

template<class Derived >
typedef Array<StorageIndex,64,1,DontAlign> Eigen::PardisoImpl< Derived >::ParameterType

Definition at line 120 of file PardisoSupport.h.

◆ RealScalar

template<class Derived >
typedef Traits::RealScalar Eigen::PardisoImpl< Derived >::RealScalar

Definition at line 114 of file PardisoSupport.h.

◆ Scalar

template<class Derived >
typedef Traits::Scalar Eigen::PardisoImpl< Derived >::Scalar

Definition at line 113 of file PardisoSupport.h.

◆ SparseMatrixType

template<class Derived >
typedef SparseMatrix<Scalar,RowMajor,StorageIndex> Eigen::PardisoImpl< Derived >::SparseMatrixType

Definition at line 116 of file PardisoSupport.h.

◆ StorageIndex

template<class Derived >
typedef Traits::StorageIndex Eigen::PardisoImpl< Derived >::StorageIndex

Definition at line 115 of file PardisoSupport.h.

◆ Traits

template<class Derived >
typedef internal::pardiso_traits<Derived> Eigen::PardisoImpl< Derived >::Traits
protected

Definition at line 108 of file PardisoSupport.h.

◆ VectorType

template<class Derived >
typedef Matrix<Scalar,Dynamic,1> Eigen::PardisoImpl< Derived >::VectorType

Definition at line 117 of file PardisoSupport.h.

Member Enumeration Documentation

◆ anonymous enum

template<class Derived >
anonymous enum
Enumerator
ScalarIsComplex 
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 121 of file PardisoSupport.h.

Constructor & Destructor Documentation

◆ PardisoImpl()

template<class Derived >
Eigen::PardisoImpl< Derived >::PardisoImpl ( )
inline

Definition at line 127 of file PardisoSupport.h.

128  : m_analysisIsOk(false), m_factorizationIsOk(false)
129  {
130  eigen_assert((sizeof(StorageIndex) >= sizeof(_INTEGER_t) && sizeof(StorageIndex) <= 8) && "Non-supported index type");
131  m_iparm.setZero();
132  m_msglvl = 0; // No output
133  m_isInitialized = false;
134  }
#define eigen_assert(x)
Definition: Macros.h:902
Traits::StorageIndex StorageIndex
StorageIndex m_msglvl
ParameterType m_iparm
Derived & setZero(Index size)

◆ ~PardisoImpl()

template<class Derived >
Eigen::PardisoImpl< Derived >::~PardisoImpl ( )
inline

Definition at line 136 of file PardisoSupport.h.

137  {
138  pardisoRelease();
139  }

Member Function Documentation

◆ _solve_impl() [1/2]

template<class Derived >
template<typename BDerived , typename XDerived >
void Eigen::PardisoImpl< Derived >::_solve_impl ( const MatrixBase< BDerived > &  b,
MatrixBase< XDerived > &  x 
) const

Definition at line 324 of file PardisoSupport.h.

325 {
326  if(m_iparm[0] == 0) // Factorization was not computed
327  {
329  return;
330  }
331 
332  //Index n = m_matrix.rows();
333  Index nrhs = Index(b.cols());
334  eigen_assert(m_size==b.rows());
335  eigen_assert(((MatrixBase<BDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) && "Row-major right hand sides are not supported");
336  eigen_assert(((MatrixBase<XDerived>::Flags & RowMajorBit) == 0 || nrhs == 1) && "Row-major matrices of unknowns are not supported");
337  eigen_assert(((nrhs == 1) || b.outerStride() == b.rows()));
338 
339 
340 // switch (transposed) {
341 // case SvNoTrans : m_iparm[11] = 0 ; break;
342 // case SvTranspose : m_iparm[11] = 2 ; break;
343 // case SvAdjoint : m_iparm[11] = 1 ; break;
344 // default:
345 // //std::cerr << "Eigen: transposition option \"" << transposed << "\" not supported by the PARDISO backend\n";
346 // m_iparm[11] = 0;
347 // }
348 
349  Scalar* rhs_ptr = const_cast<Scalar*>(b.derived().data());
350  Matrix<Scalar,Dynamic,Dynamic,ColMajor> tmp;
351 
352  // Pardiso cannot solve in-place
353  if(rhs_ptr == x.derived().data())
354  {
355  tmp = b;
356  rhs_ptr = tmp.data();
357  }
358 
359  Index error;
360  error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 33, internal::convert_index<StorageIndex>(m_size),
362  m_perm.data(), internal::convert_index<StorageIndex>(nrhs), m_iparm.data(), m_msglvl,
363  rhs_ptr, x.derived().data());
364 
365  manageErrorCode(error);
366 }
Array< int, 3, 1 > b
StorageIndex m_type
IntColVectorType m_perm
SparseMatrixType m_matrix
ComputationInfo m_info
void manageErrorCode(Index error) const
Traits::Scalar Scalar
const Scalar * data() const
const Scalar * valuePtr() const
Definition: SparseMatrix.h:177
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:195
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:186
@ InvalidInput
Definition: Constants.h:453
const unsigned int RowMajorBit
Definition: Constants.h:68
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82

◆ _solve_impl() [2/2]

template<class Derived >
template<typename Rhs , typename Dest >
void Eigen::PardisoImpl< Derived >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const

◆ analyzePattern()

template<class Derived >
Derived & Eigen::PardisoImpl< Derived >::analyzePattern ( const MatrixType matrix)

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()

Definition at line 283 of file PardisoSupport.h.

284 {
285  m_size = a.rows();
286  eigen_assert(m_size == a.cols());
287 
288  pardisoRelease();
290  derived().getMatrix(a);
291 
292  Index error;
293  error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 11, internal::convert_index<StorageIndex>(m_size),
295  m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
296 
297  manageErrorCode(error);
298  m_analysisIsOk = true;
299  m_factorizationIsOk = false;
300  m_isInitialized = true;
301  return derived();
302 }

◆ cols()

template<class Derived >
Index Eigen::PardisoImpl< Derived >::cols ( void  ) const
inline

Definition at line 141 of file PardisoSupport.h.

141 { return m_size; }

◆ compute()

template<class Derived >
Derived & Eigen::PardisoImpl< Derived >::compute ( const MatrixType matrix)

Definition at line 262 of file PardisoSupport.h.

263 {
264  m_size = a.rows();
265  eigen_assert(a.rows() == a.cols());
266 
267  pardisoRelease();
269  derived().getMatrix(a);
270 
271  Index error;
272  error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 12, internal::convert_index<StorageIndex>(m_size),
274  m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
275  manageErrorCode(error);
276  m_analysisIsOk = true;
277  m_factorizationIsOk = true;
278  m_isInitialized = true;
279  return derived();
280 }

◆ derived() [1/2]

template<class Derived >
Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

Definition at line 83 of file SparseSolverBase.h.

83 { return *static_cast<Derived*>(this); }

◆ derived() [2/2]

template<class Derived >
const Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

Definition at line 84 of file SparseSolverBase.h.

84 { return *static_cast<const Derived*>(this); }

◆ factorize()

template<class Derived >
Derived & Eigen::PardisoImpl< Derived >::factorize ( const MatrixType matrix)

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()

Definition at line 305 of file PardisoSupport.h.

306 {
307  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
308  eigen_assert(m_size == a.rows() && m_size == a.cols());
309 
310  derived().getMatrix(a);
311 
312  Index error;
313  error = internal::pardiso_run_selector<StorageIndex>::run(m_pt, 1, 1, m_type, 22, internal::convert_index<StorageIndex>(m_size),
315  m_perm.data(), 0, m_iparm.data(), m_msglvl, NULL, NULL);
316 
317  manageErrorCode(error);
318  m_factorizationIsOk = true;
319  return derived();
320 }

◆ info()

template<class Derived >
ComputationInfo Eigen::PardisoImpl< Derived >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the matrix appears to be negative.

Definition at line 149 of file PardisoSupport.h.

150  {
151  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
152  return m_info;
153  }

◆ manageErrorCode()

template<class Derived >
void Eigen::PardisoImpl< Derived >::manageErrorCode ( Index  error) const
inlineprotected

Definition at line 234 of file PardisoSupport.h.

235  {
236  switch(error)
237  {
238  case 0:
239  m_info = Success;
240  break;
241  case -4:
242  case -7:
244  break;
245  default:
247  }
248  }
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446

◆ pardisoInit()

template<class Derived >
void Eigen::PardisoImpl< Derived >::pardisoInit ( int  type)
inlineprotected

Definition at line 195 of file PardisoSupport.h.

196  {
197  m_type = type;
198  bool symmetric = std::abs(m_type) < 10;
199  m_iparm[0] = 1; // No solver default
200  m_iparm[1] = 2; // use Metis for the ordering
201  m_iparm[2] = 0; // Reserved. Set to zero. (??Numbers of processors, value of OMP_NUM_THREADS??)
202  m_iparm[3] = 0; // No iterative-direct algorithm
203  m_iparm[4] = 0; // No user fill-in reducing permutation
204  m_iparm[5] = 0; // Write solution into x, b is left unchanged
205  m_iparm[6] = 0; // Not in use
206  m_iparm[7] = 2; // Max numbers of iterative refinement steps
207  m_iparm[8] = 0; // Not in use
208  m_iparm[9] = 13; // Perturb the pivot elements with 1E-13
209  m_iparm[10] = symmetric ? 0 : 1; // Use nonsymmetric permutation and scaling MPS
210  m_iparm[11] = 0; // Not in use
211  m_iparm[12] = symmetric ? 0 : 1; // Maximum weighted matching algorithm is switched-off (default for symmetric).
212  // Try m_iparm[12] = 1 in case of inappropriate accuracy
213  m_iparm[13] = 0; // Output: Number of perturbed pivots
214  m_iparm[14] = 0; // Not in use
215  m_iparm[15] = 0; // Not in use
216  m_iparm[16] = 0; // Not in use
217  m_iparm[17] = -1; // Output: Number of nonzeros in the factor LU
218  m_iparm[18] = -1; // Output: Mflops for LU factorization
219  m_iparm[19] = 0; // Output: Numbers of CG Iterations
220 
221  m_iparm[20] = 0; // 1x1 pivoting
222  m_iparm[26] = 0; // No matrix checker
223  m_iparm[27] = (sizeof(RealScalar) == 4) ? 1 : 0;
224  m_iparm[34] = 1; // C indexing
225  m_iparm[36] = 0; // CSR
226  m_iparm[59] = 0; // 0 - In-Core ; 1 - Automatic switch between In-Core and Out-of-Core modes ; 2 - Out-of-Core
227 
228  memset(m_pt, 0, sizeof(m_pt));
229  }
const AbsReturnType abs() const
Traits::RealScalar RealScalar

◆ pardisoParameterArray()

template<class Derived >
ParameterType& Eigen::PardisoImpl< Derived >::pardisoParameterArray ( )
inline
Warning
for advanced usage only.
Returns
a reference to the parameter array controlling PARDISO. See the PARDISO manual to know how to use it.

Definition at line 158 of file PardisoSupport.h.

159  {
160  return m_iparm;
161  }

◆ pardisoRelease()

template<class Derived >
void Eigen::PardisoImpl< Derived >::pardisoRelease ( )
inlineprotected

Definition at line 185 of file PardisoSupport.h.

186  {
187  if(m_isInitialized) // Factorization ran at least once
188  {
189  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,
190  m_iparm.data(), m_msglvl, NULL, NULL);
191  m_isInitialized = false;
192  }
193  }

◆ rows()

template<class Derived >
Index Eigen::PardisoImpl< Derived >::rows ( void  ) const
inline

Definition at line 142 of file PardisoSupport.h.

142 { return m_size; }

Member Data Documentation

◆ m_analysisIsOk

template<class Derived >
bool Eigen::PardisoImpl< Derived >::m_analysisIsOk
protected

Definition at line 252 of file PardisoSupport.h.

◆ m_factorizationIsOk

template<class Derived >
bool Eigen::PardisoImpl< Derived >::m_factorizationIsOk
protected

Definition at line 252 of file PardisoSupport.h.

◆ m_info

template<class Derived >
ComputationInfo Eigen::PardisoImpl< Derived >::m_info
mutableprotected

Definition at line 251 of file PardisoSupport.h.

◆ m_iparm

template<class Derived >
ParameterType Eigen::PardisoImpl< Derived >::m_iparm
mutableprotected

Definition at line 255 of file PardisoSupport.h.

◆ m_isInitialized

template<class Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

Definition at line 123 of file SparseSolverBase.h.

◆ m_matrix

template<class Derived >
SparseMatrixType Eigen::PardisoImpl< Derived >::m_matrix
mutableprotected

Definition at line 250 of file PardisoSupport.h.

◆ m_msglvl

template<class Derived >
StorageIndex Eigen::PardisoImpl< Derived >::m_msglvl
protected

Definition at line 253 of file PardisoSupport.h.

◆ m_perm

template<class Derived >
IntColVectorType Eigen::PardisoImpl< Derived >::m_perm
mutableprotected

Definition at line 256 of file PardisoSupport.h.

◆ m_pt

template<class Derived >
void* Eigen::PardisoImpl< Derived >::m_pt[64]
mutableprotected

Definition at line 254 of file PardisoSupport.h.

◆ m_size

template<class Derived >
Index Eigen::PardisoImpl< Derived >::m_size
protected

Definition at line 257 of file PardisoSupport.h.

◆ m_type

template<class Derived >
StorageIndex Eigen::PardisoImpl< Derived >::m_type
protected

Definition at line 253 of file PardisoSupport.h.


The documentation for this class was generated from the following file: