OGS
NumLib::NonlinearSolver< NonlinearSolverTag::Newton > Class Referencefinal

Detailed Description

Find a solution to a nonlinear equation using the Newton-Raphson method.

Definition at line 76 of file NonlinearSolver.h.

#include <NonlinearSolver.h>

Inheritance diagram for NumLib::NonlinearSolver< NonlinearSolverTag::Newton >:
[legend]
Collaboration diagram for NumLib::NonlinearSolver< NonlinearSolverTag::Newton >:
[legend]

Public Types

using System = NonlinearSystem<NonlinearSolverTag::Newton>
 Type of the nonlinear equation system to be solved.
 

Public Member Functions

 NonlinearSolver (GlobalLinearSolver &linear_solver, int const maxiter, double const damping=1.0)
 
 ~NonlinearSolver ()
 
void setEquationSystem (System &eq, ConvergenceCriterion &conv_crit)
 
void calculateNonEquilibriumInitialResiduum (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
 
NonlinearSolverStatus solve (std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::function< void(int, std::vector< GlobalVector * > const &)> const &postIterationCallback, int const process_id) override
 
void compensateNonEquilibriumInitialResiduum (bool const value)
 
- Public Member Functions inherited from NumLib::NonlinearSolverBase
virtual ~NonlinearSolverBase ()=default
 

Private Attributes

GlobalLinearSolver_linear_solver
 
System_equation_system = nullptr
 
ConvergenceCriterion_convergence_criterion = nullptr
 
int const _maxiter
 maximum number of iterations
 
double const _damping
 
GlobalVector_r_neq = nullptr
 non-equilibrium initial residuum.
 
std::size_t _res_id = 0u
 ID of the residual vector.
 
std::size_t _J_id = 0u
 ID of the Jacobian matrix.
 
std::size_t _minus_delta_x_id = 0u
 ID of the \( -\Delta x\) vector.
 
std::size_t _x_new_id
 ID of the vector storing \( x - (-\Delta x) \).
 
std::size_t _r_neq_id = 0u
 
bool _compensate_non_equilibrium_initial_residuum = false
 

Member Typedef Documentation

◆ System

Type of the nonlinear equation system to be solved.

Definition at line 81 of file NonlinearSolver.h.

Constructor & Destructor Documentation

◆ NonlinearSolver()

NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::NonlinearSolver ( GlobalLinearSolver & linear_solver,
int const maxiter,
double const damping = 1.0 )
inlineexplicit

Constructs a new instance.

Parameters
linear_solverthe linear solver used by this nonlinear solver.
maxiterthe maximum number of iterations used to solve the equation.
dampingA positive damping factor.
See also
_damping

Definition at line 91 of file NonlinearSolver.h.

94 : _linear_solver(linear_solver), _maxiter(maxiter), _damping(damping)
95 {
96 }
int const _maxiter
maximum number of iterations

◆ ~NonlinearSolver()

Definition at line 556 of file NonlinearSolver.cpp.

557{
558 if (_r_neq != nullptr)
559 {
561 }
562}
GlobalVector * _r_neq
non-equilibrium initial residuum.
virtual void releaseVector(GlobalVector const &x)=0
static NUMLIB_EXPORT VectorProvider & provider

References NumLib::GlobalVectorProvider::provider, and NumLib::VectorProvider::releaseVector().

Member Function Documentation

◆ calculateNonEquilibriumInitialResiduum()

void NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::calculateNonEquilibriumInitialResiduum ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )
overridevirtual

Implements NumLib::NonlinearSolverBase.

Definition at line 292 of file NonlinearSolver.cpp.

296{
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 =
313 std::vector<double> zero_entries(selected_global_indices.size(), 0.0);
314
315 _r_neq->set(selected_global_indices, zero_entries);
316
318}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void set(IndexType rowId, double v)
set entry
Definition EigenVector.h:73
virtual std::vector< GlobalIndexType > getIndicesOfResiduumWithoutInitialCompensation() const =0
virtual void assemble(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)=0
virtual void getResidual(GlobalVector const &x, GlobalVector const &x_prev, GlobalVector &res) const =0
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170

References MathLib::LinAlg::finalizeAssembly(), NumLib::VectorProvider::getVector(), INFO(), and NumLib::GlobalVectorProvider::provider.

◆ compensateNonEquilibriumInitialResiduum()

void NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::compensateNonEquilibriumInitialResiduum ( bool const value)
inline

Definition at line 120 of file NonlinearSolver.h.

◆ setEquationSystem()

void NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::setEquationSystem ( System & eq,
ConvergenceCriterion & conv_crit )
inline

Set the nonlinear equation system that will be solved. TODO doc

Definition at line 102 of file NonlinearSolver.h.

103 {
104 _equation_system = &eq;
105 _convergence_criterion = &conv_crit;
106 }

◆ solve()

NonlinearSolverStatus NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve ( std::vector< GlobalVector * > & x,
std::vector< GlobalVector * > const & x_prev,
std::function< void(int, std::vector< GlobalVector * > const &)> const & postIterationCallback,
int const process_id )
overridevirtual

Assemble and solve the equation system.

Parameters
xin: the initial guess, out: the solution.
x_prevprevious time step solution.
postIterationCallbackcalled after each iteration if set.
process_idusually used in staggered schemes.
Return values
trueif the equation system could be solved
falseotherwise

Implements NumLib::NonlinearSolverBase.

Definition at line 320 of file NonlinearSolver.cpp.

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
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 {
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 {
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}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
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)
virtual void checkResidual(GlobalVector const &residual)=0
Check if the residual satisfies the convergence criterion.
virtual bool hasResidualCheck() const =0
virtual bool isSatisfied() const
Tell if the convergence criterion is satisfied.
virtual void checkDeltaX(GlobalVector const &minus_delta_x, GlobalVector const &x)=0
virtual bool hasDeltaXCheck() const =0
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
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.
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57
static NUMLIB_EXPORT MatrixProvider & provider

References MathLib::LinAlg::axpy(), MathLib::LinAlg::copy(), BaseLib::RunTime::elapsed(), ERR(), NumLib::FAILURE, NumLib::MatrixProvider::getMatrix(), NumLib::VectorProvider::getVector(), INFO(), NumLib::GlobalVectorProvider::provider, NumLib::GlobalMatrixProvider::provider, NumLib::MatrixProvider::releaseMatrix(), NumLib::VectorProvider::releaseVector(), NumLib::REPEAT_ITERATION, BaseLib::RunTime::start(), and NumLib::SUCCESS.

Member Data Documentation

◆ _compensate_non_equilibrium_initial_residuum

bool NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_compensate_non_equilibrium_initial_residuum = false
private

Enables computation of the non-equilibrium initial residuum \( r_{\rm neq} \) before the first time step. The forces are zero if the external forces are in equilibrium with the initial state/initial conditions. During the simulation the new residuum reads \( \tilde r = r - r_{\rm neq} \).

Definition at line 153 of file NonlinearSolver.h.

◆ _convergence_criterion

ConvergenceCriterion* NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_convergence_criterion = nullptr
private

Definition at line 130 of file NonlinearSolver.h.

◆ _damping

double const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_damping
private

A positive damping factor. The default value 1.0 gives a non-damped Newton method. Common values are in the range 0.5 to 0.7 for somewhat conservative method and seldom become smaller than 0.2 for very conservative approach.

Definition at line 137 of file NonlinearSolver.h.

◆ _equation_system

System* NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_equation_system = nullptr
private

Definition at line 127 of file NonlinearSolver.h.

◆ _J_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_J_id = 0u
private

ID of the Jacobian matrix.

Definition at line 141 of file NonlinearSolver.h.

◆ _linear_solver

Definition at line 126 of file NonlinearSolver.h.

◆ _maxiter

int const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_maxiter
private

maximum number of iterations

Definition at line 131 of file NonlinearSolver.h.

◆ _minus_delta_x_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_minus_delta_x_id = 0u
private

ID of the \( -\Delta x\) vector.

Definition at line 142 of file NonlinearSolver.h.

◆ _r_neq

non-equilibrium initial residuum.

Definition at line 139 of file NonlinearSolver.h.

◆ _r_neq_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_r_neq_id = 0u
private

ID of the non-equilibrium initial residuum vector.

Definition at line 145 of file NonlinearSolver.h.

◆ _res_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_res_id = 0u
private

ID of the residual vector.

Definition at line 140 of file NonlinearSolver.h.

◆ _x_new_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_x_new_id
private
Initial value:
=
0u

ID of the vector storing \( x - (-\Delta x) \).

Definition at line 143 of file NonlinearSolver.h.


The documentation for this class was generated from the following files: