OGS
NonlinearSolver.cpp
Go to the documentation of this file.
1
11#include "NonlinearSolver.h"
12
13#include <boost/algorithm/string.hpp>
14
15#include "BaseLib/ConfigTree.h"
16#include "BaseLib/Error.h"
17#include "BaseLib/Logging.h"
18#include "BaseLib/RunTime.h"
22#include "NumLib/Exceptions.h"
24
25namespace NumLib
26{
27namespace detail
28{
29#if !defined(USE_PETSC) && !defined(USE_LIS)
30bool solvePicard(GlobalLinearSolver& linear_solver, GlobalMatrix& A,
32 MathLib::LinearSolverBehaviour const linear_solver_behaviour)
33{
34 BaseLib::RunTime time_linear_solver;
35 time_linear_solver.start();
36
37 if (linear_solver_behaviour == MathLib::LinearSolverBehaviour::RECOMPUTE ||
38 linear_solver_behaviour ==
40 {
41 if (!linear_solver.compute(A, linear_solver_behaviour))
42 {
43 ERR("Picard: The linear solver failed in the compute() step.");
44 return false;
45 }
46 }
47
48 // REUSE the previously computed preconditioner or LU decomposition
49 bool const iteration_succeeded = linear_solver.solve(rhs, x);
50
51 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
52
53 if (iteration_succeeded)
54 {
55 return true;
56 }
57
58 ERR("Picard: The linear solver failed in the solve() step.");
59 return false;
60}
61#else
64 MathLib::LinearSolverBehaviour const linear_solver_behaviour)
65{
66 if (linear_solver_behaviour ==
68 linear_solver_behaviour == MathLib::LinearSolverBehaviour::REUSE)
69 {
70 WARN(
71 "The performance optimization to skip the linear solver compute() "
72 "step is not implemented for PETSc or LIS linear solvers.");
73 }
74
75 BaseLib::RunTime time_linear_solver;
76 time_linear_solver.start();
77
78 bool const iteration_succeeded = linear_solver.solve(A, rhs, x);
79
80 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
81
82 if (iteration_succeeded)
83 {
84 return true;
85 }
86
87 ERR("Picard: The linear solver failed in the solve() step.");
88 return false;
89}
90#endif
91} // namespace detail
92
95 std::vector<GlobalVector*> const& x,
96 std::vector<GlobalVector*> const& x_prev, int const process_id)
97{
98 if (!_compensate_non_equilibrium_initial_residuum)
99 {
100 return;
101 }
102
103 INFO("Calculate non-equilibrium initial residuum.");
104
107 _equation_system->assemble(x, x_prev, process_id);
108 _equation_system->getA(A);
109 _equation_system->getRhs(*x_prev[process_id], rhs);
110
111 // r_neq = A * x - rhs
113 MathLib::LinAlg::matMult(A, *x[process_id], *_r_neq);
114 MathLib::LinAlg::axpy(*_r_neq, -1.0, rhs); // res -= rhs
115
116 // Set the values of the selected entries of _r_neq, which are associated
117 // with the equations that do not need initial residual compensation, to
118 // zero.
119 const auto selected_global_indices =
120 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
121
122 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
123 _r_neq->set(selected_global_indices, zero_entries);
124
126
129}
130
132 std::vector<GlobalVector*>& x,
133 std::vector<GlobalVector*> const& x_prev,
134 std::function<void(int, std::vector<GlobalVector*> const&)> const&
135 postIterationCallback,
136 int const process_id)
137{
138 namespace LinAlg = MathLib::LinAlg;
139 auto& sys = *_equation_system;
140
143
144 std::vector<GlobalVector*> x_new{x};
145 x_new[process_id] =
147 LinAlg::copy(*x[process_id], *x_new[process_id]); // set initial guess
148
149 bool error_norms_met = false;
150
151 _convergence_criterion->preFirstIteration();
152
153 int iteration = 1;
154 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
155 {
156 BaseLib::RunTime timer_dirichlet;
157 double time_dirichlet = 0.0;
158
159 BaseLib::RunTime time_iteration;
160 time_iteration.start();
161
162 timer_dirichlet.start();
163 auto& x_new_process = *x_new[process_id];
165 sys.computeKnownSolutions(x_new_process, process_id);
166 sys.applyKnownSolutions(x_new_process);
167 time_dirichlet += timer_dirichlet.elapsed();
168
169 sys.preIteration(iteration, x_new_process);
170
171 BaseLib::RunTime time_assembly;
172 time_assembly.start();
173 sys.assemble(x_new, x_prev, process_id);
174 sys.getA(A);
175 sys.getRhs(*x_prev[process_id], rhs);
176 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
177
178 // Subtract non-equilibrium initial residuum if set
179 if (_r_neq != nullptr)
180 {
181 LinAlg::axpy(rhs, -1, *_r_neq);
182 }
183
184 timer_dirichlet.start();
185 sys.applyKnownSolutionsPicard(A, rhs, x_new_process);
186 time_dirichlet += timer_dirichlet.elapsed();
187 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
188
189 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
190 {
191 GlobalVector res;
192 LinAlg::matMult(A, x_new_process, res); // res = A * x_new
193 LinAlg::axpy(res, -1.0, rhs); // res -= rhs
194 _convergence_criterion->checkResidual(res);
195 }
196
197 bool iteration_succeeded =
198 detail::solvePicard(_linear_solver, A, rhs, x_new_process,
199 sys.linearSolverNeedsToCompute());
200
201 if (iteration_succeeded)
202 {
203 if (postIterationCallback)
204 {
205 postIterationCallback(iteration, x_new);
206 }
207
208 switch (sys.postIteration(x_new_process))
209 {
211 // Don't copy here. The old x might still be used further
212 // below. Although currently it is not.
213 break;
215 ERR("Picard: The postIteration() hook reported a "
216 "non-recoverable error.");
217 iteration_succeeded = false;
218 // Copy new solution to x.
219 // Thereby the failed solution can be used by the caller for
220 // debugging purposes.
221 LinAlg::copy(x_new_process, *x[process_id]);
222 break;
224 INFO(
225 "Picard: The postIteration() hook decided that this "
226 "iteration has to be repeated.");
228 *x[process_id],
229 x_new_process); // throw the iteration result away
230 continue;
231 }
232 }
233
234 if (!iteration_succeeded)
235 {
236 // Don't compute error norms, break here.
237 error_norms_met = false;
238 break;
239 }
240
241 if (sys.isLinear())
242 {
243 error_norms_met = true;
244 }
245 else
246 {
247 if (_convergence_criterion->hasDeltaXCheck())
248 {
249 GlobalVector minus_delta_x(*x[process_id]);
250 LinAlg::axpy(minus_delta_x, -1.0,
251 x_new_process); // minus_delta_x = x - x_new
252 _convergence_criterion->checkDeltaX(minus_delta_x,
253 x_new_process);
254 }
255
256 error_norms_met = _convergence_criterion->isSatisfied();
257 }
258
259 // Update x s.t. in the next iteration we will compute the right delta x
260 LinAlg::copy(x_new_process, *x[process_id]);
261
262 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
263 time_iteration.elapsed());
264
265 if (error_norms_met)
266 {
267 break;
268 }
269
270 // Avoid increment of the 'iteration' if the error norms are not met,
271 // but maximum number of iterations is reached.
272 if (iteration >= _maxiter)
273 {
274 break;
275 }
276 }
277
278 if (iteration > _maxiter)
279 {
280 ERR("Picard: Could not solve the given nonlinear system within {:d} "
281 "iterations",
282 _maxiter);
283 }
284
288
289 return {error_norms_met, iteration};
290}
291
294 std::vector<GlobalVector*> const& x,
295 std::vector<GlobalVector*> const& x_prev, int const process_id)
296{
297 if (!_compensate_non_equilibrium_initial_residuum)
298 {
299 return;
300 }
301
302 INFO("Calculate non-equilibrium initial residuum.");
303
304 _equation_system->assemble(x, x_prev, process_id);
306 _equation_system->getResidual(*x[process_id], *x_prev[process_id], *_r_neq);
307
308 // Set the values of the selected entries of _r_neq, which are associated
309 // with the equations that do not need initial residual compensation, to
310 // zero.
311 const auto selected_global_indices =
312 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
313 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
314
315 _r_neq->set(selected_global_indices, zero_entries);
316
318}
319
321 std::vector<GlobalVector*>& x,
322 std::vector<GlobalVector*> const& x_prev,
323 std::function<void(int, std::vector<GlobalVector*> const&)> const&
324 postIterationCallback,
325 int const process_id)
326{
327 namespace LinAlg = MathLib::LinAlg;
328 auto& sys = *_equation_system;
329
331 auto& minus_delta_x =
334
335 bool error_norms_met = false;
336
337 // TODO be more efficient
338 // init minus_delta_x to the right size
339 LinAlg::copy(*x[process_id], minus_delta_x);
340
341 _convergence_criterion->preFirstIteration();
342
343 int iteration = 1;
344 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
345 {
346 BaseLib::RunTime timer_dirichlet;
347 double time_dirichlet = 0.0;
348
349 BaseLib::RunTime time_iteration;
350 time_iteration.start();
351
352 timer_dirichlet.start();
353 sys.computeKnownSolutions(*x[process_id], process_id);
354 sys.applyKnownSolutions(*x[process_id]);
355 time_dirichlet += timer_dirichlet.elapsed();
356
357 sys.preIteration(iteration, *x[process_id]);
358
359 BaseLib::RunTime time_assembly;
360 time_assembly.start();
361 try
362 {
363 sys.assemble(x, x_prev, process_id);
364 }
365 catch (AssemblyException const& e)
366 {
367 ERR("Abort nonlinear iteration. Repeating timestep. Reason: {:s}",
368 e.what());
369 error_norms_met = false;
370 iteration = _maxiter;
371 break;
372 }
373 sys.getResidual(*x[process_id], *x_prev[process_id], res);
374 sys.getJacobian(J);
375 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
376
377 // Subtract non-equilibrium initial residuum if set
378 if (_r_neq != nullptr)
379 LinAlg::axpy(res, -1, *_r_neq);
380
381 minus_delta_x.setZero();
382
383 timer_dirichlet.start();
384 sys.applyKnownSolutionsNewton(J, res, minus_delta_x);
385 time_dirichlet += timer_dirichlet.elapsed();
386 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
387
388 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
389 {
390 _convergence_criterion->checkResidual(res);
391 }
392
393 BaseLib::RunTime time_linear_solver;
394 time_linear_solver.start();
395 bool iteration_succeeded = _linear_solver.solve(J, res, minus_delta_x);
396 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
397
398 if (!iteration_succeeded)
399 {
400 ERR("Newton: The linear solver failed.");
401 }
402 else
403 {
404 // TODO could be solved in a better way
405 // cf.
406 // https://petsc.org/release/manualpages/Vec/VecWAXPY
407
408 // Copy pointers, replace the one for the given process id.
409 std::vector<GlobalVector*> x_new{x};
410 x_new[process_id] =
412 *x[process_id], _x_new_id);
413 LinAlg::axpy(*x_new[process_id], -_damping, minus_delta_x);
414
415 if (postIterationCallback)
416 {
417 postIterationCallback(iteration, x_new);
418 }
419
420 switch (sys.postIteration(*x_new[process_id]))
421 {
423 break;
425 ERR("Newton: The postIteration() hook reported a "
426 "non-recoverable error.");
427 iteration_succeeded = false;
428 break;
430 INFO(
431 "Newton: The postIteration() hook decided that this "
432 "iteration has to be repeated.");
433 // TODO introduce some onDestroy hook.
435 *x_new[process_id]);
436 continue; // That throws the iteration result away.
437 }
438
439 LinAlg::copy(*x_new[process_id],
440 *x[process_id]); // copy new solution to x
442 *x_new[process_id]);
443 }
444
445 if (!iteration_succeeded)
446 {
447 // Don't compute further error norms, but break here.
448 error_norms_met = false;
449 break;
450 }
451
452 if (sys.isLinear())
453 {
454 error_norms_met = true;
455 }
456 else
457 {
458 if (_convergence_criterion->hasDeltaXCheck())
459 {
460 // Note: x contains the new solution!
461 _convergence_criterion->checkDeltaX(minus_delta_x,
462 *x[process_id]);
463 }
464
465 error_norms_met = _convergence_criterion->isSatisfied();
466 }
467
468 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
469 time_iteration.elapsed());
470
471 if (error_norms_met)
472 {
473 break;
474 }
475
476 // Avoid increment of the 'iteration' if the error norms are not met,
477 // but maximum number of iterations is reached.
478 if (iteration >= _maxiter)
479 {
480 break;
481 }
482 }
483
484 if (iteration > _maxiter)
485 {
486 ERR("Newton: Could not solve the given nonlinear system within {:d} "
487 "iterations",
488 _maxiter);
489 }
490
494
495 return {error_norms_met, iteration};
496}
497
498std::pair<std::unique_ptr<NonlinearSolverBase>, NonlinearSolverTag>
500 BaseLib::ConfigTree const& config)
501{
503 auto const type = config.getConfigParameter<std::string>("type");
505 auto const max_iter = config.getConfigParameter<int>("max_iter");
506
507 if (type == "Picard")
508 {
509 auto const tag = NonlinearSolverTag::Picard;
510 using ConcreteNLS = NonlinearSolver<tag>;
511 return std::make_pair(
512 std::make_unique<ConcreteNLS>(linear_solver, max_iter), tag);
513 }
514 if (type == "Newton")
515 {
517 auto const damping = config.getConfigParameter<double>("damping", 1.0);
518 if (damping <= 0)
519 {
520 OGS_FATAL(
521 "The damping factor for the Newon method must be positive, got "
522 "{:g}.",
523 damping);
524 }
525 auto const tag = NonlinearSolverTag::Newton;
526 using ConcreteNLS = NonlinearSolver<tag>;
527 return std::make_pair(
528 std::make_unique<ConcreteNLS>(linear_solver, max_iter, damping),
529 tag);
530 }
531#ifdef USE_PETSC
532 if (boost::iequals(type, "PETScSNES"))
533 {
534 auto prefix =
536 config.getConfigParameter<std::string>("prefix", "");
537 auto const tag = NonlinearSolverTag::Newton;
538 using ConcreteNLS = PETScNonlinearSolver;
539 return std::make_pair(std::make_unique<ConcreteNLS>(
540 linear_solver, max_iter, std::move(prefix)),
541 tag);
542 }
543
544#endif
545 OGS_FATAL("Unsupported nonlinear solver type '{:s}'.", type.c_str());
546}
547
555
563
564} // namespace NumLib
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Definition of the RunTime class.
T getConfigParameter(std::string const &param) const
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
bool solve(EigenMatrix &A, EigenVector &b, EigenVector &x)
Global vector based on Eigen vector.
Definition EigenVector.h:25
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
NonlinearSolverTag
Tag used to specify which nonlinear solver will be used.
Definition Types.h:20
std::pair< std::unique_ptr< NonlinearSolverBase >, NonlinearSolverTag > createNonlinearSolver(GlobalLinearSolver &linear_solver, BaseLib::ConfigTree const &config)
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:149
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57
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.