OGS
PETScLinearSolver.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#include "PETScLinearSolver.h"
5
6#include "BaseLib/RunTime.h"
7
8namespace MathLib
9{
10PETScLinearSolver::PETScLinearSolver(std::string const& prefix,
11 std::string const& petsc_options)
12{
13#if PETSC_VERSION_LT(3, 7, 0)
14 PetscOptionsInsertString(petsc_options.c_str());
15#else
16 PetscOptionsInsertString(nullptr, petsc_options.c_str());
17#endif
18
19 KSPCreate(PETSC_COMM_WORLD, &solver_);
20
21 KSPGetPC(solver_, &pc_);
22
23 if (!prefix.empty())
24 {
25 KSPSetOptionsPrefix(solver_, prefix.c_str());
26 }
27
28 KSPSetInitialGuessNonzero(solver_, PETSC_TRUE);
29 KSPSetFromOptions(solver_); // set run-time options
30
31 KSPType ksp_type;
32 KSPGetType(solver_, &ksp_type);
34}
35
37{
38 BaseLib::RunTime wtimer;
39 wtimer.start();
40
41// define TEST_MEM_PETSC
42#ifdef TEST_MEM_PETSC
43 PetscLogDouble mem1, mem2;
44 PetscMemoryGetCurrentUsage(&mem1);
45#endif
46
47 KSPNormType norm_type;
48 KSPGetNormType(solver_, &norm_type);
49 const char* ksp_type;
50 const char* pc_type;
51 KSPGetType(solver_, &ksp_type);
52 PCGetType(pc_, &pc_type);
53
54 PetscPrintf(PETSC_COMM_WORLD,
55 "\n================================================");
56 PetscPrintf(PETSC_COMM_WORLD,
57 "\nLinear solver %s with %s preconditioner using %s", ksp_type,
58 pc_type, KSPNormTypes[norm_type]);
59
60 KSPSetOperators(solver_, A.getRawMatrix(), A.getRawMatrix());
61
62 KSPSolve(solver_, b.getRawVector(), x.getRawVector());
63
64 KSPConvergedReason reason;
65 KSPGetConvergedReason(solver_, &reason);
66
67 bool converged = true;
68 if (reason > 0)
69 {
70 PetscInt its;
71 KSPGetIterationNumber(solver_, &its);
72 PetscPrintf(PETSC_COMM_WORLD, "\nconverged in %d iterations", its);
73 switch (reason)
74 {
75 case KSP_CONVERGED_RTOL:
76 PetscPrintf(PETSC_COMM_WORLD,
77 " (relative convergence criterion fulfilled).");
78 break;
79 case KSP_CONVERGED_ATOL:
80 PetscPrintf(PETSC_COMM_WORLD,
81 " (absolute convergence criterion fulfilled).");
82 break;
83 default:
84 PetscPrintf(PETSC_COMM_WORLD, ".");
85 }
86
87 PetscPrintf(PETSC_COMM_WORLD,
88 "\n================================================\n");
89 }
90 else if (reason == KSP_DIVERGED_ITS)
91 {
92 PetscPrintf(PETSC_COMM_WORLD,
93 "\nWarning: maximum number of iterations reached.\n");
94 }
95 else
96 {
97 converged = false;
98 if (reason == KSP_DIVERGED_INDEFINITE_PC)
99 {
100 PetscPrintf(PETSC_COMM_WORLD,
101 "\nDivergence because of indefinite preconditioner,");
102 PetscPrintf(PETSC_COMM_WORLD,
103 "\nTry to run again with "
104 "-pc_factor_shift_positive_definite option.\n");
105 }
106 else if (reason == KSP_DIVERGED_BREAKDOWN_BICG)
107 {
108 PetscPrintf(PETSC_COMM_WORLD,
109 "\nKSPBICG method was detected so the method could not "
110 "continue to enlarge the Krylov space.");
111 PetscPrintf(PETSC_COMM_WORLD,
112 "\nTry to run again with another solver.\n");
113 }
114 else if (reason == KSP_DIVERGED_NONSYMMETRIC)
115 {
116 PetscPrintf(PETSC_COMM_WORLD,
117 "\nMatrix or preconditioner is unsymmetric but KSP "
118 "requires symmetric.\n");
119 }
120 else
121 {
122 PetscPrintf(PETSC_COMM_WORLD,
123 "\nDivergence detected, use command option "
124 "-ksp_monitor or -log_summary to check the details.\n");
125 }
126 }
127
128#ifdef TEST_MEM_PETSC
129 PetscMemoryGetCurrentUsage(&mem2);
130 PetscPrintf(
131 PETSC_COMM_WORLD,
132 "###Memory usage by solver. Before: %f After: %f Increase: %d\n", mem1,
133 mem2, (int)(mem2 - mem1));
134#endif
135
136 elapsed_ctime_ += wtimer.elapsed();
137
138 return converged;
139}
140
141} // namespace MathLib
Count the running time.
Definition RunTime.h:18
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
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:21
Mat & getRawMatrix()
Get matrix reference.
Definition PETScMatrix.h:75
Wrapper class for PETSc vector.
Definition PETScVector.h:28
PETSc_Vec & getRawVector()
Exposes the underlying PETSc vector.