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

Detailed Description

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

Definition at line 69 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, int const recompute_jacobian=1, double const damping=1.0, std::optional< double > const damping_reduction=std::nullopt)
 ~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, bool, 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
int const _recompute_jacobian = 1
double const _damping
std::optional< double > const _damping_reduction
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 74 of file NonlinearSolver.h.

Constructor & Destructor Documentation

◆ NonlinearSolver()

NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::NonlinearSolver ( GlobalLinearSolver & linear_solver,
int const maxiter,
int const recompute_jacobian = 1,
double const damping = 1.0,
std::optional< double > const damping_reduction = std::nullopt )
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.
recompute_jacobianrecompute jacobian for the follow-up steps
dampingA positive damping factor.
See also
_damping
Parameters
damping_reductionA parameter that reduces damping
See also
_damping_reduction

Definition at line 87 of file NonlinearSolver.h.

References _damping, _damping_reduction, _linear_solver, _maxiter, and _recompute_jacobian.

Referenced by ~NonlinearSolver(), calculateNonEquilibriumInitialResiduum(), and solve().

◆ ~NonlinearSolver()

Definition at line 630 of file NonlinearSolver.cpp.

631{
632 if (_r_neq != nullptr)
633 {
635 }
636}
GlobalVector * _r_neq
non-equilibrium initial residuum.

References NonlinearSolver(), _r_neq, and NumLib::GlobalVectorProvider::provider.

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 315 of file NonlinearSolver.cpp.

319{
321 {
322 return;
323 }
324
325 INFO("Calculate non-equilibrium initial residuum.");
326
327 _equation_system->assemble(x, x_prev, process_id);
330
331 // Set the values of the selected entries of _r_neq, which are associated
332 // with the equations that do not need initial residual compensation, to
333 // zero.
334 const auto selected_global_indices =
335 _equation_system->getIndicesOfResiduumWithoutInitialCompensation();
337
339 _equation_system->setReleaseNodalForces(_r_neq, process_id);
340
342}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:191

References NonlinearSolver(), _compensate_non_equilibrium_initial_residuum, _equation_system, _r_neq, _r_neq_id, MathLib::LinAlg::finalizeAssembly(), INFO(), and NumLib::GlobalVectorProvider::provider.

◆ compensateNonEquilibriumInitialResiduum()

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

◆ 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 105 of file NonlinearSolver.h.

References _convergence_criterion, and _equation_system.

◆ solve()

NonlinearSolverStatus NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve ( std::vector< GlobalVector * > & x,
std::vector< GlobalVector * > const & x_prev,
std::function< void(int, bool, 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 344 of file NonlinearSolver.cpp.

350{
351 namespace LinAlg = MathLib::LinAlg;
352 auto& sys = *_equation_system;
353
355 auto& minus_delta_x =
358
359 bool error_norms_met = false;
360
361 // TODO be more efficient
362 // init minus_delta_x to the right size
364
365 _convergence_criterion->preFirstIteration();
366
367 int iteration = 1;
368#if !defined(USE_PETSC) && !defined(USE_LIS)
370#endif
371 for (; iteration <= _maxiter; ++iteration, _convergence_criterion->reset())
372 {
374 double time_dirichlet = 0.0;
375
377 INFO("Iteration #{:d} started.", iteration);
378 time_iteration.start();
379
380 timer_dirichlet.start();
381 sys.computeKnownSolutions(*x[process_id], process_id);
382 time_dirichlet += timer_dirichlet.elapsed();
383
384 sys.preIteration(iteration, *x[process_id]);
385
387 time_assembly.start();
388 bool mpi_rank_assembly_ok = true;
389 try
390 {
391 sys.assemble(x, x_prev, process_id);
392 }
393 catch (AssemblyException const& e)
394 {
395 ERR("Abort nonlinear iteration. Repeating timestep. Reason: {:s}",
396 e.what());
397 error_norms_met = false;
399 mpi_rank_assembly_ok = false;
400 }
402 {
403 break;
404 }
405 sys.getResidual(*x[process_id], *x_prev[process_id], res);
406 sys.getJacobian(J);
407 INFO("[time] Assembly took {:g} s.", time_assembly.elapsed());
408
409 // Subtract non-equilibrium initial residuum if set
410 if (_r_neq != nullptr)
411 LinAlg::axpy(res, -1, *_r_neq);
412
413 minus_delta_x.setZero();
414
415 timer_dirichlet.start();
416 sys.applyKnownSolutionsNewton(J, res, *x[process_id], minus_delta_x);
417 time_dirichlet += timer_dirichlet.elapsed();
418 INFO("[time] Applying Dirichlet BCs took {:g} s.", time_dirichlet);
419
420 if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
421 {
422 _convergence_criterion->checkResidual(res);
423 }
424
426 time_linear_solver.start();
427#if !defined(USE_PETSC) && !defined(USE_LIS)
430 {
435 }
436
437 bool iteration_succeeded = false;
439 {
440 ERR("Newton: The linear solver failed in the compute() step.");
441 }
442 else
443 {
445 }
446#else
448#endif
449 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
450
452 {
453 ERR("Newton: The linear solver failed.");
454 }
455 else
456 {
457 // TODO could be solved in a better way
458 // cf.
459 // https://petsc.org/release/manualpages/Vec/VecWAXPY
460
461 // Copy pointers, replace the one for the given process id.
466 double damping = _damping;
467 if (_convergence_criterion->hasNonNegativeDamping())
468 {
469 damping = _convergence_criterion->getDampingFactor(
471 }
473 {
474 damping =
475 damping +
476 (1 - damping) *
478 }
480
482 {
484 }
485
486 switch (sys.postIteration(*x_new[process_id]))
487 {
489 break;
491 ERR("Newton: The postIteration() hook reported a "
492 "non-recoverable error.");
493 iteration_succeeded = false;
494 break;
496 INFO(
497 "Newton: The postIteration() hook decided that this "
498 "iteration has to be repeated.");
499 // TODO introduce some onDestroy hook.
501 *x_new[process_id]);
502 continue; // That throws the iteration result away.
503 }
504
506 *x[process_id]); // copy new solution to x
508 *x_new[process_id]);
509 }
510
512 {
513 // Don't compute further error norms, but break here.
514 error_norms_met = false;
515 break;
516 }
517
518 if (sys.isLinear())
519 {
520 error_norms_met = true;
521 }
522 else
523 {
524 if (_convergence_criterion->hasDeltaXCheck())
525 {
526 // Note: x contains the new solution!
528 *x[process_id]);
529 }
530
532 }
533
534 INFO("[time] Iteration #{:d} took {:g} s.", iteration,
535 time_iteration.elapsed());
536
537 if (error_norms_met)
538 {
539 break;
540 }
541
542 // Avoid increment of the 'iteration' if the error norms are not met,
543 // but maximum number of iterations is reached.
544 if (iteration >= _maxiter)
545 {
546 break;
547 }
548 }
549
550 if (iteration > _maxiter)
551 {
552 ERR("Newton: Could not solve the given nonlinear system within {:d} "
553 "iterations",
554 _maxiter);
555 }
556
560
561 return {error_norms_met, iteration};
562}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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:30
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50

References NonlinearSolver(), _convergence_criterion, _damping, _damping_reduction, _equation_system, _J_id, _linear_solver, _maxiter, _minus_delta_x_id, _r_neq, _recompute_jacobian, _res_id, _x_new_id, BaseLib::MPI::anyOf(), MathLib::LinAlg::axpy(), MathLib::LinAlg::copy(), BaseLib::RunTime::elapsed(), ERR(), NumLib::FAILURE, INFO(), NumLib::GlobalMatrixProvider::provider, NumLib::GlobalVectorProvider::provider, MathLib::RECOMPUTE_AND_STORE, NumLib::REPEAT_ITERATION, MathLib::REUSE, 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 163 of file NonlinearSolver.h.

Referenced by calculateNonEquilibriumInitialResiduum(), and compensateNonEquilibriumInitialResiduum().

◆ _convergence_criterion

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

Definition at line 133 of file NonlinearSolver.h.

Referenced by setEquationSystem(), and solve().

◆ _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 142 of file NonlinearSolver.h.

Referenced by NonlinearSolver(), and solve().

◆ _damping_reduction

std::optional<double> const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_damping_reduction
private

A parameter that reduces the damping by iteration / _damping_reduction If the expression is smaller than zero, then damping is the same for all iterations. If the expression is bigger than one it is clamped by one.

Definition at line 147 of file NonlinearSolver.h.

Referenced by NonlinearSolver(), and solve().

◆ _equation_system

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

◆ _J_id

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

ID of the Jacobian matrix.

Definition at line 151 of file NonlinearSolver.h.

Referenced by solve().

◆ _linear_solver

Definition at line 129 of file NonlinearSolver.h.

Referenced by NonlinearSolver(), and solve().

◆ _maxiter

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

maximum number of iterations

Definition at line 134 of file NonlinearSolver.h.

Referenced by NonlinearSolver(), and solve().

◆ _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 152 of file NonlinearSolver.h.

Referenced by solve().

◆ _r_neq

non-equilibrium initial residuum.

Definition at line 149 of file NonlinearSolver.h.

Referenced by ~NonlinearSolver(), calculateNonEquilibriumInitialResiduum(), and solve().

◆ _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 155 of file NonlinearSolver.h.

Referenced by calculateNonEquilibriumInitialResiduum().

◆ _recompute_jacobian

int const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_recompute_jacobian = 1
private

Definition at line 136 of file NonlinearSolver.h.

Referenced by NonlinearSolver(), and solve().

◆ _res_id

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

ID of the residual vector.

Definition at line 150 of file NonlinearSolver.h.

Referenced by solve().

◆ _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 153 of file NonlinearSolver.h.

Referenced by solve().


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