INMOST
A toolkit for distributed mathematical modeling
INMOST::Solver Class Reference

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...
 
Solveroperator= (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...
 

Detailed Description

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.

See also
Sparse::Matrix
Sparse::Vector
Sparse::Solve

Definition at line 29 of file inmost_solver.h.

Constructor & Destructor Documentation

◆ Solver() [1/2]

INMOST::Solver::Solver ( std::string  solverName,
std::string  prefix = "",
INMOST_MPI_Comm  _comm = INMOST_MPI_COMM_WORLD 
)

Main constructor of the solver.

Parameters
solverNameThe solver name to be used for solution.
prefixThe user specified name of the current solver.
commCommunicator for parallel data exchanges, MPI_COMM_WORLD by default.
See also
Solver::Initialize
Solver::SetMatrix
Solver::Solve
Solver::Finalize

◆ Solver() [2/2]

INMOST::Solver::Solver ( const Solver other)

Copy a solver.

Warning
Not all solvers support assignment. This operation may be very expensive.

Member Function Documentation

◆ Condest()

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.

Parameters
tolTolerance used for power series.
maxitsMaximum number of iterations allowed.
Returns
Condition number or 1.0e100 if not converged.

◆ Finalize()

static void INMOST::Solver::Finalize ( )
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.

See also
Solver::Initialize
Solver::isFinalized

◆ getAvailableSolvers()

static std::vector<std::string> INMOST::Solver::getAvailableSolvers ( )
static

Return the list of all available solvers.

See also
Solver::isSolverAvailable

◆ GetParameter()

std::string INMOST::Solver::GetParameter ( std::string  name) const

Get the solver output parameter.

Parameters
nameThe name of solver's output parameter
See also
Solver::SetParameter

◆ Initialize()

static void INMOST::Solver::Initialize ( int *  argc,
char ***  argv,
const char *  database = NULL 
)
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.

Parameters
argcThe number of arguments transmitted to the function main.
argvThe pointer to arguments transmitted to the function main.
databaseUsually 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,"");

See also
Solver::Finalize
Solver::isInitialized

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

◆ IsSolved()

bool INMOST::Solver::IsSolved ( ) const

Return the last solution successfullness.

See also
Sparse::Solve

◆ isSolverAvailable()

static bool INMOST::Solver::isSolverAvailable ( std::string  name)
static

Checks if solver available.

Parameters
nameSolver name
See also
Solver::getAvailableSolvers

◆ Iterations()

INMOST_DATA_ENUM_TYPE INMOST::Solver::Iterations ( ) const

Return the number of iterations performed by the last solution.

See also
Sparse::Solve

◆ IterationsTime()

double INMOST::Solver::IterationsTime ( ) const

Return the time of iteration stage by the last solution.

See also
Sparse::Solve

◆ operator=()

Solver& INMOST::Solver::operator= ( const Solver other)

Assign a solver.

Warning
Not all solvers support assignment. This operation may be very expensive.

◆ PreconditionerTime()

double INMOST::Solver::PreconditionerTime ( ) const

Return the time of preconditioner evaluation by the last solution.

See also
Sparse::Solve

◆ Residual()

INMOST_DATA_REAL_TYPE INMOST::Solver::Residual ( ) const

Return the final precondioned residual achieved by the last solution.

See also
Sparse::Solve

◆ ReturnReason()

const std::string INMOST::Solver::ReturnReason ( ) const

Get the reason of convergence or divergence of the last solution.

See also
Sparse::Solve

◆ SetMatrix()

void INMOST::Solver::SetMatrix ( Sparse::Matrix A,
bool  ModifiedPattern = true,
bool  OldPreconditioner = false 
)

Set the matrix and construct the preconditioner.

Parameters
AMatrix A in linear problem Ax = b
ModifiedPatternIndicates whether the structure of the matrix have changed since last call to Solver::SetMatrix.
OldPreconditionerIf 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

  • for INNER_*, PETSc and ANI packages
  • for Trilinos preconditioner will be constructed each time Sparse::Solve is called

Any changes to preconditioner parameters should happen before that point. If you increase gmres_substep after this point, inner methods most likely will fail

◆ SetParameter()

void INMOST::Solver::SetParameter ( std::string  name,
std::string  value 
)
Parameters
nameThe name of parameter
valueThe value of parameter Set the solver parameter of the integer type.

Parameters:

  • "maximum_iterations" - total number of iterations
  • "schwartz_overlap" - number of overlapping levels for additive schwartz method, works for: INNER_ILU2, INNER_MPTILU2, INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC Trilinos_Aztec, Trilinos_Belos, Trilinos_ML, Trilinos_Ifpack PETSc
  • "gmres_substeps" - number of gmres steps performed after each bicgstab step, works for: INNER_ILU2, INNER_MPTILU2, INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC
  • "reorder_nonzeros" - place sparser rows at the beggining of matrix during reordering, works for: INNER_DDPQILUC
  • "rescale_iterations" - number of iterations for two-side matrix rescaling, works for: INNER_ILU2, INNER_MPTILU2, INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC
  • "condition_estimation" - exploit condition estimation of inversed factors to adapt drop and reuse tolerances, works for: INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC
  • "adapt_ddpq_tolerance" - adapt ddpq tolerance depending from the complexity of calculation of Schur complement, works for: INNER_DDPQILUC Set the solver parameter of the real type.

Parameters:

  • "absolute_tolerance" - iterative method will stop on i-th iteration if ||A x(i)-b|| < absolute_tolerance
  • "relative_tolerance" - iterative method will stop on i-th iteration if ||A x(i)-b||/||A x(0) - b||
  • "divergence_tolerance" - iterative method will fail if ||A x(i) - b|| > divergence_tolerance
  • "drop_tolerance" - tolerance for dropping values during incomplete factorization, works for: INNER_ILU2, INNER_MPTILU2, INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC Trilinos_Aztec, Trilinos_Ifpack PETSc
  • "reuse_tolerance" - tolerance for reusing values during incomplete factorization, these values are used only during calculation of L and U factors and/or Schur complement and discarded once factorization is done, value should be less than "drop_tolerance", typical value is drop_tolerance^2, works for: INNER_ILU2, INNER_MPTILU2, INNER_MPTILUC, INNER_MLMPTILUC, INNER_DDPQILUC
  • "ddpq_tolerance" - by this tolerance most diagonnaly-dominant elements will be selected to form the next level of factorization, the closer the tolerance is to one the smaller will be the level. Actual rule is: A(i,j)/(sum(A(i,:))+sum(A(:,j))-A(i,j)) > ddpq_tolerance * A(imax,jmax)/(sum(A(imax,:))+sum(A(:,jmax))-A(imax,jmax)) where on imax, jmax maximum is reached. works for: INNER_DDPQILUC
  • "fill_level" - level of fill for ILU-type preconditioners, works for: INNER_ILU2 (if LFILL is defined in solver_ilu2.hpp) Trilinos, Trilinos_Ifpack
  • "pivot_condition" - delay factoriation of the row of the matrix to the next level, if the estimated condition number is above prescribed value. works for: INNER_MLMPTILUC
  • "pivot_diag" - delay factoriation of the row of the matrix to the next level, if the inverted diagonal value is above the prescribed value. works for: INNER_MLMPTILUC
    See also
    Solver::GetParameter

◆ SetParameterEnum()

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.

Parameters
nameName of the parameter.
valueValue for the parameter.
See also
Solver::SetParameter

◆ SetParameterReal()

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.

Parameters
nameName of the parameter.
valueValue for the parameter.
See also
Solver::SetParameter

◆ SolutionMetadataLine()

const std::string INMOST::Solver::SolutionMetadataLine ( const std::string &  delimiter) const

Get the useful info of the last solution.

See also
Sparse::Solve

◆ Solve()

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

Parameters
RHSThe right-hand side Vector b.
SOLThe 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.

See also
Sparse::SetMatrix

◆ SolverName()

std::string INMOST::Solver::SolverName ( ) const

Return the solver name.

See also
Sparse::Solve

◆ SolverPrefix()

std::string INMOST::Solver::SolverPrefix ( ) const

Return the solver user specified name of the current solver.

See also
Sparse::Solve

◆ VerbosityLevel()

SolverVerbosityLevel INMOST::Solver::VerbosityLevel ( ) const

Return the solver verbosity specified level of the current solver.

See also
Sparse::Solve

Member Data Documentation

◆ ANI

const Type INMOST::Solver::ANI
static

external Solver from ANI3D based on ILU2 (sequential Fortran version),

See also
http://ani3d.sourceforge.net/

Definition at line 43 of file inmost_solver.h.

◆ PETSc

const Type INMOST::Solver::PETSc
static

external Solver PETSc,

See also
http://www.mcs.anl.gov/petsc/

Definition at line 42 of file inmost_solver.h.

◆ SUPERLU

const Type INMOST::Solver::SUPERLU
static

external Solver SuperLU

See also
https://github.com/starseeker/SuperLU

Definition at line 46 of file inmost_solver.h.


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