INMOST
Mathematical Modelling Toolkit
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
INMOST::Solver Class Reference

#include <inmost_solver.h>

Collaboration diagram for INMOST::Solver:

Classes

class  OrderInfo
 Base class for low level operations with objects of Solver class. More...
 

Public Types

enum  Type {
  INNER_ILU2, INNER_DDPQILUC, INNER_MPTILUC, INNER_MPTILU2,
  Trilinos_Aztec, Trilinos_Belos, Trilinos_ML, Trilinos_Ifpack,
  PETSc, ANI, FCBIILU2, K3BIILU2
}
 Type of the Solver can be currently used in this version of INMOST. More...
 

Public Member Functions

INMOST_DATA_REAL_TYPE GetConditionNumberL ()
 
INMOST_DATA_REAL_TYPE GetConditionNumberU ()
 
void SetParameterEnum (std::string name, INMOST_DATA_ENUM_TYPE value)
 
void SetParameterReal (std::string name, INMOST_DATA_REAL_TYPE value)
 
std::string GetName ()
 Get the used defined name of the Solver. More...
 
Type GetPackage () const
 Get the package Type. More...
 
void SetMatrix (Sparse::Matrix &A, bool OldPreconditioner=false)
 
bool Solve (Sparse::Vector &RHS, Sparse::Vector &SOL)
 
INMOST_DATA_ENUM_TYPE Iterations ()
 
INMOST_DATA_REAL_TYPE Residual ()
 
std::string GetReason ()
 
INMOST_DATA_REAL_TYPE Condest (INMOST_DATA_REAL_TYPE tol, INMOST_DATA_ENUM_TYPE maxits=100)
 
 Solver (Type pack, std::string _name="", INMOST_MPI_Comm comm=INMOST_MPI_COMM_WORLD)
 
 ~Solver ()
 
void Clear ()
 Clear all internal data of the current solver including matrix, preconditioner etc. More...
 

Static Public Member Functions

static std::string TypeName (Type t)
 
static void Initialize (int *argc, char ***argv, const char *database="")
 
static void Finalize ()
 
static bool isInitialized ()
 
static bool isFinalized ()
 

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 38 of file inmost_solver.h.

Member Enumeration Documentation

Type of the Solver can be currently used in this version of INMOST.

Enumerator
INNER_ILU2 

inner Solver based on BiCGStab(L) solver with second order ILU factorization as preconditioner.

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.

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.

INNER_MPTILU2 

inner Solver based on BiCGStab(L) solver with second order ILU and maximum product transversal reordering as preconditioner.

Trilinos_Aztec 

external Solver AztecOO from Trilinos package.

Trilinos_Belos 

external Solver Belos from Trilinos package, currently without preconditioner.

Trilinos_ML 

external Solver AztecOO with ML preconditioner.

Trilinos_Ifpack 

external Solver AztecOO with Ifpack preconditioner.

PETSc 

external Solver PETSc,

See Also
http://www.mcs.anl.gov/petsc/
ANI 

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

See Also
http://ani3d.sourceforge.net/
FCBIILU2 

external FCBIILU2 Solver (BIILU2 parallel F2C version).

K3BIILU2 

inner K3BIILU2 Solver (BIILU2 parallel version).

Definition at line 44 of file inmost_solver.h.

Constructor & Destructor Documentation

INMOST::Solver::Solver ( Type  pack,
std::string  _name = "",
INMOST_MPI_Comm  comm = INMOST_MPI_COMM_WORLD 
)

Main constructor of the solver. Solver name provided here is used to extract options from database file for PETSc and Trilinos packages.

Parameters
packThe package Type to be used for solution.
_nameThe 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
INMOST::Solver::~Solver ( )

Member Function Documentation

void INMOST::Solver::Clear ( )

Clear all internal data of the current solver including matrix, preconditioner etc.

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.
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
INMOST_DATA_REAL_TYPE INMOST::Solver::GetConditionNumberL ( )

Retrive approximate condition number produced by INNER_MPTILUC. The number is cond(L^-1).

INMOST_DATA_REAL_TYPE INMOST::Solver::GetConditionNumberU ( )

Retrive approximate condition number produced by INNER_MPTILUC. The number is cond(U^-1).

std::string INMOST::Solver::GetName ( )
inline

Get the used defined name of the Solver.

Definition at line 255 of file inmost_solver.h.

Type INMOST::Solver::GetPackage ( ) const
inline

Get the package Type.

Definition at line 258 of file inmost_solver.h.

std::string INMOST::Solver::GetReason ( )

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

See Also
Sparse::Solve
static void INMOST::Solver::Initialize ( int *  argc,
char ***  argv,
const char *  database = "" 
)
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 PETSc and Trilinos packages. if database file for is provided any changes through SetParameterEnum, SetParameterReal would not be effective for PETSc and Trilinos packages. Currently this database file provides directions for package-specific files. In future it is supposed to set up parameters for internal solvers.

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:

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  
static bool INMOST::Solver::isFinalized ( )
inlinestatic

Definition at line 349 of file inmost_solver.h.

static bool INMOST::Solver::isInitialized ( )
inlinestatic

Definition at line 348 of file inmost_solver.h.

INMOST_DATA_ENUM_TYPE INMOST::Solver::Iterations ( )

Return the number of iterations performed by the last solution.

See Also
Sparse::Solve
INMOST_DATA_REAL_TYPE INMOST::Solver::Residual ( )

Return the final residual achieved by the last solution.

See Also
Sparse::Solve
void INMOST::Solver::SetMatrix ( Sparse::Matrix A,
bool  OldPreconditioner = false 
)

Set the matrix and construct the preconditioner.

Parameters
AMatrix A in linear problem Ax = b
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

void INMOST::Solver::SetParameterEnum ( std::string  name,
INMOST_DATA_ENUM_TYPE  value 
)

Set the solver parameter of the integer type. You can find defaults for parameters in the top of the file inmost_solver.h.

Parameters:

  • "maximum_iterations" - total number of iterations
  • "schwartz_overlap" - number of overlapping levels for additive schwartz method, works for: INNER_ILU2, INNER_MLILUC 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_MLILUC
  • "reorder_nonzeros" - place sparser rows at the beggining of matrix during reordering, works for: INNER_MLILUC
  • "rescale_iterations" - number of iterations for two-side matrix rescaling, works for: INNER_ILU2, INNER_MLILUC
  • "condition_estimation" - exploit condition estimation of inversed factors to adapt drop and reuse tolerances, works for: INNER_MLILUC
  • "adapt_ddpq_tolerance" - adapt ddpq tolerance depending from the complexity of calculation of Schur complement, works for: INNER_MLILUC
void INMOST::Solver::SetParameterReal ( std::string  name,
INMOST_DATA_REAL_TYPE  value 
)

Set the solver parameter of the real type. You can find defaults for parameters in the top of the file inmost_solver.h.

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_MLILUC 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 then "drop_tolerance", typical value is drop_tolerance^2, works for: INNER_ILU2, INNER_MLILUC
  • "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_MLILUC
  • "fill_level" - level of fill for ILU-type preconditioners, works for: INNER_ILU2 (if LFILL is defined in solver_ilu2.hpp) Trilinos, Trilinos_Ifpack
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
static std::string INMOST::Solver::TypeName ( Type  t)
static

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