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

#include <PETScLinearSolver.h>

Public Member Functions

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

Private Attributes

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

Constructor & Destructor Documentation

◆ PETScLinearSolver()

MathLib::PETScLinearSolver::PETScLinearSolver ( const std::string  prefix,
BaseLib::ConfigTree const *const  option 
)

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.
optionPetsc options, which will be inserted into the global petsc options database.
Input File Parameter:
prj__linear_solvers__linear_solver__petsc
Input File Parameter:
prj__linear_solvers__linear_solver__petsc__parameters
Input File Parameter:
prj__linear_solvers__linear_solver__petsc__prefix

Definition at line 24 of file PETScLinearSolver.cpp.

26 {
27  // Insert options into petsc database. Default options are given in the
28  // string below.
29  std::string petsc_options =
30  "-ksp_type cg -pc_type bjacobi -ksp_rtol 1e-16 -ksp_max_it 10000";
31 
32  std::string prefix;
33 
34  if (option)
35  {
36  ignoreOtherLinearSolvers(*option, "petsc");
37 
39  if (auto const subtree = option->getConfigSubtreeOptional("petsc"))
40  {
41  if (auto const parameters =
43  subtree->getConfigParameterOptional<std::string>("parameters"))
44  {
45  petsc_options = *parameters;
46  }
47 
48  if (auto const pre =
50  subtree->getConfigParameterOptional<std::string>("prefix"))
51  {
52  if (!pre->empty())
53  prefix = *pre + "_";
54  }
55  }
56  }
57 #if PETSC_VERSION_LT(3, 7, 0)
58  PetscOptionsInsertString(petsc_options.c_str());
59 #else
60  PetscOptionsInsertString(nullptr, petsc_options.c_str());
61 #endif
62 
63  KSPCreate(PETSC_COMM_WORLD, &solver_);
64 
65  KSPGetPC(solver_, &pc_);
66 
67  if (!prefix.empty())
68  {
69  KSPSetOptionsPrefix(solver_, prefix.c_str());
70  }
71 
72  KSPSetInitialGuessNonzero(solver_, PETSC_TRUE);
73  KSPSetFromOptions(solver_); // set run-time options
74 }
PC pc_
Preconditioner type.
void ignoreOtherLinearSolvers(const BaseLib::ConfigTree &config, const std::string &solver_name)

References BaseLib::ConfigTree::getConfigSubtreeOptional(), MathLib::ignoreOtherLinearSolvers(), pc_, and solver_.

◆ ~PETScLinearSolver()

MathLib::PETScLinearSolver::~PETScLinearSolver ( )
inline

Definition at line 52 of file PETScLinearSolver.h.

52 { KSPDestroy(&solver_); }

References solver_.

Member Function Documentation

◆ getElapsedTime()

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

Get elapsed wall clock time.

Definition at line 65 of file PETScLinearSolver.h.

65 { 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 57 of file PETScLinearSolver.h.

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

References solver_.

◆ solve()

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

Definition at line 76 of file PETScLinearSolver.cpp.

77 {
78  BaseLib::RunTime wtimer;
79  wtimer.start();
80 
81 // define TEST_MEM_PETSC
82 #ifdef TEST_MEM_PETSC
83  PetscLogDouble mem1, mem2;
84  PetscMemoryGetCurrentUsage(&mem1);
85 #endif
86 
87  KSPSetOperators(solver_, A.getRawMatrix(), A.getRawMatrix());
88 
89  KSPSolve(solver_, b.getRawVector(), x.getRawVector());
90 
91  KSPConvergedReason reason;
92  KSPGetConvergedReason(solver_, &reason);
93 
94  bool converged = true;
95  if (reason > 0)
96  {
97  const char* ksp_type;
98  const char* pc_type;
99  KSPGetType(solver_, &ksp_type);
100  PCGetType(pc_, &pc_type);
101 
102  PetscPrintf(PETSC_COMM_WORLD,
103  "\n================================================");
104  PetscPrintf(PETSC_COMM_WORLD,
105  "\nLinear solver %s with %s preconditioner", ksp_type,
106  pc_type);
107 
108  PetscInt its;
109  KSPGetIterationNumber(solver_, &its);
110  PetscPrintf(PETSC_COMM_WORLD, "\nconverged in %d iterations", its);
111  switch (reason)
112  {
113  case KSP_CONVERGED_RTOL:
114  PetscPrintf(PETSC_COMM_WORLD,
115  " (relative convergence criterion fulfilled).");
116  break;
117  case KSP_CONVERGED_ATOL:
118  PetscPrintf(PETSC_COMM_WORLD,
119  " (absolute convergence criterion fulfilled).");
120  break;
121  default:
122  PetscPrintf(PETSC_COMM_WORLD, ".");
123  }
124 
125  PetscPrintf(PETSC_COMM_WORLD,
126  "\n================================================\n");
127  }
128  else if (reason == KSP_DIVERGED_ITS)
129  {
130  const char* ksp_type;
131  const char* pc_type;
132  KSPGetType(solver_, &ksp_type);
133  PCGetType(pc_, &pc_type);
134  PetscPrintf(PETSC_COMM_WORLD,
135  "\nLinear solver %s with %s preconditioner", ksp_type,
136  pc_type);
137  PetscPrintf(PETSC_COMM_WORLD,
138  "\nWarning: maximum number of iterations reached.\n");
139  }
140  else
141  {
142  converged = false;
143  if (reason == KSP_DIVERGED_INDEFINITE_PC)
144  {
145  PetscPrintf(PETSC_COMM_WORLD,
146  "\nDivergence because of indefinite preconditioner,");
147  PetscPrintf(PETSC_COMM_WORLD,
148  "\nTry to run again with "
149  "-pc_factor_shift_positive_definite option.\n");
150  }
151  else if (reason == KSP_DIVERGED_BREAKDOWN_BICG)
152  {
153  PetscPrintf(PETSC_COMM_WORLD,
154  "\nKSPBICG method was detected so the method could not "
155  "continue to enlarge the Krylov space.");
156  PetscPrintf(PETSC_COMM_WORLD,
157  "\nTry to run again with another solver.\n");
158  }
159  else if (reason == KSP_DIVERGED_NONSYMMETRIC)
160  {
161  PetscPrintf(PETSC_COMM_WORLD,
162  "\nMatrix or preconditioner is unsymmetric but KSP "
163  "requires symmetric.\n");
164  }
165  else
166  {
167  PetscPrintf(PETSC_COMM_WORLD,
168  "\nDivergence detected, use command option "
169  "-ksp_monitor or -log_summary to check the details.\n");
170  }
171  }
172 
173 #ifdef TEST_MEM_PETSC
174  PetscMemoryGetCurrentUsage(&mem2);
175  PetscPrintf(
176  PETSC_COMM_WORLD,
177  "###Memory usage by solver. Before: %f After: %f Increase: %d\n", mem1,
178  mem2, (int)(mem2 - mem1));
179 #endif
180 
181  elapsed_ctime_ += wtimer.elapsed();
182 
183  return converged;
184 }
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 71 of file PETScLinearSolver.h.

Referenced by getElapsedTime(), and solve().

◆ pc_

PC MathLib::PETScLinearSolver::pc_
private

Preconditioner type.

Definition at line 69 of file PETScLinearSolver.h.

Referenced by PETScLinearSolver(), and solve().

◆ solver_

KSP MathLib::PETScLinearSolver::solver_
private

Solver type.

Definition at line 68 of file PETScLinearSolver.h.

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


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