SpECTRE
v2024.12.16
|
A Newton-Raphson correction scheme for nonlinear systems of equations \(A_\mathrm{nonlinear}(x)=b\). More...
#include <NewtonRaphson.hpp>
Public Types | |
using | fields_tag = FieldsTag |
using | options_group = OptionsGroup |
using | source_tag = SourceTag |
using | operand_tag = fields_tag |
using | linear_solver_fields_tag = db::add_tag_prefix< NonlinearSolver::Tags::Correction, fields_tag > |
using | linear_solver_source_tag = db::add_tag_prefix< NonlinearSolver::Tags::Residual, fields_tag > |
using | component_list = tmpl::list< detail::ResidualMonitor< Metavariables, FieldsTag, OptionsGroup > > |
using | initialize_element = detail::InitializeElement< FieldsTag, OptionsGroup, SourceTag > |
using | register_element = tmpl::list<> |
using | amr_projectors = initialize_element |
template<typename ApplyNonlinearOperator , typename SolveLinearization , typename ObserveActions = tmpl::list<>, typename Label = OptionsGroup> | |
using | solve = tmpl::list< detail::PrepareSolve< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, ApplyNonlinearOperator, ObserveActions, detail::SendInitialResidualMagnitude< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, detail::ReceiveInitialHasConverged< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, detail::PrepareStep< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, SolveLinearization, detail::PerformStep< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, ApplyNonlinearOperator, detail::ContributeToResidualMagnitudeReduction< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, detail::Globalize< FieldsTag, OptionsGroup, Label, ArraySectionIdTag >, ObserveActions, detail::CompleteStep< FieldsTag, OptionsGroup, Label, ArraySectionIdTag > > |
A Newton-Raphson correction scheme for nonlinear systems of equations \(A_\mathrm{nonlinear}(x)=b\).
We can use a correction scheme to solve a nonlinear problem \(A_\mathrm{nonlinear}(x)=b\) by repeatedly solving a linearization of it. A Newton-Raphson scheme iteratively refines an initial guess \(x_0\) by repeatedly solving the linearized problem
\begin{equation} \frac{\delta A_\mathrm{nonlinear}}{\delta x}(x_k) \delta x_k = b-A_\mathrm{nonlinear}(x_k) \equiv r_k \end{equation}
for the correction \(\delta x_k\) and then updating the solution as \(x_{k+1}=x_k + \delta x_k\).
The operations we need to supply to the algorithm are the nonlinear operator \(A_\mathrm{nonlinear}(x)\) and a linear solver for the linearized problem \(\frac{\delta A_\mathrm{nonlinear}}{\delta x}(x_k) \delta x_k = r_k\). Each step of the algorithm expects that \(A_\mathrm{nonlinear}(x)\) is computed and stored in the DataBox as db::add_tag_prefix<NonlinearSolver::Tags::OperatorAppliedTo, FieldsTag>
. To perform a solve, add the solve
action list to an array parallel component. Pass the actions that compute \(A_\mathrm{nonlinear}(x)\) as the first template parameter to solve
. As the second template parameter, pass the action list that performs a linear solve of the linearized operator \(\frac{\delta A_\mathrm{nonlinear}}{\delta x}(x)\) for the field \(\delta x\) (the linear_solver_fields_tag
) sourced by \(r\) (the linear_solver_source_tag
). You will find suitable iterative linear solvers in the LinearSolver
namespace. The third template parameter allows you to pass any further actions you wish to run in each step of the algorithm (such as observations). If you add the solve
action list multiple times, use the fourth template parameter to label each solve with a different type.
NonlinearSolver::OptionTags::SufficientDecrease
for details on the sufficient decrease condition), the step length is reduced until the residual sufficiently decreases or the maximum number of globalization steps is reached. The reduction in step length is determined by the minimum of a quadratic (first globalization step) or cubic (subsequent globalization steps) polynomial interpolation according to Algorithm A6.1.3 in [50] (p. 325) (see NonlinearSolver::newton_raphson::next_step_length
). Alternatives to a line-search globalization, such as a trust-region globalization or more sophisticated nonlinear preconditioning techniques (see e.g. [32] for an overview), are not currently implemented.Parallel::Section
). Set the ArraySectionIdTag
template parameter to restrict the solver to elements in that section. Only a single section must be associated with the ArraySectionIdTag
. The default is void
, which means running over all elements in the array. Note that the actions in the SolveLinearization
list passed to solve
will not be restricted to run only on section elements, so all elements in the array may participate in preconditioning (see LinearSolver::multigrid::Multigrid).