DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
NewtonSolver Class Reference

#include <NewtonSolver.h>

Public Member Functions

 NewtonSolver (MPI_Comm comm)
 
 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 and the vector to the assemble the residual into.
 
void setJ (std::function< void(const Vec, Mat)> J, Mat Jmat)
 Set the function for computing the Jacobian (dF/dx) and the matrix to assemble the residual into.
 
void setP (std::function< void(const Vec, Mat)> P, Mat Pmat)
 Set the function for computing the preconditioner matrix (optional).
 
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)
 Set function that is called after each Newton iteration to update the solution.
 
std::pair< int, bool > solve (Vec x)
 Solve the nonlinear problem \(`F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\).
 
int iteration () const
 The number of Newton iterations. It can can called by functions that check for convergence during a solve.
 
int krylov_iterations () const
 Get number of Krylov iterations elapsed since solve started.
 
double residual () const
 Get current residual.
 
double residual0 () const
 Return initial residual.
 
MPI_Comm comm () const
 Get MPI communicator.
 

Public Attributes

int max_it = 50
 Maximum number of iterations.
 
double rtol = 1e-9
 Relative tolerance.
 
double atol = 1e-10
 Absolute 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

This class defines a Newton solver for nonlinear systems of equations of the form \(F(x) = 0\).

Constructor & Destructor Documentation

◆ NewtonSolver()

NewtonSolver ( MPI_Comm comm)
explicit

Create nonlinear solver

Parameters
[in]commThe MPI 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

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

Returns
The number of Newton iterations performed

◆ krylov_iterations()

int krylov_iterations ( ) const

Get number of Krylov iterations elapsed since solve started.

Returns
Number of iterations.

◆ residual()

double residual ( ) const

Get current residual.

Returns
Current residual

◆ residual0()

double residual0 ( ) const

Return 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]cThe function 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]formThe function 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)

Set function that is called after each Newton iteration to update the solution.

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 and the vector to the assemble the residual into.

Parameters
[in]FFunction to compute the residual vector b (x, b)
[in]bThe vector to assemble to residual into

◆ setJ()

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

Set the function for computing the Jacobian (dF/dx) and the matrix to assemble the residual into.

Parameters
[in]JFunction to compute the Jacobian matrix b (x, A)
[in]JmatThe matrix 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 (optional).

Parameters
[in]PFunction to compute the preconditioner matrix b (x, P)
[in]PmatThe matrix to assemble the preconditioner into

◆ solve()

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

Solve the nonlinear problem \(`F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\).

Parameters
[in,out]xThe vector
Returns
(number of Newton iterations, whether iteration converged)

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