OGS
NonlinearSolver.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 "NonlinearSolver.h"
5
6#include "BaseLib/Error.h"
7#include "BaseLib/Logging.h"
8#include "BaseLib/MPI.h"
9#include "BaseLib/RunTime.h"
13#include "NumLib/Exceptions.h"
14
15#ifdef USE_PETSC
17#endif // USE_PETSC
18
19namespace NumLib
20{
21namespace detail
22{
23#if !defined(USE_PETSC) && !defined(USE_LIS)
24bool solvePicard(GlobalLinearSolver& linear_solver, GlobalMatrix& A,
26 MathLib::LinearSolverBehaviour const linear_solver_behaviour)
27{
28 BaseLib::RunTime time_linear_solver;
29 time_linear_solver.start();
30
31 if (!linear_solver.compute(A, linear_solver_behaviour))
32 {
33 ERR("Picard: The linear solver failed in the compute() step.");
34 return false;
35 }
36
37 bool const iteration_succeeded = linear_solver.solve(rhs, x);
38
39 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
40
41 if (iteration_succeeded)
42 {
43 return true;
44 }
45
46 ERR("Picard: The linear solver failed in the solve() step.");
47 return false;
48}
49#else
52 MathLib::LinearSolverBehaviour const linear_solver_behaviour)
53{
54 if (linear_solver_behaviour ==
56 linear_solver_behaviour == MathLib::LinearSolverBehaviour::REUSE)
57 {
58 WARN(
59 "The performance optimization to skip the linear solver compute() "
60 "step is not implemented for PETSc or LIS linear solvers.");
61 }
62
63 BaseLib::RunTime time_linear_solver;
64 time_linear_solver.start();
65
66 bool const iteration_succeeded = linear_solver.solve(A, rhs, x);
67
68 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
69
70 if (iteration_succeeded)
71 {
72 return true;
73 }
74
75 ERR("Picard: The linear solver failed in the solve() step.");
76 return false;
77}
78#endif
79} // namespace detail
80
83 std::vector<GlobalVector*> const& x,
84 std::vector<GlobalVector*> const& x_prev, int const process_id)
85{
87 {
88 return;
89 }
90
91 INFO("Calculate non-equilibrium initial residuum.");
92
95 _equation_system->assemble(x, x_prev, process_id);
96 _equation_system->getA(A);
97 _equation_system->getRhs(*x_prev[process_id], rhs);
98
99 // r_neq = A * x - rhs
101 MathLib::LinAlg::matMult(A, *x[process_id], *_r_neq);
102 MathLib::LinAlg::axpy(*_r_neq, -1.0, rhs); // res -= rhs
103
104 // Set the values of the selected entries of _r_neq, which are associated
105 // with the equations that do not need initial residual compensation, to
106 // zero.
107 auto selected_global_indices =
108 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
109
110#ifdef USE_PETSC
111 // Ghost entry with global index 0 is encoded as -global_size
112 // After abs(), it appears as global_size and must be converted back to 0
113 auto const global_size = _r_neq->size();
114 for (auto& idx : selected_global_indices)
115 {
116 if (idx == global_size)
117 {
118 idx = 0;
119 }
120 }
121#endif
122
123 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
124 _r_neq->set(selected_global_indices, zero_entries);
125 _equation_system->setReleaseNodalForces(_r_neq, process_id);
126
128
131}
132
134 std::vector<GlobalVector*>& x,
135 std::vector<GlobalVector*> const& x_prev,
136 std::function<void(int, bool, std::vector<GlobalVector*> const&)> const&
137 postIterationCallback,
138 int const process_id)
139{
140 namespace LinAlg = MathLib::LinAlg;
141 auto& sys = *_equation_system;
142
144 auto& rhs = NumLib::GlobalVectorProvider::provider.getVector(_rhs_id);
145
146 std::vector<GlobalVector*> x_new{x};
147 x_new[process_id] =
149 LinAlg::copy(*x[process_id], *x_new[process_id]); // set initial guess
150
151 bool error_norms_met = false;
152
153 _convergence_criterion->preFirstIteration();
154
155 int iteration = 1;
156 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
157 {
158 BaseLib::RunTime timer_dirichlet;
159 double time_dirichlet = 0.0;
160
161 BaseLib::RunTime time_iteration;
162 time_iteration.start();
163
164 INFO("Iteration #{:d} started.", iteration);
165 timer_dirichlet.start();
166 auto& x_new_process = *x_new[process_id];
168 sys.computeKnownSolutions(x_new_process, process_id);
169 sys.applyKnownSolutions(x_new_process);
170 time_dirichlet += timer_dirichlet.elapsed();
171
172 sys.preIteration(iteration, x_new_process);
173
174 BaseLib::RunTime time_assembly;
175 time_assembly.start();
176 sys.assemble(x_new, x_prev, process_id);
177 sys.getA(A);
178 sys.getRhs(*x_prev[process_id], rhs);
179
180 // Normalize the linear equation system, if required
181 if (sys.requiresNormalization() &&
182 !_linear_solver.canSolveRectangular())
183 {
184 sys.getAandRhsNormalized(A, rhs);
185 WARN(
186 "The equation system is rectangular, but the current linear "
187 "solver only supports square systems. "
188 "The system will be normalized, which lead to a squared "
189 "condition number and potential numerical issues. "
190 "It is recommended to use a solver that supports rectangular "
191 "equation systems for better numerical stability.");
192 }
193
194 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
195
196 // Subtract non-equilibrium initial residuum if set
197 if (_r_neq != nullptr)
198 {
199 LinAlg::axpy(rhs, -1, *_r_neq);
200 }
201
202 auto const solver_needs_to_compute = sys.linearSolverNeedsToCompute();
203 bool const solver_will_compute =
204 _linear_solver.willCompute(solver_needs_to_compute);
205
206 timer_dirichlet.start();
207 sys.applyKnownSolutionsPicard(
208 A, rhs, x_new_process,
209 solver_will_compute
212 FAST_INCOMPLETE_MATRIX_UPDATE);
213 time_dirichlet += timer_dirichlet.elapsed();
214 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
215
216 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
217 {
218 if (!solver_will_compute)
219 {
220 // !solver_will_compute means that the Dirichlet BC application
221 // is incomplete (i.e., A not properly modified) and the
222 // computed residual is wrong.
223 OGS_FATAL(
224 "Logic error. The solver skips the compute step for a "
225 "non-linear equation system.");
226 }
227 GlobalVector res;
228 LinAlg::matMult(A, x_new_process, res); // res = A * x_new
229 LinAlg::axpy(res, -1.0, rhs); // res -= rhs
230 _convergence_criterion->checkResidual(res);
231 }
232
233 bool iteration_succeeded = detail::solvePicard(
234 _linear_solver, A, rhs, x_new_process, solver_needs_to_compute);
235
236 if (iteration_succeeded)
237 {
238 if (postIterationCallback)
239 {
240 postIterationCallback(iteration, error_norms_met, x_new);
241 }
242
243 switch (sys.postIteration(x_new_process))
244 {
246 // Don't copy here. The old x might still be used further
247 // below. Although currently it is not.
248 break;
250 ERR("Picard: The postIteration() hook reported a "
251 "non-recoverable error.");
252 iteration_succeeded = false;
253 // Copy new solution to x.
254 // Thereby the failed solution can be used by the caller for
255 // debugging purposes.
256 LinAlg::copy(x_new_process, *x[process_id]);
257 break;
259 INFO(
260 "Picard: The postIteration() hook decided that this "
261 "iteration has to be repeated.");
263 *x[process_id],
264 x_new_process); // throw the iteration result away
265 continue;
266 }
267 }
268
269 if (!iteration_succeeded)
270 {
271 // Don't compute error norms, break here.
272 error_norms_met = false;
273 break;
274 }
275
276 if (sys.isLinear())
277 {
278 error_norms_met = true;
279 }
280 else
281 {
282 if (_convergence_criterion->hasDeltaXCheck())
283 {
284 GlobalVector minus_delta_x(*x[process_id]);
285 LinAlg::axpy(minus_delta_x, -1.0,
286 x_new_process); // minus_delta_x = x - x_new
287 _convergence_criterion->checkDeltaX(minus_delta_x,
288 x_new_process);
289 }
290
291 error_norms_met = _convergence_criterion->isSatisfied();
292 }
293
294 // Update x s.t. in the next iteration we will compute the right delta x
295 LinAlg::copy(x_new_process, *x[process_id]);
296
297 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
298 time_iteration.elapsed());
299
300 if (error_norms_met)
301 {
302 break;
303 }
304
305 // Avoid increment of the 'iteration' if the error norms are not met,
306 // but maximum number of iterations is reached.
307 if (iteration >= _maxiter)
308 {
309 break;
310 }
311 }
312
313 if (iteration > _maxiter)
314 {
315 ERR("Picard: Could not solve the given nonlinear system within {:d} "
316 "iterations",
317 _maxiter);
318 }
319
322 NumLib::GlobalVectorProvider::provider.releaseVector(*x_new[process_id]);
323
324 return {error_norms_met, iteration};
325}
326
329 std::vector<GlobalVector*> const& x,
330 std::vector<GlobalVector*> const& x_prev, int const process_id)
331{
333 {
334 return;
335 }
336
337 INFO("Calculate non-equilibrium initial residuum.");
338
339 _equation_system->assemble(x, x_prev, process_id);
341 _equation_system->getResidual(*x[process_id], *x_prev[process_id], *_r_neq);
342
343 // Set the values of the selected entries of _r_neq, which are associated
344 // with the equations that do not need initial residual compensation, to
345 // zero.
346 auto selected_global_indices =
347 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
348
349#ifdef USE_PETSC
350 // Ghost entry with global index 0 is encoded as -global_size
351 // After abs(), it appears as global_size and must be converted back to 0
352 auto const global_size = _r_neq->size();
353 for (auto& idx : selected_global_indices)
354 {
355 if (idx == global_size)
356 {
357 idx = 0;
358 }
359 }
360#endif
361
362 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
363 _r_neq->set(selected_global_indices, zero_entries);
364 _equation_system->setReleaseNodalForces(_r_neq, process_id);
365
367}
368
370 std::vector<GlobalVector*>& x,
371 std::vector<GlobalVector*> const& x_prev,
372 std::function<void(int, bool, std::vector<GlobalVector*> const&)> const&
373 postIterationCallback,
374 int const process_id)
375{
376 namespace LinAlg = MathLib::LinAlg;
377 auto& sys = *_equation_system;
378
379 auto& res = NumLib::GlobalVectorProvider::provider.getVector(_res_id);
380 auto& minus_delta_x =
383
384 bool error_norms_met = false;
385
386 // TODO be more efficient
387 // init minus_delta_x to the right size
388 LinAlg::copy(*x[process_id], minus_delta_x);
389
390 _convergence_criterion->preFirstIteration();
391
392 NewtonStepContext step_ctx{sys, x_prev, process_id};
393
394 int iteration = 1;
395#if !defined(USE_PETSC) && !defined(USE_LIS)
396 int next_iteration_inv_jacobian_recompute = 1;
397#endif
398 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
399 {
400 BaseLib::RunTime timer_dirichlet;
401 double time_dirichlet = 0.0;
402
403 BaseLib::RunTime time_iteration;
404 INFO("Iteration #{:d} started.", iteration);
405 time_iteration.start();
406
407 timer_dirichlet.start();
408 sys.computeKnownSolutions(*x[process_id], process_id);
409 time_dirichlet += timer_dirichlet.elapsed();
410
411 sys.preIteration(iteration, *x[process_id]);
412
413 BaseLib::RunTime time_assembly;
414 time_assembly.start();
415 bool mpi_rank_assembly_ok = true;
416 try
417 {
418 sys.assemble(x, x_prev, process_id);
419 }
420 catch (AssemblyException const& e)
421 {
422 ERR("Abort nonlinear iteration. Repeating timestep. Reason: {:s}",
423 e.what());
424 error_norms_met = false;
425 iteration = _maxiter;
426 mpi_rank_assembly_ok = false;
427 }
428 if (BaseLib::MPI::anyOf(!mpi_rank_assembly_ok))
429 {
430 break;
431 }
432 sys.getResidual(*x[process_id], *x_prev[process_id], res);
433 sys.getJacobian(J);
434 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
435
436 // Subtract non-equilibrium initial residuum if set
437 if (_r_neq != nullptr)
438 {
439 LinAlg::axpy(res, -1, *_r_neq);
440 }
441
442 minus_delta_x.setZero();
443
444 timer_dirichlet.start();
445 sys.applyKnownSolutionsNewton(J, res, *x[process_id], minus_delta_x);
446 time_dirichlet += timer_dirichlet.elapsed();
447 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
448
449 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
450 {
451 _convergence_criterion->checkResidual(res);
452 }
453
454 BaseLib::RunTime time_linear_solver;
455 time_linear_solver.start();
456#if !defined(USE_PETSC) && !defined(USE_LIS)
457 auto linear_solver_behaviour = MathLib::LinearSolverBehaviour::REUSE;
458 if (iteration == next_iteration_inv_jacobian_recompute)
459 {
460 linear_solver_behaviour =
462 next_iteration_inv_jacobian_recompute =
463 next_iteration_inv_jacobian_recompute + _recompute_jacobian;
464 }
465
466 bool iteration_succeeded = false;
467 if (!_linear_solver.compute(J, linear_solver_behaviour))
468 {
469 ERR("Newton: The linear solver failed in the compute() step.");
470 }
471 else
472 {
473 iteration_succeeded = _linear_solver.solve(res, minus_delta_x);
474 }
475#else
476 bool iteration_succeeded = _linear_solver.solve(J, res, minus_delta_x);
477#endif
478 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
479
480 if (!iteration_succeeded)
481 {
482 ERR("Newton: The linear solver failed.");
483 }
484 else
485 {
486 // TODO could be solved in a better way
487 // cf.
488 // https://petsc.org/release/manualpages/Vec/VecWAXPY
489
490 // Copy pointers, replace the one for the given process id.
491 std::vector<GlobalVector*> x_new{x};
492 x_new[process_id] =
494 *x[process_id], _x_new_id);
495 auto const step_result = _step_strategy->applyStep(
496 *x[process_id], minus_delta_x, res, J, *x_new[process_id],
497 step_ctx, iteration);
498
499 if (step_result.step_length != 1.0)
500 {
501 INFO("Step length: {:g}", step_result.step_length);
502 }
503
504 if (!step_result.success)
505 {
506 ERR("Newton: step strategy failed.");
507 iteration_succeeded = false;
508 }
509 else if (!step_result.x_new_is_set)
510 {
511 LinAlg::axpy(*x_new[process_id], -1.0, minus_delta_x);
512 }
513
514 if (postIterationCallback)
515 {
516 postIterationCallback(iteration, error_norms_met, x_new);
517 }
518
519 switch (sys.postIteration(*x_new[process_id]))
520 {
522 break;
524 ERR("Newton: The postIteration() hook reported a "
525 "non-recoverable error.");
526 iteration_succeeded = false;
527 break;
529 INFO(
530 "Newton: The postIteration() hook decided that this "
531 "iteration has to be repeated.");
532 // TODO introduce some onDestroy hook.
534 *x_new[process_id]);
535 continue; // That throws the iteration result away.
536 }
537
538 LinAlg::copy(*x_new[process_id],
539 *x[process_id]); // copy new solution to x
541 *x_new[process_id]);
542 }
543
544 if (!iteration_succeeded)
545 {
546 // Don't compute further error norms, but break here.
547 error_norms_met = false;
548 break;
549 }
550
551 if (sys.isLinear())
552 {
553 error_norms_met = true;
554 }
555 else
556 {
557 if (_convergence_criterion->hasDeltaXCheck())
558 {
559 // Note: x contains the new solution!
560 _convergence_criterion->checkDeltaX(minus_delta_x,
561 *x[process_id]);
562 }
563
564 error_norms_met = _convergence_criterion->isSatisfied();
565 }
566
567 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
568 time_iteration.elapsed());
569
570 if (error_norms_met)
571 {
572 break;
573 }
574
575 // Avoid increment of the 'iteration' if the error norms are not met,
576 // but maximum number of iterations is reached.
577 if (iteration >= _maxiter)
578 {
579 break;
580 }
581 }
582
583 if (iteration > _maxiter)
584 {
585 ERR("Newton: Could not solve the given nonlinear system within {:d} "
586 "iterations",
587 _maxiter);
588 }
589
592 NumLib::GlobalVectorProvider::provider.releaseVector(minus_delta_x);
593
594 return {error_norms_met, iteration};
595}
596
604
612
613} // namespace NumLib
#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 ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
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
bool solve(EigenMatrix &A, EigenVector &b, EigenVector &x)
ConvergenceCriterion * _convergence_criterion
Convergence criterion used to terminate the Newton iteration.
std::size_t _J_id
ID of the Jacobian matrix.
std::size_t _x_new_id
ID of the vector storing .
std::size_t _res_id
ID of the residual vector.
GlobalVector * _r_neq
non-equilibrium initial residuum.
int const _maxiter
maximum number of iterations
int const _recompute_jacobian
Recompute Jacobian every this many steps.
std::unique_ptr< NewtonStepStrategy > _step_strategy
Globalization / step-acceptance strategy (e.g. fixed damping).
NonlinearSolver(GlobalLinearSolver &linear_solver, int const maxiter, std::unique_ptr< NewtonStepStrategy > newton_strategy, int const recompute_jacobian=1)
std::size_t _rhs_id
ID of the right-hand side vector.
const int _maxiter
maximum number of iterations
GlobalVector * _r_neq
non-equilibrium initial residuum.
NonlinearSolver(GlobalLinearSolver &linear_solver, const int maxiter)
static bool anyOf(bool const val, Mpi const &mpi=Mpi{OGS_COMM_WORLD})
Definition MPI.h:174
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:191
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:142
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
DirichletBCApplicationMode
Definition LinAlgEnums.h:33
@ COMPLETE_MATRIX_UPDATE
Both A and b fully updated.
Definition LinAlgEnums.h:34
bool solvePicard(GlobalLinearSolver &linear_solver, GlobalMatrix &A, GlobalVector &rhs, GlobalVector &x, MathLib::LinearSolverBehaviour const linear_solver_behaviour)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider
Status of the non-linear solver.