OGS
PETScNonlinearSolver.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#ifdef USE_PETSC
5
7
8#include <petscmat.h>
9#include <petscvec.h>
10
11#include "BaseLib/RunTime.h"
12
13namespace
14{
16{
19 std::vector<GlobalVector*>& x;
20 std::vector<GlobalVector*> const& x_prev;
23 int const process_id;
24};
25
26void printConvergenceReason(SNES const& snes_solver,
27 SNESConvergedReason const& reason)
28{
29 const char* strreason;
30 SNESGetConvergedReasonString(snes_solver, &strreason);
31 if (reason < 0)
32 {
33 INFO("PETScSNES diverged reason: {:s}.", strreason);
34 }
35 else
36 {
37 INFO("PETScSNES converged reason: {:s}.", strreason);
38 }
39}
40
41PetscErrorCode updateResidual(SNES /*snes*/, Vec x, Vec petsc_r,
42 void* petsc_context)
43{
44 auto context = static_cast<PetscContext*>(petsc_context);
45
46 DBUG("PETScNonlinearSolver: residual callback called.");
47
48 VecCopy(x, context->x[context->process_id]->getRawVector());
49
50 // Assemble in ogs context.
51 BaseLib::RunTime time_assembly;
52 time_assembly.start();
53 context->system->assemble(context->x, context->x_prev, context->process_id);
54
55 INFO("[time] Assembly took {} s.", time_assembly.elapsed());
56 context->system->getResidual(*context->x[context->process_id],
57 *context->x_prev[context->process_id],
58 *context->r);
59 context->r->finalizeAssembly();
60 context->J->finalizeAssembly();
61
62 context->system->getJacobian(*context->J);
63 context->system->applyKnownSolutionsPETScSNES(
64 *context->J, *context->r, *context->x[context->process_id]);
65
66 VecCopy(context->r->getRawVector(), petsc_r);
67
68 return 0;
69}
70
71PetscErrorCode updateJacobian(SNES /*snes*/, Vec /*x*/, Mat J,
72 Mat /*same as J*/, void* petsc_context)
73{
74 DBUG("PETScNonlinearSolver: jacobian callback called.");
75 // Assume the system is already assembled.
76 // Copy the results into petsc vectors.
77
78 auto context = static_cast<PetscContext*>(petsc_context);
79 MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
80 MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
81 MatCopy(context->J->getRawMatrix(), J, DIFFERENT_NONZERO_PATTERN);
82 return 0;
83}
84} // namespace
85
86namespace NumLib
87{
89 GlobalLinearSolver& /*linear_solver*/, int const maxiter,
90 std::string prefix)
91{
92 SNESCreate(PETSC_COMM_WORLD, &_snes_solver);
93 if (!prefix.empty())
94 {
95 prefix = prefix + "_";
96 SNESSetOptionsPrefix(_snes_solver, prefix.c_str());
97 }
98 // force SNESSolve() to take at least one iteration regardless of the
99 // initial residual norm
100 SNESSetForceIteration(_snes_solver, PETSC_TRUE);
101
102 // Set the maximum iterations.
103 PetscReal atol, rtol, stol;
104 PetscInt maxf;
105 SNESGetTolerances(_snes_solver, &atol, &rtol, &stol, nullptr, &maxf);
106 SNESSetTolerances(_snes_solver, atol, rtol, stol, maxiter, maxf);
107
108 SNESSetFromOptions(_snes_solver);
109#ifndef NDEBUG
110 PetscOptionsView(nullptr, PETSC_VIEWER_STDOUT_WORLD);
111#endif // NDEBUG
112}
113
120
126
128 std::vector<GlobalVector*> const& /*x*/,
129 std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/)
130{
132 {
133 return;
134 }
135
136 OGS_FATAL(
137 "Non-equilibrium initial residuum is not implemented for the "
138 "PETScNonlinearSolver.");
139}
140
142 std::vector<GlobalVector*>& x,
143 std::vector<GlobalVector*> const& x_prev,
144 std::function<void(
145 int, bool,
146 std::vector<GlobalVector*> const&)> const& /*postIterationCallback*/,
147 int const process_id)
148{
149 DBUG("PETScNonlinearSolver: solve()");
150 using TimeDiscretizedSystem =
153
154 auto* system = static_cast<TimeDiscretizedSystem*>(_equation_system);
155
156 DBUG("PETScNonlinearSolver: create vectors");
157 // r and J on which the ogs assembly operates.
158 auto& J = NumLib::GlobalMatrixProvider::provider.getMatrix(
159 system->getMatrixSpecifications(process_id), _jacobian_id);
160 auto& r = NumLib::GlobalVectorProvider::provider.getVector(
161 system->getMatrixSpecifications(process_id), _residual_id);
162
163 // Vectors and matrices used by SNES for solutions. These will be locked
164 // during the SNES' solve call.
165 auto& J_snes = NumLib::GlobalMatrixProvider::provider.getMatrix(
166 system->getMatrixSpecifications(process_id), _petsc_jacobian_id);
167 BaseLib::RunTime timer_dirichlet;
168 double time_dirichlet = 0.0;
169
170 timer_dirichlet.start();
171 system->computeKnownSolutions(*x[process_id], process_id);
172 system->applyKnownSolutions(*x[process_id]);
173 time_dirichlet += timer_dirichlet.elapsed();
174 INFO("[time] Applying Dirichlet BCs took {} s.", time_dirichlet);
175
176 auto& r_snes = NumLib::GlobalVectorProvider::provider.getVector(
177 system->getMatrixSpecifications(process_id), _petsc_residual_id);
178 auto& x_snes = NumLib::GlobalVectorProvider::provider.getVector(
179 system->getMatrixSpecifications(process_id), _petsc_x_id);
180 MathLib::LinAlg::copy(*x[process_id], x_snes); // Initial guess.
181
182 PetscContext petsc_context{_equation_system, x, x_prev, &r, &J, process_id};
183
184 DBUG("PETScNonlinearSolver: set function");
185 SNESSetFunction(_snes_solver, r_snes.getRawVector(), updateResidual,
186 &petsc_context);
187
188 DBUG("PETScNonlinearSolver: set jacobian");
189 // The jacobian and the preconditioner matrices point to the same location.
190 SNESSetJacobian(_snes_solver, J_snes.getRawMatrix(), J_snes.getRawMatrix(),
191 updateJacobian, &petsc_context);
192
193 std::unique_ptr<GlobalVector> xl = nullptr;
194 std::unique_ptr<GlobalVector> xu = nullptr;
195
196 SNESType snes_type;
197 SNESGetType(_snes_solver, &snes_type);
198 if ((std::strcmp(snes_type, SNESVINEWTONRSLS) == 0) ||
199 (std::strcmp(snes_type, SNESVINEWTONSSLS) == 0))
200 {
201 // Set optional constraints via callback.
202 DBUG("PETScNonlinearSolver: set constraints");
204 system->getMatrixSpecifications(process_id));
206 system->getMatrixSpecifications(process_id));
207
208 system->updateConstraints(*xl, *xu, process_id);
211
212 SNESVISetVariableBounds(_snes_solver, xl->getRawVector(),
213 xu->getRawVector());
214 }
215
216 DBUG("PETScNonlinearSolver: call SNESSolve");
217 SNESSolve(_snes_solver, nullptr, x_snes.getRawVector());
218
219 SNESConvergedReason reason;
220 SNESGetConvergedReason(_snes_solver, &reason);
221 printConvergenceReason(_snes_solver, reason);
222
223 PetscInt iterations;
224 SNESGetIterationNumber(_snes_solver, &iterations);
225 INFO("PETScSNES used {} iterations.", iterations);
226
227 // Copy back the solution.
228 MathLib::LinAlg::copy(x_snes, *x[process_id]);
229
230 NumLib::GlobalVectorProvider::provider.releaseVector(x_snes);
231 NumLib::GlobalVectorProvider::provider.releaseVector(r_snes);
232 NumLib::GlobalMatrixProvider::provider.releaseMatrix(J_snes);
235
236 return {reason >= 0, iterations};
237}
238
239} // namespace NumLib
240#endif // USE_PETSC
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenLisLinearSolver GlobalLinearSolver
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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
NonlinearSolverStatus solve(std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::function< void(int, bool, std::vector< GlobalVector * > const &)> const &postIterationCallback, int const process_id) override
ConvergenceCriterion * _convergence_criterion
void setEquationSystem(System &eq, ConvergenceCriterion &conv_crit)
Set the nonlinear equation system that will be solved.
PETScNonlinearSolver(GlobalLinearSolver &linear_solver, int maxiter, std::string prefix)
std::size_t _jacobian_id
ID of the Jacobian matrix.
void compensateNonEquilibriumInitialResiduum(bool const value)
NonlinearSystem< NonlinearSolverTag::Newton > System
Type of the nonlinear equation system to be solved.
std::size_t _residual_id
ID of the residual vector.
void calculateNonEquilibriumInitialResiduum(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
@ FirstOrderImplicitQuasilinear
Definition Types.h:27
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
PetscErrorCode updateResidual(SNES, Vec x, Vec petsc_r, void *petsc_context)
void printConvergenceReason(SNES const &snes_solver, SNESConvergedReason const &reason)
PetscErrorCode updateJacobian(SNES, Vec, Mat J, Mat, void *petsc_context)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider
Status of the non-linear solver.
NumLib::NonlinearSystem< NumLib::NonlinearSolverTag::Newton > System