OGS
PETScLinearSolver.cpp
Go to the documentation of this file.
1
17#include "PETScLinearSolver.h"
18
19#include "BaseLib/RunTime.h"
20
21namespace MathLib
22{
23PETScLinearSolver::PETScLinearSolver(std::string const& prefix,
24 std::string const& petsc_options)
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
44 KSPType ksp_type;
45 KSPGetType(solver_, &ksp_type);
47}
48
50{
51 BaseLib::RunTime wtimer;
52 wtimer.start();
53
54// define TEST_MEM_PETSC
55#ifdef TEST_MEM_PETSC
56 PetscLogDouble mem1, mem2;
57 PetscMemoryGetCurrentUsage(&mem1);
58#endif
59
60 KSPNormType norm_type;
61 KSPGetNormType(solver_, &norm_type);
62 const char* ksp_type;
63 const char* pc_type;
64 KSPGetType(solver_, &ksp_type);
65 PCGetType(pc_, &pc_type);
66
67 PetscPrintf(PETSC_COMM_WORLD,
68 "\n================================================");
69 PetscPrintf(PETSC_COMM_WORLD,
70 "\nLinear solver %s with %s preconditioner using %s", ksp_type,
71 pc_type, KSPNormTypes[norm_type]);
72
73 KSPSetOperators(solver_, A.getRawMatrix(), A.getRawMatrix());
74
75 KSPSolve(solver_, b.getRawVector(), x.getRawVector());
76
77 KSPConvergedReason reason;
78 KSPGetConvergedReason(solver_, &reason);
79
80 bool converged = true;
81 if (reason > 0)
82 {
83 PetscInt its;
84 KSPGetIterationNumber(solver_, &its);
85 PetscPrintf(PETSC_COMM_WORLD, "\nconverged in %d iterations", its);
86 switch (reason)
87 {
88 case KSP_CONVERGED_RTOL:
89 PetscPrintf(PETSC_COMM_WORLD,
90 " (relative convergence criterion fulfilled).");
91 break;
92 case KSP_CONVERGED_ATOL:
93 PetscPrintf(PETSC_COMM_WORLD,
94 " (absolute convergence criterion fulfilled).");
95 break;
96 default:
97 PetscPrintf(PETSC_COMM_WORLD, ".");
98 }
99
100 PetscPrintf(PETSC_COMM_WORLD,
101 "\n================================================\n");
102 }
103 else if (reason == KSP_DIVERGED_ITS)
104 {
105 PetscPrintf(PETSC_COMM_WORLD,
106 "\nWarning: maximum number of iterations reached.\n");
107 }
108 else
109 {
110 converged = false;
111 if (reason == KSP_DIVERGED_INDEFINITE_PC)
112 {
113 PetscPrintf(PETSC_COMM_WORLD,
114 "\nDivergence because of indefinite preconditioner,");
115 PetscPrintf(PETSC_COMM_WORLD,
116 "\nTry to run again with "
117 "-pc_factor_shift_positive_definite option.\n");
118 }
119 else if (reason == KSP_DIVERGED_BREAKDOWN_BICG)
120 {
121 PetscPrintf(PETSC_COMM_WORLD,
122 "\nKSPBICG method was detected so the method could not "
123 "continue to enlarge the Krylov space.");
124 PetscPrintf(PETSC_COMM_WORLD,
125 "\nTry to run again with another solver.\n");
126 }
127 else if (reason == KSP_DIVERGED_NONSYMMETRIC)
128 {
129 PetscPrintf(PETSC_COMM_WORLD,
130 "\nMatrix or preconditioner is unsymmetric but KSP "
131 "requires symmetric.\n");
132 }
133 else
134 {
135 PetscPrintf(PETSC_COMM_WORLD,
136 "\nDivergence detected, use command option "
137 "-ksp_monitor or -log_summary to check the details.\n");
138 }
139 }
140
141#ifdef TEST_MEM_PETSC
142 PetscMemoryGetCurrentUsage(&mem2);
143 PetscPrintf(
144 PETSC_COMM_WORLD,
145 "###Memory usage by solver. Before: %f After: %f Increase: %d\n", mem1,
146 mem2, (int)(mem2 - mem1));
147#endif
148
149 elapsed_ctime_ += wtimer.elapsed();
150
151 return converged;
152}
153
154} // namespace MathLib
Declaration of class PETScLinearSolver, which defines a solver object based on PETSc routines.
Definition of the RunTime class.
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
double elapsed_ctime_
Clock time.
PETScLinearSolver(std::string const &prefix, std::string const &petsc_options)
PC pc_
Preconditioner type.
bool solve(PETScMatrix &A, PETScVector &b, PETScVector &x)
bool canSolverHandleRectangular(std::string_view ksp_type)
Wrapper class for PETSc matrix routines for matrix.
Definition PETScMatrix.h:32
Mat & getRawMatrix()
Get matrix reference.
Definition PETScMatrix.h:86
Wrapper class for PETSc vector.
Definition PETScVector.h:40
PETSc_Vec & getRawVector()
Exposes the underlying PETSc vector.