INMOST
A toolkit for distributed mathematical modeling
|
Main class to set and solve linear system. More...
#include <inmost_solver.h>
Classes | |
class | OrderInfo |
Base class for low level parallel operations with objects of Solver class. More... | |
Public Types | |
typedef std::string | Type |
Public Member Functions | |
void | SetParameterEnum (std::string name, INMOST_DATA_ENUM_TYPE value) |
Backward-compatibility access to integer-typed parameters. More... | |
void | SetParameterReal (std::string name, INMOST_DATA_REAL_TYPE value) |
Backward-compatibility access to real-typed parameters. More... | |
std::string | GetReason () |
Backward-compatibility acces to return reason. | |
Solver (std::string solverName, std::string prefix="", INMOST_MPI_Comm _comm=INMOST_MPI_COMM_WORLD) | |
Main constructor of the solver. More... | |
Solver (const Solver &other) | |
Copy a solver. More... | |
Solver & | operator= (const Solver &other) |
Assign a solver. More... | |
std::string | SolverName () const |
Return the solver name. More... | |
std::string | SolverPrefix () const |
Return the solver user specified name of the current solver. More... | |
SolverVerbosityLevel | VerbosityLevel () const |
Return the solver verbosity specified level of the current solver. More... | |
void | SetVerbosityLevel (SolverVerbosityLevel level) |
void | SetMatrix (Sparse::Matrix &A, bool ModifiedPattern=true, bool OldPreconditioner=false) |
Set the matrix and construct the preconditioner. More... | |
bool | Solve (Sparse::Vector &RHS, Sparse::Vector &SOL) |
Solver the linear system: A*x = b. More... | |
bool | Clear () |
Clear all internal data of the current solver including matrix, preconditioner etc. | |
std::string | GetParameter (std::string name) const |
Get the solver output parameter. More... | |
void | SetParameter (std::string name, std::string value) |
INMOST_DATA_ENUM_TYPE | Iterations () const |
Return the number of iterations performed by the last solution. More... | |
INMOST_DATA_REAL_TYPE | Residual () const |
Return the final precondioned residual achieved by the last solution. More... | |
const std::string | ReturnReason () const |
Get the reason of convergence or divergence of the last solution. More... | |
double | PreconditionerTime () const |
Return the time of preconditioner evaluation by the last solution. More... | |
double | IterationsTime () const |
Return the time of iteration stage by the last solution. More... | |
bool | IsSolved () const |
Return the last solution successfullness. More... | |
const std::string | SolutionMetadataLine (const std::string &delimiter) const |
Get the useful info of the last solution. More... | |
INMOST_DATA_REAL_TYPE | Condest (INMOST_DATA_REAL_TYPE tol, INMOST_DATA_ENUM_TYPE maxits=100) |
Computes the smallest and the largest eigenvalue with the power method. More... | |
~Solver () | |
Delete solver and associated data. | |
Static Public Member Functions | |
static void | Initialize (int *argc, char ***argv, const char *database=NULL) |
Initialize the stage of parallel solution. More... | |
static void | Finalize () |
Finalize the stage of parallel solution. More... | |
static bool | isInitialized () |
Checks the stage of parallel solution is initialized. | |
static bool | isFinalized () |
Checks the stage of parallel solution is finalized. | |
static bool | isSolverAvailable (std::string name) |
Checks if solver available. More... | |
static std::vector< std::string > | getAvailableSolvers () |
Return the list of all available solvers. More... | |
Static Public Attributes | |
static const Type | INNER_ILU2 |
inner Solver based on BiCGStab(L) solver with second order ILU factorization as preconditioner. | |
static const Type | INNER_DDPQILUC |
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition estimation and unsymmetric reordering for diagonal dominance as preconditioner. | |
static const Type | INNER_MPTILUC |
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition estimation and maximum product transversal reordering as preconditioner. | |
static const Type | INNER_MLMPTILUC |
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition estimation and maximum product transversal reordering as preconditioner. | |
static const Type | INNER_MPTILU2 |
inner Solver based on BiCGStab(L) solver with second order ILU and maximum product transversal reordering as preconditioner. | |
static const Type | Trilinos_Aztec |
external Solver AztecOO from Trilinos package. | |
static const Type | Trilinos_Belos |
external Solver Belos from Trilinos package, currently without preconditioner. | |
static const Type | Trilinos_ML |
external Solver AztecOO with ML preconditioner. | |
static const Type | Trilinos_Ifpack |
external Solver AztecOO with Ifpack preconditioner. | |
static const Type | PETSc |
external Solver PETSc, More... | |
static const Type | ANI |
external Solver from ANI3D based on ILU2 (sequential Fortran version), More... | |
static const Type | FCBIILU2 |
external FCBIILU2 Solver (BIILU2 parallel F2C version). | |
static const Type | K3BIILU2 |
inner K3BIILU2 Solver (BIILU2 parallel version). | |
static const Type | SUPERLU |
external Solver SuperLU More... | |
Main class to set and solve linear system.
Solver class is used to set the coefficient Matrix, the right-hand side Vector and the initial guess Vector, construct the preconditioner and Solve the linear system.
Formally, Solver class is independent of INMOST::Mesh class.
Definition at line 29 of file inmost_solver.h.
INMOST::Solver::Solver | ( | std::string | solverName, |
std::string | prefix = "" , |
||
INMOST_MPI_Comm | _comm = INMOST_MPI_COMM_WORLD |
||
) |
Main constructor of the solver.
solverName | The solver name to be used for solution. |
prefix | The user specified name of the current solver. |
comm | Communicator for parallel data exchanges, MPI_COMM_WORLD by default. |
INMOST::Solver::Solver | ( | const Solver & | other | ) |
Copy a solver.
INMOST_DATA_REAL_TYPE INMOST::Solver::Condest | ( | INMOST_DATA_REAL_TYPE | tol, |
INMOST_DATA_ENUM_TYPE | maxits = 100 |
||
) |
Computes the smallest and the largest eigenvalue with the power method.
Requires SetMatrix to be called to compute the preconditioner. Currently works for internal methods only, since it uses internal matrix-vector multiplication. Largest eigenvalue: vprev = 0; v = rand(); while( |v|-|vprev| > tol ) {vprev = v; v = A*v; v /= |v|;} lambda_max = |v|; Smallest eigenvalue: vprev = 0; v = rand(); while( |v|-|vprev| > tol ){vprev = v; solve(A*v = v); v /= |v|;} lambda_min = 1.0/|v|; See answer by Blair Perot in: https://www.researchgate.net/post/What_is_the_best_way_to_estimate_the_condition_number_of_a_sparse_matrix.
tol | Tolerance used for power series. |
maxits | Maximum number of iterations allowed. |
|
static |
Finalize the stage of parallel solution.
If MPI was initialized in Solver::Initialize, then it will be finalized. By this reason, do not use any MPI function after call to this function.
|
static |
Return the list of all available solvers.
std::string INMOST::Solver::GetParameter | ( | std::string | name | ) | const |
Get the solver output parameter.
name | The name of solver's output parameter |
|
static |
Initialize the stage of parallel solution.
If MPI is not initialized yet, then it will be initialized.
database file is used to pass parameters to Inner solvers, PETSc and Trilinos packages. if database file for is provided any changes through SetParameter would not be effective for PETSc and Trilinos packages.
argc | The number of arguments transmitted to the function main. |
argv | The pointer to arguments transmitted to the function main. |
database | Usually the name of the file with the Solver parameters. |
The shortest call to this function with the default solver parameters is the following: Initialize(NULL,NULL,"");
Example of contents of the database file: Main: database.xml PETSc: petsc_options.txt Trilinos_Ifpack: trilinos_ifpack_options.xml Trilinos_ML: trilinos_ml_options.xml Trilinos_Aztec: trilinos_aztec_options.xml Trilinos_Belos: trilinos_belos_options.xml
bool INMOST::Solver::IsSolved | ( | ) | const |
Return the last solution successfullness.
|
static |
INMOST_DATA_ENUM_TYPE INMOST::Solver::Iterations | ( | ) | const |
Return the number of iterations performed by the last solution.
double INMOST::Solver::IterationsTime | ( | ) | const |
Return the time of iteration stage by the last solution.
Assign a solver.
double INMOST::Solver::PreconditionerTime | ( | ) | const |
Return the time of preconditioner evaluation by the last solution.
INMOST_DATA_REAL_TYPE INMOST::Solver::Residual | ( | ) | const |
Return the final precondioned residual achieved by the last solution.
const std::string INMOST::Solver::ReturnReason | ( | ) | const |
Get the reason of convergence or divergence of the last solution.
void INMOST::Solver::SetMatrix | ( | Sparse::Matrix & | A, |
bool | ModifiedPattern = true , |
||
bool | OldPreconditioner = false |
||
) |
Set the matrix and construct the preconditioner.
A | Matrix A in linear problem Ax = b |
ModifiedPattern | Indicates whether the structure of the matrix have changed since last call to Solver::SetMatrix. |
OldPreconditioner | If this parameter is set to true, then the previous preconditioner will be used, otherwise the new preconditioner will be constructed. |
Preconditioner will be constructed on call to this function
Any changes to preconditioner parameters should happen before that point. If you increase gmres_substep after this point, inner methods most likely will fail
void INMOST::Solver::SetParameter | ( | std::string | name, |
std::string | value | ||
) |
name | The name of parameter |
value | The value of parameter Set the solver parameter of the integer type. |
Parameters:
Parameters:
void INMOST::Solver::SetParameterEnum | ( | std::string | name, |
INMOST_DATA_ENUM_TYPE | value | ||
) |
Backward-compatibility access to integer-typed parameters.
Check Solver::SetParameter for available parameter names.
name | Name of the parameter. |
value | Value for the parameter. |
void INMOST::Solver::SetParameterReal | ( | std::string | name, |
INMOST_DATA_REAL_TYPE | value | ||
) |
Backward-compatibility access to real-typed parameters.
Check Solver::SetParameter for available parameter names.
name | Name of the parameter. |
value | Value for the parameter. |
const std::string INMOST::Solver::SolutionMetadataLine | ( | const std::string & | delimiter | ) | const |
Get the useful info of the last solution.
bool INMOST::Solver::Solve | ( | Sparse::Vector & | RHS, |
Sparse::Vector & | SOL | ||
) |
Solver the linear system: A*x = b.
Prior to this call you should call SetMatrix
RHS | The right-hand side Vector b. |
SOL | The initial guess to the solution on input and the solution Vector x on return. |
It is assumed that the coefficient matrix A have been set and the preconditioner have been already constructed.
std::string INMOST::Solver::SolverName | ( | ) | const |
Return the solver name.
std::string INMOST::Solver::SolverPrefix | ( | ) | const |
Return the solver user specified name of the current solver.
SolverVerbosityLevel INMOST::Solver::VerbosityLevel | ( | ) | const |
Return the solver verbosity specified level of the current solver.
|
static |
external Solver from ANI3D based on ILU2 (sequential Fortran version),
Definition at line 43 of file inmost_solver.h.
|
static |
external Solver PETSc,
Definition at line 42 of file inmost_solver.h.
|
static |
external Solver SuperLU
Definition at line 46 of file inmost_solver.h.