|
DOLFINx 0.11.0.0
DOLFINx C++ interface
|
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 | |
| NewtonSolver & | operator= (NewtonSolver &&solver)=default |
| Move assignment constructor. | |
| NewtonSolver & | operator= (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::KrylovSolver & | get_krylov_solver () const |
| Get the internal Krylov solver used to solve for the Newton updates (const version). | |
| la::petsc::KrylovSolver & | get_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. | |
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.
|
explicit |
| 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_.
| 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_.
| 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.
| 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.
| [in] | c | Function that tests for convergence |
| 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.
| [in] | form | Function to call. It takes the (latest) solution vector x as an argument. |
| 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
| [in] | update | The function that updates the solution |
| 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.
| [in] | F | Function 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] | b | Vector to assemble to 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.
| [in] | J | Function 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] | Jmat | 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.
It is optional to set the preconditioner matrix. By default the solver will use the Jacobian matrix.
| [in] | P | Function 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] | Pmat | Matrix to assemble the preconditioner into. |
| std::pair< int, bool > solve | ( | Vec | x | ) |
Solve the nonlinear problem.
| [in,out] | x | The solution vector. It should be set the initial solution guess. |
| std::string convergence_criterion = "residual" |
Convergence criterion.