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 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
171
172 // Subtract non-equilibrium initial residuum if set
173 if (_r_neq != nullptr)
174 {
175 LinAlg::axpy(rhs, -1, *_r_neq);
176 }
177
178 timer_dirichlet.start();
179 sys.applyKnownSolutionsPicard(A, rhs, x_new_process);
180 time_dirichlet += timer_dirichlet.elapsed();
181 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
182
183 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
184 {
185 GlobalVector res;
186 LinAlg::matMult(A, x_new_process, res); // res = A * x_new
187 LinAlg::axpy(res, -1.0, rhs); // res -= rhs
188 _convergence_criterion->checkResidual(res);
189 }
190
191 bool iteration_succeeded =
192 detail::solvePicard(_linear_solver, A, rhs, x_new_process,
193 sys.linearSolverNeedsToCompute());
194
195 if (iteration_succeeded)
196 {
197 if (postIterationCallback)
198 {
199 postIterationCallback(iteration, x_new);
200 }
201
202 switch (sys.postIteration(x_new_process))
203 {
205 // Don't copy here. The old x might still be used further
206 // below. Although currently it is not.
207 break;
209 ERR("Picard: The postIteration() hook reported a "
210 "non-recoverable error.");
211 iteration_succeeded = false;
212 // Copy new solution to x.
213 // Thereby the failed solution can be used by the caller for
214 // debugging purposes.
215 LinAlg::copy(x_new_process, *x[process_id]);
216 break;
218 INFO(
219 "Picard: The postIteration() hook decided that this "
220 "iteration has to be repeated.");
222 *x[process_id],
223 x_new_process); // throw the iteration result away
224 continue;
225 }
226 }
227
228 if (!iteration_succeeded)
229 {
230 // Don't compute error norms, break here.
231 error_norms_met = false;
232 break;
233 }
234
235 if (sys.isLinear())
236 {
237 error_norms_met = true;
238 }
239 else
240 {
241 if (_convergence_criterion->hasDeltaXCheck())
242 {
243 GlobalVector minus_delta_x(*x[process_id]);
244 LinAlg::axpy(minus_delta_x, -1.0,
245 x_new_process); // minus_delta_x = x - x_new
246 _convergence_criterion->checkDeltaX(minus_delta_x,
247 x_new_process);
248 }
249
250 error_norms_met = _convergence_criterion->isSatisfied();
251 }
252
253 // Update x s.t. in the next iteration we will compute the right delta x
254 LinAlg::copy(x_new_process, *x[process_id]);
255
256 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
257 time_iteration.elapsed());
258
259 if (error_norms_met)
260 {
261 break;
262 }
263
264 // Avoid increment of the 'iteration' if the error norms are not met,
265 // but maximum number of iterations is reached.
266 if (iteration >= _maxiter)
267 {
268 break;
269 }
270 }
271
272 if (iteration > _maxiter)
273 {
274 ERR("Picard: Could not solve the given nonlinear system within {:d} "
275 "iterations",
276 _maxiter);
277 }
278
282
283 return {error_norms_met, iteration};
284}
285
288 std::vector<GlobalVector*> const& x,
289 std::vector<GlobalVector*> const& x_prev, int const process_id)
290{
291 if (!_compensate_non_equilibrium_initial_residuum)
292 {
293 return;
294 }
295
296 INFO("Calculate non-equilibrium initial residuum.");
297
298 _equation_system->assemble(x, x_prev, process_id);
300 _equation_system->getResidual(*x[process_id], *x_prev[process_id], *_r_neq);
301
302 // Set the values of the selected entries of _r_neq, which are associated
303 // with the equations that do not need initial residual compensation, to
304 // zero.
305 const auto selected_global_indices =
306 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
307 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
308
309 _r_neq->set(selected_global_indices, zero_entries);
310
312}
313
315 std::vector<GlobalVector*>& x,
316 std::vector<GlobalVector*> const& x_prev,
317 std::function<void(int, std::vector<GlobalVector*> const&)> const&
318 postIterationCallback,
319 int const process_id)
320{
321 namespace LinAlg = MathLib::LinAlg;
322 auto& sys = *_equation_system;
323
325 auto& minus_delta_x =
328
329 bool error_norms_met = false;
330
331 // TODO be more efficient
332 // init minus_delta_x to the right size
333 LinAlg::copy(*x[process_id], minus_delta_x);
334
335 _convergence_criterion->preFirstIteration();
336
337 int iteration = 1;
338 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
339 {
340 BaseLib::RunTime timer_dirichlet;
341 double time_dirichlet = 0.0;
342
343 BaseLib::RunTime time_iteration;
344 time_iteration.start();
345
346 timer_dirichlet.start();
347 sys.computeKnownSolutions(*x[process_id], process_id);
348 sys.applyKnownSolutions(*x[process_id]);
349 time_dirichlet += timer_dirichlet.elapsed();
350
351 sys.preIteration(iteration, *x[process_id]);
352
353 BaseLib::RunTime time_assembly;
354 time_assembly.start();
355 try
356 {
357 sys.assemble(x, x_prev, process_id);
358 }
359 catch (AssemblyException const& e)
360 {
361 ERR("Abort nonlinear iteration. Repeating timestep. Reason: {:s}",
362 e.what());
363 error_norms_met = false;
364 iteration = _maxiter;
365 break;
366 }
367 sys.getResidual(*x[process_id], *x_prev[process_id], res);
368 sys.getJacobian(J);
369 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
370
371 // Subtract non-equilibrium initial residuum if set
372 if (_r_neq != nullptr)
373 LinAlg::axpy(res, -1, *_r_neq);
374
375 minus_delta_x.setZero();
376
377 timer_dirichlet.start();
378 sys.applyKnownSolutionsNewton(J, res, minus_delta_x);
379 time_dirichlet += timer_dirichlet.elapsed();
380 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
381
382 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
383 {
384 _convergence_criterion->checkResidual(res);
385 }
386
387 BaseLib::RunTime time_linear_solver;
388 time_linear_solver.start();
389 bool iteration_succeeded = _linear_solver.solve(J, res, minus_delta_x);
390 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
391
392 if (!iteration_succeeded)
393 {
394 ERR("Newton: The linear solver failed.");
395 }
396 else
397 {
398 // TODO could be solved in a better way
399 // cf.
400 // https://petsc.org/release/manualpages/Vec/VecWAXPY
401
402 // Copy pointers, replace the one for the given process id.
403 std::vector<GlobalVector*> x_new{x};
404 x_new[process_id] =
406 *x[process_id], _x_new_id);
407 LinAlg::axpy(*x_new[process_id], -_damping, minus_delta_x);
408
409 if (postIterationCallback)
410 {
411 postIterationCallback(iteration, x_new);
412 }
413
414 switch (sys.postIteration(*x_new[process_id]))
415 {
417 break;
419 ERR("Newton: The postIteration() hook reported a "
420 "non-recoverable error.");
421 iteration_succeeded = false;
422 break;
424 INFO(
425 "Newton: The postIteration() hook decided that this "
426 "iteration has to be repeated.");
427 // TODO introduce some onDestroy hook.
429 *x_new[process_id]);
430 continue; // That throws the iteration result away.
431 }
432
433 LinAlg::copy(*x_new[process_id],
434 *x[process_id]); // copy new solution to x
436 *x_new[process_id]);
437 }
438
439 if (!iteration_succeeded)
440 {
441 // Don't compute further error norms, but break here.
442 error_norms_met = false;
443 break;
444 }
445
446 if (sys.isLinear())
447 {
448 error_norms_met = true;
449 }
450 else
451 {
452 if (_convergence_criterion->hasDeltaXCheck())
453 {
454 // Note: x contains the new solution!
455 _convergence_criterion->checkDeltaX(minus_delta_x,
456 *x[process_id]);
457 }
458
459 error_norms_met = _convergence_criterion->isSatisfied();
460 }
461
462 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
463 time_iteration.elapsed());
464
465 if (error_norms_met)
466 {
467 break;
468 }
469
470 // Avoid increment of the 'iteration' if the error norms are not met,
471 // but maximum number of iterations is reached.
472 if (iteration >= _maxiter)
473 {
474 break;
475 }
476 }
477
478 if (iteration > _maxiter)
479 {
480 ERR("Newton: Could not solve the given nonlinear system within {:d} "
481 "iterations",
482 _maxiter);
483 }
484
488
489 return {error_norms_met, iteration};
490}
491
492std::pair<std::unique_ptr<NonlinearSolverBase>, NonlinearSolverTag>
494 BaseLib::ConfigTree const& config)
495{
497 auto const type = config.getConfigParameter<std::string>("type");
499 auto const max_iter = config.getConfigParameter<int>("max_iter");
500
501 if (type == "Picard")
502 {
503 auto const tag = NonlinearSolverTag::Picard;
504 using ConcreteNLS = NonlinearSolver<tag>;
505 return std::make_pair(
506 std::make_unique<ConcreteNLS>(linear_solver, max_iter), tag);
507 }
508 if (type == "Newton")
509 {
511 auto const damping = config.getConfigParameter<double>("damping", 1.0);
512 if (damping <= 0)
513 {
514 OGS_FATAL(
515 "The damping factor for the Newon method must be positive, got "
516 "{:g}.",
517 damping);
518 }
519 auto const tag = NonlinearSolverTag::Newton;
520 using ConcreteNLS = NonlinearSolver<tag>;
521 return std::make_pair(
522 std::make_unique<ConcreteNLS>(linear_solver, max_iter, damping),
523 tag);
524 }
525#ifdef USE_PETSC
526 if (boost::iequals(type, "PETScSNES"))
527 {
528 auto prefix =
530 config.getConfigParameter<std::string>("prefix", "");
531 auto const tag = NonlinearSolverTag::Newton;
532 using ConcreteNLS = PETScNonlinearSolver;
533 return std::make_pair(std::make_unique<ConcreteNLS>(
534 linear_solver, max_iter, std::move(prefix)),
535 tag);
536 }
537
538#endif
539 OGS_FATAL("Unsupported nonlinear solver type '{:s}'.", type.c_str());
540}
541
549
557
558} // 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.