OGS
MathLib::PETScLinearSolver Class Reference

Detailed Description

A class of linear solver based on PETSc routines.

All command-line options that are not recognized by OGS are passed on to PETSc, i.e., they potentially affect the linear solver. The linear solver options in the project file take precedence over the command-line options, because the former are processed at a later stage.

Definition at line 35 of file PETScLinearSolver.h.

#include <PETScLinearSolver.h>

Public Member Functions

 PETScLinearSolver (std::string const &prefix, std::string const &petsc_options)
 
 ~PETScLinearSolver ()
 
bool solve (PETScMatrix &A, PETScVector &b, PETScVector &x)
 
PetscInt getNumberOfIterations () const
 Get number of iterations.
 
double getElapsedTime () const
 Get elapsed wall clock time.
 

Private Attributes

KSP solver_
 Solver type.
 
PC pc_
 Preconditioner type.
 
double elapsed_ctime_ = 0.0
 Clock time.
 

Constructor & Destructor Documentation

◆ PETScLinearSolver()

MathLib::PETScLinearSolver::PETScLinearSolver ( std::string const & prefix,
std::string const & petsc_options )

Constructor

Parameters
prefixName used to distinguish the options in the command line for this solver. It can be the name of the PDE that owns an instance of this class.
petsc_optionsPETSc options string which is passed to PETSc lib and inserted in the PETSc option database (see https://petsc.org/release/manualpages/Sys/PetscOptionsInsertString/).

Definition at line 23 of file PETScLinearSolver.cpp.

25{
26#if PETSC_VERSION_LT(3, 7, 0)
27 PetscOptionsInsertString(petsc_options.c_str());
28#else
29 PetscOptionsInsertString(nullptr, petsc_options.c_str());
30#endif
31
32 KSPCreate(PETSC_COMM_WORLD, &solver_);
33
34 KSPGetPC(solver_, &pc_);
35
36 if (!prefix.empty())
37 {
38 KSPSetOptionsPrefix(solver_, prefix.c_str());
39 }
40
41 KSPSetInitialGuessNonzero(solver_, PETSC_TRUE);
42 KSPSetFromOptions(solver_); // set run-time options
43}
PC pc_
Preconditioner type.

References pc_, and solver_.

◆ ~PETScLinearSolver()

MathLib::PETScLinearSolver::~PETScLinearSolver ( )
inline

Definition at line 50 of file PETScLinearSolver.h.

50{ KSPDestroy(&solver_); }

References solver_.

Member Function Documentation

◆ getElapsedTime()

double MathLib::PETScLinearSolver::getElapsedTime ( ) const
inline

Get elapsed wall clock time.

Definition at line 63 of file PETScLinearSolver.h.

63{ return elapsed_ctime_; }
double elapsed_ctime_
Clock time.

References elapsed_ctime_.

◆ getNumberOfIterations()

PetscInt MathLib::PETScLinearSolver::getNumberOfIterations ( ) const
inline

Get number of iterations.

Definition at line 55 of file PETScLinearSolver.h.

56 {
57 PetscInt its = 0;
58 KSPGetIterationNumber(solver_, &its);
59 return its;
60 }

References solver_.

◆ solve()

bool MathLib::PETScLinearSolver::solve ( PETScMatrix & A,
PETScVector & b,
PETScVector & x )

Definition at line 45 of file PETScLinearSolver.cpp.

46{
47 BaseLib::RunTime wtimer;
48 wtimer.start();
49
50// define TEST_MEM_PETSC
51#ifdef TEST_MEM_PETSC
52 PetscLogDouble mem1, mem2;
53 PetscMemoryGetCurrentUsage(&mem1);
54#endif
55
56 KSPNormType norm_type;
57 KSPGetNormType(solver_, &norm_type);
58 const char* ksp_type;
59 const char* pc_type;
60 KSPGetType(solver_, &ksp_type);
61 PCGetType(pc_, &pc_type);
62
63 PetscPrintf(PETSC_COMM_WORLD,
64 "\n================================================");
65 PetscPrintf(PETSC_COMM_WORLD,
66 "\nLinear solver %s with %s preconditioner using %s", ksp_type,
67 pc_type, KSPNormTypes[norm_type]);
68
69 KSPSetOperators(solver_, A.getRawMatrix(), A.getRawMatrix());
70
71 KSPSolve(solver_, b.getRawVector(), x.getRawVector());
72
73 KSPConvergedReason reason;
74 KSPGetConvergedReason(solver_, &reason);
75
76 bool converged = true;
77 if (reason > 0)
78 {
79 PetscInt its;
80 KSPGetIterationNumber(solver_, &its);
81 PetscPrintf(PETSC_COMM_WORLD, "\nconverged in %d iterations", its);
82 switch (reason)
83 {
84 case KSP_CONVERGED_RTOL:
85 PetscPrintf(PETSC_COMM_WORLD,
86 " (relative convergence criterion fulfilled).");
87 break;
88 case KSP_CONVERGED_ATOL:
89 PetscPrintf(PETSC_COMM_WORLD,
90 " (absolute convergence criterion fulfilled).");
91 break;
92 default:
93 PetscPrintf(PETSC_COMM_WORLD, ".");
94 }
95
96 PetscPrintf(PETSC_COMM_WORLD,
97 "\n================================================\n");
98 }
99 else if (reason == KSP_DIVERGED_ITS)
100 {
101 PetscPrintf(PETSC_COMM_WORLD,
102 "\nWarning: maximum number of iterations reached.\n");
103 }
104 else
105 {
106 converged = false;
107 if (reason == KSP_DIVERGED_INDEFINITE_PC)
108 {
109 PetscPrintf(PETSC_COMM_WORLD,
110 "\nDivergence because of indefinite preconditioner,");
111 PetscPrintf(PETSC_COMM_WORLD,
112 "\nTry to run again with "
113 "-pc_factor_shift_positive_definite option.\n");
114 }
115 else if (reason == KSP_DIVERGED_BREAKDOWN_BICG)
116 {
117 PetscPrintf(PETSC_COMM_WORLD,
118 "\nKSPBICG method was detected so the method could not "
119 "continue to enlarge the Krylov space.");
120 PetscPrintf(PETSC_COMM_WORLD,
121 "\nTry to run again with another solver.\n");
122 }
123 else if (reason == KSP_DIVERGED_NONSYMMETRIC)
124 {
125 PetscPrintf(PETSC_COMM_WORLD,
126 "\nMatrix or preconditioner is unsymmetric but KSP "
127 "requires symmetric.\n");
128 }
129 else
130 {
131 PetscPrintf(PETSC_COMM_WORLD,
132 "\nDivergence detected, use command option "
133 "-ksp_monitor or -log_summary to check the details.\n");
134 }
135 }
136
137#ifdef TEST_MEM_PETSC
138 PetscMemoryGetCurrentUsage(&mem2);
139 PetscPrintf(
140 PETSC_COMM_WORLD,
141 "###Memory usage by solver. Before: %f After: %f Increase: %d\n", mem1,
142 mem2, (int)(mem2 - mem1));
143#endif
144
145 elapsed_ctime_ += wtimer.elapsed();
146
147 return converged;
148}
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32

References BaseLib::RunTime::elapsed(), elapsed_ctime_, MathLib::PETScMatrix::getRawMatrix(), MathLib::PETScVector::getRawVector(), pc_, solver_, and BaseLib::RunTime::start().

Member Data Documentation

◆ elapsed_ctime_

double MathLib::PETScLinearSolver::elapsed_ctime_ = 0.0
private

Clock time.

Definition at line 69 of file PETScLinearSolver.h.

Referenced by getElapsedTime(), and solve().

◆ pc_

PC MathLib::PETScLinearSolver::pc_
private

Preconditioner type.

Definition at line 67 of file PETScLinearSolver.h.

Referenced by PETScLinearSolver(), and solve().

◆ solver_

KSP MathLib::PETScLinearSolver::solver_
private

Solver type.

Definition at line 66 of file PETScLinearSolver.h.

Referenced by PETScLinearSolver(), ~PETScLinearSolver(), getNumberOfIterations(), and solve().


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