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