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
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}
149
150} // 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)
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.