KASKADE 7 development version
|
Factorization of sparse linear systems with UMFPack. More...
#include <umfpack_solve.hh>
Factorization of sparse linear systems with UMFPack.
Scalar | the scalar matrix entry type |
UMFPackIndex | the integer type to be used by UMFPack. Either long or int. |
Definition at line 88 of file umfpack_solve.hh.
Public Types | |
typedef Scalar | field_type |
The type of matrix elements (a field type). More... | |
using | Options = FactorizationOptions |
Public Member Functions | |
UMFFactorization (UMFPackIndex n_, std::vector< UMFPackIndex > const &ridx, std::vector< UMFPackIndex > const &cidx, std::vector< Scalar > const &values, FactorizationOptions opt) | |
Constructor keeping input data in triplet format (aka coordinate format) constant. Construction is factorization! More... | |
UMFFactorization (UMFPackIndex n_, std::unique_ptr< std::vector< UMFPackIndex > > const ridx, std::unique_ptr< std::vector< UMFPackIndex > > const cidx, std::unique_ptr< std::vector< Scalar > > const values, FactorizationOptions opt) | |
Constructor destroying input data before factorization: more memory efficient. More... | |
void | solve (std::vector< Scalar > const &b, std::vector< Scalar > &x, bool transposed=false) const |
void | solve (Scalar const *b, Scalar *x, bool transposed=false) const |
Solves the system \( Ax=b \) for the given right hand side \( b \). More... | |
void | solve (std::vector< Scalar > const &b, std::vector< Scalar > &x, int nr, bool transposed=false) const |
void | solve (std::vector< Scalar > &b) const |
Solves the system for the given right hand side. More... | |
void | solve (Scalar *b) const |
Solves the system for the given right hand side. More... | |
virtual size_t | size () const |
reports the dimension of the system More... | |
void | setVerbose (int verbose_) |
int | getVerbose () const |
Info const & | info () const |
Protected Attributes | |
Info | info_ |
|
inherited |
The type of matrix elements (a field type).
Definition at line 48 of file factorization.hh.
|
inherited |
Definition at line 50 of file factorization.hh.
|
inline |
Constructor keeping input data in triplet format (aka coordinate format) constant. Construction is factorization!
n | size of the (square) matrix, i.e. the number of rows |
ridx | row indices |
cidx | column indices |
values | entry values |
Definition at line 101 of file umfpack_solve.hh.
|
inline |
Constructor destroying input data before factorization: more memory efficient.
Construction is factorization!
n | size of the (square) matrix, i.e. the number of rows |
ridx | row indices |
cidx | column indices |
values | entry values |
Definition at line 129 of file umfpack_solve.hh.
|
inlineinherited |
Definition at line 115 of file factorization.hh.
Referenced by Kaskade::MUMPSFactorization< Scalar >::MUMPSFactorization(), and Kaskade::SUPERLUFactorization< Scalar >::SUPERLUFactorization().
|
inlineinherited |
Definition at line 117 of file factorization.hh.
|
inlineinherited |
Definition at line 114 of file factorization.hh.
Referenced by Kaskade::MUMPSFactorization< Scalar >::MUMPSFactorization().
|
inlinevirtual |
reports the dimension of the system
Implements Kaskade::Factorization< Scalar >.
Definition at line 208 of file umfpack_solve.hh.
|
inlinevirtual |
Solves the system for the given right hand side.
This is less efficient than providing the solution vector separately.
Implements Kaskade::Factorization< Scalar >.
Definition at line 198 of file umfpack_solve.hh.
|
inlinevirtual |
Solves the system \( Ax=b \) for the given right hand side \( b \).
[in] | b | the right hand side |
[out] | x | the solution. x must point to a memory region of length of system dimension. |
Implements Kaskade::Factorization< Scalar >.
Definition at line 165 of file umfpack_solve.hh.
|
inlinevirtual |
Solves the system for the given right hand side.
This is less efficient than providing the solution vector separately.
Reimplemented from Kaskade::Factorization< Scalar >.
Definition at line 184 of file umfpack_solve.hh.
|
inlinevirtual |
Solves the system for the given right hand side
Reimplemented from Kaskade::Factorization< Scalar >.
Definition at line 158 of file umfpack_solve.hh.
Referenced by Kaskade::UMFFactorization< Scalar, UMFPackIndex >::solve(), and Kaskade::Limex< Eq >::step().
|
inline |
Definition at line 170 of file umfpack_solve.hh.
|
protectedinherited |
Definition at line 126 of file factorization.hh.
Referenced by Kaskade::MUMPSFactorization< Scalar >::factorize(), and Kaskade::Factorization< Scalar >::info().