OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
petsc_solver< double > Class Template Reference

This class is able to do Matrix inversion in parallel with PETSC solvers. More...

Detailed Description

template<>
class petsc_solver< double >

This class is able to do Matrix inversion in parallel with PETSC solvers.

Example of use of this class

// velocity in the grid is the property 0, pressure is the property 1
constexpr int velocity = 0;
constexpr int pressure = 1;
// Domain, a rectangle
Box<2,float> domain({0.0,0.0},{3.0,1.0});
// Ghost (Not important in this case but required)
Ghost<2,float> g(0.01);
// Grid points on x=256 and y=64
long int sz[] = {256,64};
size_t szu[2];
szu[0] = (size_t)sz[0];
szu[1] = (size_t)sz[1];
// We need one more point on the left and down part of the domain
// This is given by the boundary conditions that we impose, the
// reason is mathematical in order to have a well defined system
// and cannot be discussed here
Padding<2> pd({1,1},{0,0});
// Distributed grid that store the solution
grid_dist_id<2,float,aggregate<float[2],float>,CartDecomposition<2,float>> g_dist(szu,domain,g);
// It is the maximum extension of the stencil
Ghost<2,long int> stencil_max(1);
// Finite difference scheme
FDScheme<lid_nn> fd(pd, stencil_max, domain,g_dist);
// Here we impose the equation, we start from the incompressibility Eq imposed in the bulk with the
// exception of the first point {0,0} and than we set P = 0 in {0,0}, why we are doing this is again
// mathematical to have a well defined system, an intuitive explanation is that P and P + c are both
// solution for the incompressibility equation, this produce an ill-posed problem to make it well posed
// we set one point in this case {0,0} the pressure to a fixed constant for convenience P = 0
fd.impose(ic_eq(),0.0, EQ_3, {0,0},{sz[0]-2,sz[1]-2},true);
fd.impose(Prs(), 0.0, EQ_3, {0,0},{0,0});
// Here we impose the Eq1 and Eq2
fd.impose(vx_eq(),0.0, EQ_1, {1,0},{sz[0]-2,sz[1]-2});
fd.impose(vy_eq(),0.0, EQ_2, {0,1},{sz[0]-2,sz[1]-2});
// v_x and v_y
// Imposing B1
fd.impose(v_x(),0.0, EQ_1, {0,0},{0,sz[1]-2});
fd.impose(avg_vy_f(),0.0, EQ_2 , {-1,0},{-1,sz[1]-1});
// Imposing B2
fd.impose(v_x(),0.0, EQ_1, {sz[0]-1,0},{sz[0]-1,sz[1]-2});
fd.impose(avg_vy(),1.0, EQ_2, {sz[0]-1,0},{sz[0]-1,sz[1]-1});
// Imposing B3
fd.impose(avg_vx_f(),0.0, EQ_1, {0,-1},{sz[0]-1,-1});
fd.impose(v_y(), 0.0, EQ_2, {0,0},{sz[0]-2,0});
// Imposing B4
fd.impose(avg_vx(),0.0, EQ_1, {0,sz[1]-1},{sz[0]-1,sz[1]-1});
fd.impose(v_y(), 0.0, EQ_2, {0,sz[1]-1},{sz[0]-2,sz[1]-1});
// When we pad the grid, there are points of the grid that are not
// touched by the previous condition. Mathematically this lead
// to have too many variables for the conditions that we are imposing.
// Here we are imposing variables that we do not touch to zero
//
// Padding pressure
fd.impose(Prs(), 0.0, EQ_3, {-1,-1},{sz[0]-1,-1});
fd.impose(Prs(), 0.0, EQ_3, {-1,sz[1]-1},{sz[0]-1,sz[1]-1});
fd.impose(Prs(), 0.0, EQ_3, {-1,0},{-1,sz[1]-2});
fd.impose(Prs(), 0.0, EQ_3, {sz[0]-1,0},{sz[0]-1,sz[1]-2});
// Impose v_x Padding Impose v_y padding
fd.impose(v_x(), 0.0, EQ_1, {-1,-1},{-1,sz[1]-1});
fd.impose(v_y(), 0.0, EQ_2, {-1,-1},{sz[0]-1,-1});
solver_type solver;
auto x = solver.solve(fd.getA(),fd.getB());

Definition at line 101 of file petsc_solver.hpp.

#include <petsc_solver.hpp>

Data Structures

struct  itError
 contain the infinity norm of the residual at each iteration More...
 
struct  solv_bench_info
 It contain the benchmark information for each solver. More...
 

Public Types

typedef Vector< double,
PETSC_BASE > 
return_type
 Type of the solution object.
 

Public Member Functions

void addTestSolver (std::string &solver)
 Add a test solver. More...
 
void removeTestSolver (const std::string &solver)
 Remove a test solver. More...
 
void log_monitor ()
 Set the Petsc solver. More...
 
void setSolver (KSPType type)
 Set the Petsc solver. More...
 
void setRelTol (PetscReal rtol_)
 Set the relative tolerance as stop criteria. More...
 
void setAbsTol (PetscReal abstol_)
 Set the absolute tolerance as stop criteria. More...
 
void setDivTol (PetscReal dtol_)
 Set the divergence tolerance. More...
 
void setMaxIter (PetscInt n)
 Set the maximum number of iteration for Krylov solvers. More...
 
void searchDirections (PetscInt l)
 
void setRestart (PetscInt n)
 For GMRES based method, the number of Krylov directions to orthogonalize against. More...
 
void setPreconditioner (PCType type)
 Set the preconditioner of the linear solver. More...
 
void setPreconditionerAMG_nl (int nl)
 Set the number of levels for the algebraic-multigrid preconditioner. More...
 
void setPreconditionerAMG_maxit (int nit)
 Set the maximum number of V or W cycle for algebraic-multi-grid. More...
 
void setPreconditionerAMG_relax (const std::string &type, int k=REL_ALL)
 Set the relaxation method for the algebraic-multi-grid preconditioner. More...
 
void setPreconditionerAMG_cycleType (const std::string &cycle_type, int sweep_up=-1, int sweep_dw=-1, int sweep_crs=-1)
 It set the type of cycle and optionally the number of sweep. More...
 
void setPreconditionerAMG_coarsen (const std::string &type)
 Set the coarsening method for the algebraic-multi-grid preconditioner. More...
 
void setPreconditionerAMG_interp (const std::string &type)
 Set the interpolation method for the algebraic-multi-grid preconditioner. More...
 
void setPreconditionerAMG_coarsenNodalType (int norm)
 Set the block coarsening norm type. More...
 
void setPreconditionerAMG_interp_eu_level (int k)
 Indicate the number of levels in the ILU(k) for the Euclid smoother. More...
 
void setBlockSize (int block_sz)
 Set how many degree of freedom each node has. More...
 
Vector< double, PETSC_BASE > solve (SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b, bool initial_guess=false)
 Here we invert the matrix and solve the system. More...
 
KSP getKSP ()
 Return the KSP solver. More...
 
solError get_residual_error (SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &x, const Vector< double, PETSC_BASE > &b)
 It return the resiual error. More...
 
solError get_residual_error (const Vector< double, PETSC_BASE > &x, const Vector< double, PETSC_BASE > &b)
 It return the resiual error. More...
 
bool solve (SparseMatrix< double, int, PETSC_BASE > &A, Vector< double, PETSC_BASE > &x, const Vector< double, PETSC_BASE > &b)
 Here we invert the matrix and solve the system. More...
 
Vector< double, PETSC_BASE > solve (const Vector< double, PETSC_BASE > &b)
 Here we invert the matrix and solve the system. More...
 
bool solve (Vector< double, PETSC_BASE > &x, const Vector< double, PETSC_BASE > &b)
 Here we invert the matrix and solve the system. More...
 
void setPetscOption (const char *name, const char *value)
 this function give you the possibility to set PETSC options More...
 
Vector< double, PETSC_BASE > try_solve (SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b)
 Try to solve the system using all the solvers and generate a report. More...
 

Private Member Functions

void write_bench_report ()
 Here we write the benchmark report. More...
 
void pre_solve_impl (const Mat &A_, const Vec &b_, Vec &x_)
 It set up the solver based on the provided options. More...
 
void print_progress_bar ()
 Print a progress bar on standard out. More...
 
void progress (PetscInt it)
 This function print an "*" showing the progress of the solvers. More...
 
void print_stat (solError &err)
 Print statistic about the solution error and method used. More...
 
std::string to_string_method (const std::string &solv)
 It convert the KSP type into a human read-able string. More...
 
void new_bench (const std::string &str)
 Allocate a new benchmark slot for a method. More...
 
void copy_if_better (double res, Vec &sol, double &best_res, Vec &best_sol)
 Copy the solution if better. More...
 
void try_solve_simple (const Mat &A_, const Vec &b_, Vec &x_)
 Try to solve the system x=inv(A)*b using all the Krylov solvers with simple Jacobi pre-conditioner. More...
 
void bench_solve_simple (const Mat &A_, const Vec &b_, Vec &x_, solv_bench_info &bench)
 Benchmark solve simple solving x=inv(A)*b. More...
 
void solve_simple (const Mat &A_, const Vec &b_, Vec &x_)
 solve simple use a Krylov solver + Simple selected Parallel Pre-conditioner More...
 
void solve_simple (const Vec &b_, Vec &x_)
 solve simple use a Krylov solver + Simple selected Parallel Pre-conditioner More...
 
void initKSP ()
 initialize the KSP object More...
 
void initKSPForTest ()
 initialize the KSP object for solver testing More...
 
void destroyKSP ()
 Destroy the KSP object. More...
 

Static Private Member Functions

static double calculate_it (double t, solv_bench_info &slv)
 Calculate the residual error at time t for one method. More...
 
static PetscErrorCode monitor_progress_residual (KSP ksp, PetscInt it, PetscReal res, void *data)
 procedure print the progress of the solver in benchmark mode More...
 
static solError statSolutionError (const Mat &A_, const Vec &b_, Vec &x_, KSP ksp)
 Calculate statistic on the error solution. More...
 
static solError getSolNormError (const Vec &b_, const Vec &x_, KSP ksp)
 Return the norm error of the solution. More...
 
static solError getSolNormError (const Mat &A_, const Vec &b_, const Vec &x_)
 Return the norm error of the solution. More...
 

Private Attributes

bool is_preconditioner_set = false
 indicate if the preconditioner is set
 
PetscInt maxits
 KSP Maximum number of iterations.
 
KSP ksp
 Main parallel solver.
 
size_t tmp
 Temporal variable used for calculation of static members.
 
openfpm::vector< std::string > solvs
 The full set of solvers.
 
openfpm::vector< solv_bench_info > bench
 It contain the solver benchmark results.
 
AMG_type atype = NONE_AMG
 Type of the algebraic multi-grid preconditioner.
 
int block_sz = 0
 Block size.
 

Member Function Documentation

void petsc_solver< double >::addTestSolver ( std::string &  solver)
inline

Add a test solver.

The try solve function use the most robust solvers in PETSC, if you want to add additionally solver like KSPIBCGS,KSPFBCGSR,KSPPGMRES, use addTestSolver(std::string(KSPIBCGS))

Parameters
solveradditional solver solver to test

Definition at line 834 of file petsc_solver.hpp.

void petsc_solver< double >::bench_solve_simple ( const Mat &  A_,
const Vec &  b_,
Vec &  x_,
solv_bench_info &  bench 
)
inlineprivate

Benchmark solve simple solving x=inv(A)*b.

Parameters
A_Matrix A
b_vector b
x_solution x
benchstructure that store the benchmark information

Definition at line 604 of file petsc_solver.hpp.

static double petsc_solver< double >::calculate_it ( double  t,
solv_bench_info &  slv 
)
inlinestaticprivate

Calculate the residual error at time t for one method.

Parameters
ttime
slvsolver
Returns
the residual

Definition at line 168 of file petsc_solver.hpp.

void petsc_solver< double >::copy_if_better ( double  res,
Vec &  sol,
double &  best_res,
Vec &  best_sol 
)
inlineprivate

Copy the solution if better.

Parameters
resResidual of the solution
solsolution
best_resresidual of the best solution
best_solbest solution

Definition at line 493 of file petsc_solver.hpp.

void petsc_solver< double >::destroyKSP ( )
inlineprivate

Destroy the KSP object.

Definition at line 726 of file petsc_solver.hpp.

solError petsc_solver< double >::get_residual_error ( SparseMatrix< double, int, PETSC_BASE > &  A,
const Vector< double, PETSC_BASE > &  x,
const Vector< double, PETSC_BASE > &  b 
)
inline

It return the resiual error.

Parameters
ASparse matrix
xsolution
bright-hand-side
Returns
the solution error norms

Definition at line 1312 of file petsc_solver.hpp.

solError petsc_solver< double >::get_residual_error ( const Vector< double, PETSC_BASE > &  x,
const Vector< double, PETSC_BASE > &  b 
)
inline

It return the resiual error.

Parameters
xsolution
bright-hand-side
Returns
the solution error norms

Definition at line 1325 of file petsc_solver.hpp.

KSP petsc_solver< double >::getKSP ( )
inline

Return the KSP solver.

In case you want to do fine tuning of the KSP solver before solve your system whith this function you can retrieve the KSP object

Returns
the Krylov solver

Definition at line 1298 of file petsc_solver.hpp.

static solError petsc_solver< double >::getSolNormError ( const Vec &  b_,
const Vec &  x_,
KSP  ksp 
)
inlinestaticprivate

Return the norm error of the solution.

Parameters
x_the solution
b_the right-hand-side
kspKrylov solver
Returns
the solution error

Definition at line 740 of file petsc_solver.hpp.

static solError petsc_solver< double >::getSolNormError ( const Mat &  A_,
const Vec &  b_,
const Vec &  x_ 
)
inlinestaticprivate

Return the norm error of the solution.

Parameters
A_the matrix that identity the linear system
x_the solution
b_the right-hand-side
Returns
the solution error

Definition at line 758 of file petsc_solver.hpp.

void petsc_solver< double >::initKSP ( )
inlineprivate

initialize the KSP object

Definition at line 706 of file petsc_solver.hpp.

void petsc_solver< double >::initKSPForTest ( )
inlineprivate

initialize the KSP object for solver testing

Definition at line 715 of file petsc_solver.hpp.

void petsc_solver< double >::log_monitor ( )
inline

Set the Petsc solver.

See Also
KSPType in PETSC manual for a list of all PETSC solvers

Definition at line 864 of file petsc_solver.hpp.

static PetscErrorCode petsc_solver< double >::monitor_progress_residual ( KSP  ksp,
PetscInt  it,
PetscReal  res,
void *  data 
)
inlinestaticprivate

procedure print the progress of the solver in benchmark mode

Parameters
kspSolver
itIteration number
resresudual
datacustom pointer to data
Returns
always zero

Definition at line 339 of file petsc_solver.hpp.

void petsc_solver< double >::new_bench ( const std::string &  str)
inlineprivate

Allocate a new benchmark slot for a method.

Parameters
strMethod name

Definition at line 477 of file petsc_solver.hpp.

void petsc_solver< double >::pre_solve_impl ( const Mat &  A_,
const Vec &  b_,
Vec &  x_ 
)
inlineprivate

It set up the solver based on the provided options.

Parameters
A_Matrix
x_solution
b_right-hand-side

Definition at line 300 of file petsc_solver.hpp.

void petsc_solver< double >::print_progress_bar ( )
inlineprivate

Print a progress bar on standard out.

Definition at line 320 of file petsc_solver.hpp.

void petsc_solver< double >::print_stat ( solError err)
inlineprivate

Print statistic about the solution error and method used.

Parameters
errstructure that contain the solution errors

Definition at line 405 of file petsc_solver.hpp.

void petsc_solver< double >::progress ( PetscInt  it)
inlineprivate

This function print an "*" showing the progress of the solvers.

Parameters
ititeration number

Definition at line 380 of file petsc_solver.hpp.

void petsc_solver< double >::removeTestSolver ( const std::string &  solver)
inline

Remove a test solver.

The try solve function use the most robust solvers in PETSC, if you want to remove a solver like use removeTestSolver(std::string(KSPIBCGS))

Parameters
solverremove solver to test

Definition at line 847 of file petsc_solver.hpp.

void petsc_solver< double >::searchDirections ( PetscInt  l)
inline

For the BiCGStab(L) it define the number of search directions

BiCG Methods can fail for base break-down (or near break-down). BiCGStab(L) try to avoid this problem. Such method has a parameter L (2 by default) Bigger is L more the system will try to avoid the breakdown

Parameters
lIncreasing L should reduce the probability of failure of the solver because of break-down of the base

Definition at line 963 of file petsc_solver.hpp.

void petsc_solver< double >::setAbsTol ( PetscReal  abstol_)
inline

Set the absolute tolerance as stop criteria.

See Also
PETSC manual KSPSetTolerances for an explanation
Parameters
abstol_Absolute tolerance

Definition at line 907 of file petsc_solver.hpp.

void petsc_solver< double >::setBlockSize ( int  block_sz)
inline

Set how many degree of freedom each node has.

In case you are solving a system of equations this function, help in setting the degree of freedom each grid point has. Setting this parameter to the number of variables for each grid point it should improve the convergenve of the solvers in particular using algebraic-multi-grid

Parameters
block_sznumber of degree of freedom

Definition at line 1243 of file petsc_solver.hpp.

void petsc_solver< double >::setDivTol ( PetscReal  dtol_)
inline

Set the divergence tolerance.

See Also
PETSC manual KSPSetTolerances for an explanation
Parameters
dtol_tolerance

Definition at line 925 of file petsc_solver.hpp.

void petsc_solver< double >::setMaxIter ( PetscInt  n)
inline

Set the maximum number of iteration for Krylov solvers.

Parameters
nmaximum number of iterations

Definition at line 941 of file petsc_solver.hpp.

void petsc_solver< double >::setPetscOption ( const char *  name,
const char *  value 
)
inline

this function give you the possibility to set PETSC options

this function call PetscOptionsSetValue

Parameters
namethe name of the option
valuethe value of the option

Definition at line 1415 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditioner ( PCType  type)
inline

Set the preconditioner of the linear solver.

The preconditoner that can be set are the same as PETSC

For a full list please visit Please visit: http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCType.html#PCType

An exception is PCHYPRE with BOOMERAMG in this case use PCHYPRE_BOOMERAMG. Many preconditioners has default values, but many times the default values are not good. Here we list some interesting case

Algebraic-multi-grid based preconditioners

Parameters for this type of preconditioner can be set using setPreconditionerAMG_* functions.

  • Number of levels set by setPreconditionerAMG_nl
  • Maximum number of cycles setPreconditionerAMG_maxit
  • Smooth or relax method setPreconditionerAMG_smooth,setPreconditionerAMG_relax
  • Cycle type and number of sweep (relaxation steps) when going up or down setPreconditionerAMG_cycleType
  • Coarsening method setPreconditionerAMG_coarsen
  • interpolation schemes setPreconditionerAMG_interp
Parameters
typeof the preconditioner

Definition at line 1005 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_coarsen ( const std::string &  type)
inline

Set the coarsening method for the algebraic-multi-grid preconditioner.

Possible values can be "CLJP","Ruge-Stueben","modifiedRuge-Stueben","Falgout", "PMIS", "HMIS"

Warning
in case of big problem use PMIS or HMIS otherwise the AMG can hang (take a lot of time) in setup
Parameters
typeof the preconditioner smoothing operator

Definition at line 1150 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_coarsenNodalType ( int  norm)
inline

Set the block coarsening norm type.

The use of this function make sanse if you specify the degree of freedom for each node using setBlockSize

  • 0 (default each variable treat independently)
  • 1 Frobenius norm
  • 2 sum of absolute values of elements in each block
  • 3 largest element in each block (not absolute value)
  • 4 row-sum norm
  • 6 sum of all values in each block

In case the matrix represent a system of equations in general

Parameters
normtype

Definition at line 1201 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_cycleType ( const std::string &  cycle_type,
int  sweep_up = -1,
int  sweep_dw = -1,
int  sweep_crs = -1 
)
inline

It set the type of cycle and optionally the number of sweep.

This function set the cycle type for the multigrid methods. Possible values are:

  • V cycle
  • W cycle

Optionally you can set the number of sweep or relaxation steps on each grid when going up and when going down

Parameters
cycle_typecycle type
sweep_up
sweep_dw
sweep_crsspeep at the coarse level

Definition at line 1119 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_interp ( const std::string &  type)
inline

Set the interpolation method for the algebraic-multi-grid preconditioner.

Possible values can be "classical", "direct", "multipass", "multipass-wts", "ext+i", "ext+i-cc", "standard", "standard-wts", "FF", "FF1"

Parameters
typeof the interpolation scheme

Definition at line 1172 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_interp_eu_level ( int  k)
inline

Indicate the number of levels in the ILU(k) for the Euclid smoother.

Parameters
knumber of levels for the Euclid smoother

Definition at line 1218 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_maxit ( int  nit)
inline

Set the maximum number of V or W cycle for algebraic-multi-grid.

Parameters
nitnumber of levels

Definition at line 1054 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_nl ( int  nl)
inline

Set the number of levels for the algebraic-multigrid preconditioner.

In case you select an algebraic preconditioner like PCHYPRE or PCGAMG you can set the number of levels using this function

Parameters
nlnumber of levels

Definition at line 1037 of file petsc_solver.hpp.

void petsc_solver< double >::setPreconditionerAMG_relax ( const std::string &  type,
int  k = REL_ALL 
)
inline

Set the relaxation method for the algebraic-multi-grid preconditioner.

Possible values for relazation can be "Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel", "SOR/Jacobi","backward-SOR/Jacobi",hybrid chaotic Gauss-Seidel (works only with OpenMP), "symmetric-SOR/Jacobi","l1scaled-SOR/Jacobi","Gaussian-elimination","CG","Chebyshev", "FCF-Jacobi","l1scaled-Jacobi"

Every smooth operator can have additional parameters to be set in order to correctly work.

Parameters
typeof relax method
kwhere is applied REL_UP,REL_DOWN,REL_ALL (default is all)

Definition at line 1086 of file petsc_solver.hpp.

void petsc_solver< double >::setRelTol ( PetscReal  rtol_)
inline

Set the relative tolerance as stop criteria.

See Also
PETSC manual KSPSetTolerances for an explanation
Parameters
rtol_Relative tolerance

Definition at line 889 of file petsc_solver.hpp.

void petsc_solver< double >::setRestart ( PetscInt  n)
inline

For GMRES based method, the number of Krylov directions to orthogonalize against.

Parameters
nnumber of directions

Definition at line 973 of file petsc_solver.hpp.

void petsc_solver< double >::setSolver ( KSPType  type)
inline

Set the Petsc solver.

See Also
KSPType in PETSC manual for a list of all PETSC solvers
Parameters
typepetsc solver type

Definition at line 877 of file petsc_solver.hpp.

Vector<double,PETSC_BASE> petsc_solver< double >::solve ( SparseMatrix< double, int, PETSC_BASE > &  A,
const Vector< double, PETSC_BASE > &  b,
bool  initial_guess = false 
)
inline

Here we invert the matrix and solve the system.

Warning
umfpack is not a parallel solver, this function work only with one processor
Note
if you want to use umfpack in a NON parallel, but on a distributed data, use solve with triplet
Template Parameters
implImplementation of the SparseMatrix
Parameters
Asparse matrix
bvector
initial_guesstrue if x has the initial guess
Returns
the solution

Definition at line 1263 of file petsc_solver.hpp.

bool petsc_solver< double >::solve ( SparseMatrix< double, int, PETSC_BASE > &  A,
Vector< double, PETSC_BASE > &  x,
const Vector< double, PETSC_BASE > &  b 
)
inline

Here we invert the matrix and solve the system.

Parameters
Asparse matrix
bvector
xsolution and initial guess
Returns
true if succeed

Definition at line 1339 of file petsc_solver.hpp.

Vector<double,PETSC_BASE> petsc_solver< double >::solve ( const Vector< double, PETSC_BASE > &  b)
inline

Here we invert the matrix and solve the system.

Parameters
bvector
Returns
true if succeed

Definition at line 1361 of file petsc_solver.hpp.

bool petsc_solver< double >::solve ( Vector< double, PETSC_BASE > &  x,
const Vector< double, PETSC_BASE > &  b 
)
inline

Here we invert the matrix and solve the system.

In this call we are interested in solving the system with multiple right-hand-side and the same Matrix. We do not set the Matrix again and this give us the possibility to-skip the preconditioning setting that in some case like Algebraic-multi-grid can be expensive

Parameters
bvector
xsolution and initial guess
Returns
true if succeed

Definition at line 1396 of file petsc_solver.hpp.

void petsc_solver< double >::solve_simple ( const Mat &  A_,
const Vec &  b_,
Vec &  x_ 
)
inlineprivate

solve simple use a Krylov solver + Simple selected Parallel Pre-conditioner

Parameters
A_SparseMatrix
b_right-hand-side
x_solution

Definition at line 641 of file petsc_solver.hpp.

void petsc_solver< double >::solve_simple ( const Vec &  b_,
Vec &  x_ 
)
inlineprivate

solve simple use a Krylov solver + Simple selected Parallel Pre-conditioner

Parameters
b_right hand side
x_solution

Definition at line 672 of file petsc_solver.hpp.

static solError petsc_solver< double >::statSolutionError ( const Mat &  A_,
const Vec &  b_,
Vec &  x_,
KSP  ksp 
)
inlinestaticprivate

Calculate statistic on the error solution.

Parameters
A_Matrix of the system
b_Right hand side of the matrix
x_Solution
kspKrylov solver
Returns
the solution error

Definition at line 688 of file petsc_solver.hpp.

std::string petsc_solver< double >::to_string_method ( const std::string &  solv)
inlineprivate

It convert the KSP type into a human read-able string.

Parameters
solvsolver (short form)
Returns
the name of the solver in long form

Definition at line 426 of file petsc_solver.hpp.

Vector<double,PETSC_BASE> petsc_solver< double >::try_solve ( SparseMatrix< double, int, PETSC_BASE > &  A,
const Vector< double, PETSC_BASE > &  b 
)
inline

Try to solve the system using all the solvers and generate a report.

In this mode the system will try different Solvers, Preconditioner and combination of solvers in order to find the best solver in speed, and precision. As output it will produce a performance report

Parameters
AMatrix to invert
bright hand side
Returns
the solution

Definition at line 1431 of file petsc_solver.hpp.

void petsc_solver< double >::try_solve_simple ( const Mat &  A_,
const Vec &  b_,
Vec &  x_ 
)
inlineprivate

Try to solve the system x=inv(A)*b using all the Krylov solvers with simple Jacobi pre-conditioner.

It try to solve the system using JACOBI pre-conditioner and all the Krylov solvers available at the end it write a report

Parameters
A_Matrix
b_vector of coefficents
x_solution

Definition at line 512 of file petsc_solver.hpp.

void petsc_solver< double >::write_bench_report ( )
inlineprivate

Here we write the benchmark report.

Definition at line 185 of file petsc_solver.hpp.


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