Sparse direct LU solver based on PaStiX library.
More...
#include <PaStiXSupport.h>
Inherits PastixBase< Derived >.
|
void | analyzePattern (const MatrixType &matrix) |
|
void | compute (const MatrixType &matrix) |
|
Array< RealScalar, IPARM_SIZE, 1 > & | dparm () |
|
double & | dparm (int idxparam) |
|
void | factorize (const MatrixType &matrix) |
|
ComputationInfo | info () const |
| Reports whether previous computation was successful.
|
|
Array< Index, IPARM_SIZE, 1 > & | iparm () |
|
int & | iparm (int idxparam) |
|
template<typename Rhs > |
const internal::solve_retval
< PastixBase, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
template<typename Rhs > |
const
internal::sparse_solve_retval
< PastixBase, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
template<typename _MatrixType, bool IsStrSym>
class Eigen::PastixLU< _MatrixType, IsStrSym >
Sparse direct LU solver based on PaStiX library.
This class is used to solve the linear systems A.X = B with a supernodal LU factorization in the PaStiX library. The matrix A should be squared and nonsingular PaStiX requires that the matrix A has a symmetric structural pattern. This interface can symmetrize the input matrix otherwise. The vectors or matrices X and B can be either dense or sparse.
- Template Parameters
-
_MatrixType | the type of the sparse matrix A, it must be a SparseMatrix<> |
IsStrSym | Indicates if the input matrix has a symmetric pattern, default is false NOTE : Note that if the analysis and factorization phase are called separately, the input matrix will be symmetrized at each call, hence it is advised to symmetrize the matrix in a end-user program and set IsStrSym to true |
- See Also
- Solving linear problems
void analyzePattern |
( |
const MatrixType & |
matrix | ) |
|
|
inline |
Compute the LU symbolic factorization of matrix
using its sparsity pattern. Several ordering methods can be used at this step. See the PaStiX user's manual. The result of this operation can be used with successive matrices having the same pattern as matrix
- See Also
- factorize()
void compute |
( |
const MatrixType & |
matrix | ) |
|
|
inline |
Compute the LU supernodal factorization of matrix
. iparm and dparm can be used to tune the PaStiX parameters. see the PaStiX user's manual
- See Also
- analyzePattern() factorize()
Array<RealScalar,IPARM_SIZE,1>& dparm |
( |
| ) |
|
|
inlineinherited |
Returns a reference to the double vector DPARM of PaStiX parameters The statistics related to the different phases of factorization and solve are saved here as well
- See Also
- analyzePattern() factorize()
double& dparm |
( |
int |
idxparam | ) |
|
|
inlineinherited |
Return a reference to a particular index parameter of the DPARM vector
- See Also
- dparm()
void factorize |
( |
const MatrixType & |
matrix | ) |
|
|
inline |
Compute the LU supernodal factorization of matrix
WARNING The matrix matrix
should have the same structural pattern as the same used in the analysis phase.
- See Also
- analyzePattern()
Reports whether previous computation was successful.
- Returns
Success
if computation was succesful, NumericalIssue
if the PaStiX reports a problem InvalidInput
if the input matrix is invalid
- See Also
- iparm()
Array<Index,IPARM_SIZE,1>& iparm |
( |
| ) |
|
|
inlineinherited |
Returns a reference to the integer vector IPARM of PaStiX parameters to modify the default parameters. The statistics related to the different phases of factorization and solve are saved here as well
- See Also
- analyzePattern() factorize()
int& iparm |
( |
int |
idxparam | ) |
|
|
inlineinherited |
Return a reference to a particular index parameter of the IPARM vector
- See Also
- iparm()
const internal::solve_retval<PastixBase, Rhs> solve |
( |
const MatrixBase< Rhs > & |
b | ) |
const |
|
inlineinherited |
- Returns
- the solution x of
using the current decomposition of A.
- See Also
- compute()
const internal::sparse_solve_retval<PastixBase, Rhs> solve |
( |
const SparseMatrixBase< Rhs > & |
b | ) |
const |
|
inlineinherited |
- Returns
- the solution x of
using the current decomposition of A.
- See Also
- compute()
The documentation for this class was generated from the following file: