DOLFINx 0.11.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
NewtonSolver Class Reference

A Newton solver for nonlinear systems of equations of the form \(F(x) = 0\). More...

#include <NewtonSolver.h>

Public Member Functions

 NewtonSolver (MPI_Comm comm)
 Create nonlinear solver.
 NewtonSolver (NewtonSolver &&solver)=default
 Move constructor.
 NewtonSolver (const NewtonSolver &solver)=delete
NewtonSolveroperator= (NewtonSolver &&solver)=default
 Move assignment constructor.
NewtonSolveroperator= (const NewtonSolver &solver)=delete
 ~NewtonSolver ()
 Destructor.
void setF (std::function< void(const Vec, Vec)> F, Vec b)
 Set the function for computing the residual \(F(x) = 0\) and the vector to assemble the residual into.
void setJ (std::function< void(const Vec, Mat)> J, Mat Jmat)
 Set the function for computing the Jacobian \(J:=dF/dx\) and the matrix to assemble the Jacobian into.
void setP (std::function< void(const Vec, Mat)> P, Mat Pmat)
 Set the function for computing the preconditioner matrix.
const la::petsc::KrylovSolverget_krylov_solver () const
 Get the internal Krylov solver used to solve for the Newton updates (const version).
la::petsc::KrylovSolverget_krylov_solver ()
 Get the internal Krylov solver used to solve for the Newton updates (non-const version).
void set_form (std::function< void(Vec)> form)
 Set the function that is called before the residual or Jacobian are computed. It is commonly used to update ghost values.
void set_convergence_check (std::function< std::pair< double, bool >(const NewtonSolver &, const Vec)> c)
 Set function that is called at the end of each Newton iteration to test for convergence.
void set_update (std::function< void(const NewtonSolver &solver, const Vec, Vec)> update)
 Optional set function that is called after each inner solve for the Newton increment to update the solution.
std::pair< int, bool > solve (Vec x)
 Solve the nonlinear problem.
int iteration () const
 Get number of Newton iterations. It can be called by functions that check for convergence during a solve.
int krylov_iterations () const
 Number of Krylov iterations elapsed since solve started.
double residual () const
 Get current residual.
double residual0 () const
 Get initial residual.
MPI_Comm comm () const
 Get MPI communicator.

Public Attributes

int max_it = 50
 Maximum number of iterations.
double rtol = 1e-9
 Relative convergence tolerance.
double atol = 1e-10
 Absolute convergence tolerance.
std::string convergence_criterion = "residual"
 Convergence criterion.
bool report = true
 Monitor convergence.
bool error_on_nonconvergence = true
 Throw error if solver fails to converge.
double relaxation_parameter = 1.0
 Relaxation parameter.

Detailed Description

A Newton solver for nonlinear systems of equations of the form \(F(x) = 0\).

It solves

\[ \left. \frac{dF}{dx} \right|_{x} \Delta x = F(x) \]

with default update \(x \leftarrow x - \Delta x\).

It relies on PETSc for linear algebra backends.

Deprecated
The generic NewtonSolver will be removed in a future release. It is recommended to build your own Newton Solver, or use SNES from PETSc.

Constructor & Destructor Documentation

◆ NewtonSolver()

NewtonSolver ( MPI_Comm comm)
explicit

Create nonlinear solver.

Parameters
[in]commMPI communicator for the solver

Member Function Documentation

◆ get_krylov_solver() [1/2]

la::petsc::KrylovSolver & get_krylov_solver ( )

Get the internal Krylov solver used to solve for the Newton updates (non-const version).

The Krylov solver prefix is nls_solve_.

Returns
The Krylov solver

◆ get_krylov_solver() [2/2]

const la::petsc::KrylovSolver & get_krylov_solver ( ) const

Get the internal Krylov solver used to solve for the Newton updates (const version).

The Krylov solver prefix is nls_solve_.

Returns
The Krylov solver

◆ iteration()

int iteration ( ) const

Get number of Newton iterations. It can be called by functions that check for convergence during a solve.

Returns
Number of Newton iterations performed.

◆ krylov_iterations()

int krylov_iterations ( ) const

Number of Krylov iterations elapsed since solve started.

Returns
Number of Krylov iterations.

◆ residual()

double residual ( ) const

Get current residual.

Returns
Current residual.

◆ residual0()

double residual0 ( ) const

Get initial residual.

Returns
Initial residual.

◆ set_convergence_check()

void set_convergence_check ( std::function< std::pair< double, bool >(const NewtonSolver &, const Vec)> c)

Set function that is called at the end of each Newton iteration to test for convergence.

Parameters
[in]cFunction that tests for convergence

◆ set_form()

void set_form ( std::function< void(Vec)> form)

Set the function that is called before the residual or Jacobian are computed. It is commonly used to update ghost values.

Parameters
[in]formFunction to call. It takes the (latest) solution vector x as an argument.

◆ set_update()

void set_update ( std::function< void(const NewtonSolver &solver, const Vec, Vec)> update)

Optional set function that is called after each inner solve for the Newton increment to update the solution.

The function update takes this, the Newton increment dx, and the vector x from the start of the Newton iteration.

By default, the update is x <- x - dx

Parameters
[in]updateThe function that updates the solution

◆ setF()

void setF ( std::function< void(const Vec, Vec)> F,
Vec b )

Set the function for computing the residual \(F(x) = 0\) and the vector to assemble the residual into.

Parameters
[in]FFunction to compute/assemble the residual vector b. The first argument to the function is the solution vector x and the second is the vector b to assemble into.
[in]bVector to assemble to residual into.

◆ setJ()

void setJ ( std::function< void(const Vec, Mat)> J,
Mat Jmat )

Set the function for computing the Jacobian \(J:=dF/dx\) and the matrix to assemble the Jacobian into.

Parameters
[in]JFunction to compute the Jacobian matrix Jmat. The first argument to the function is the solution vector x and the second is the matrix to assemble into.
[in]JmatMatrix to assemble the Jacobian into.

◆ setP()

void setP ( std::function< void(const Vec, Mat)> P,
Mat Pmat )

Set the function for computing the preconditioner matrix.

It is optional to set the preconditioner matrix. By default the solver will use the Jacobian matrix.

Parameters
[in]PFunction to compute the preconditioner matrix Pmat. The first argument to the function is the solution vector x and the second is the matrix to assemble into.
[in]PmatMatrix to assemble the preconditioner into.

◆ solve()

std::pair< int, bool > solve ( Vec x)

Solve the nonlinear problem.

Parameters
[in,out]xThe solution vector. It should be set the initial solution guess.
Returns
(number of Newton iterations, whether iteration converged)

Member Data Documentation

◆ convergence_criterion

std::string convergence_criterion = "residual"

Convergence criterion.

Todo
change to string to enum.

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