|
DOLFINx 0.11.0.0
DOLFINx C++ interface
|
SuperLU_DIST linear solver interface. More...
#include <superlu_dist.h>
Public Member Functions | |
| SuperLUDistSolver (std::shared_ptr< const SuperLUDistMatrix< T > > A) | |
| Create solver for a SuperLU_DIST matrix operator. | |
| SuperLUDistSolver (const SuperLUDistSolver &)=delete | |
| Copy constructor. | |
| SuperLUDistSolver & | operator= (const SuperLUDistSolver &)=delete |
| Copy assignment. | |
| void | set_option (std::string name, std::string value) |
| Set solver option name to value. | |
| void | set_options (SuperLUDistStructs::superlu_dist_options_t options) |
| Set all solver options (native struct). | |
| void | set_A (std::shared_ptr< const SuperLUDistMatrix< T > > A) |
| Set assembled left-hand side matrix A. | |
| int | solve (const Vector< T > &b, Vector< T > &u) const |
| Solve linear system Au = b. | |
SuperLU_DIST linear solver interface.
| SuperLUDistSolver | ( | std::shared_ptr< const SuperLUDistMatrix< T > > | A | ) |
Create solver for a SuperLU_DIST matrix operator.
Solves linear system Au = b via LU decomposition.
The SuperLU_DIST solver has options set to upstream defaults, except PrintStat (verbose solver output) set to NO.
| T | Scalar type. |
| A | Assembled left-hand side matrix. |
| void set_A | ( | std::shared_ptr< const SuperLUDistMatrix< T > > | A | ) |
Set assembled left-hand side matrix A.
For advanced use with SuperLU_DIST option Factor allowing use of previously computed permutations when solving with new matrix A.
| A | Assembled left-hand side matrix. |
| void set_option | ( | std::string | name, |
| std::string | value ) |
Set solver option name to value.
See SuperLU_DIST User's Guide for option names and values.
| name | Option name. |
| value | Option value. |
| void set_options | ( | SuperLUDistStructs::superlu_dist_options_t | options | ) |
Set all solver options (native struct).
See SuperLU_DIST User's Guide for option names and values.
Callers must complete the forward declared struct, e.g.:
| options | SuperLU_DIST option struct. |
| int solve | ( | const Vector< T > & | b, |
| la::Vector< T > & | u ) const |
Solve linear system Au = b.
| b | Right-hand side vector. |
| u | Solution vector, overwritten during solve. |