This class defines a Newton solver for nonlinear systems of equations of the form \(F(x) = 0\).
More...
#include <NewtonSolver.h>
|
| NewtonSolver (MPI_Comm comm) |
| Create nonlinear solver.
|
|
| NewtonSolver (NewtonSolver &&solver)=delete |
|
| NewtonSolver (const NewtonSolver &solver)=delete |
|
NewtonSolver & | operator= (const NewtonSolver &solver)=delete |
|
NewtonSolver & | operator= (const NewtonSolver &&solver)=delete |
|
| ~NewtonSolver () |
| Destructor.
|
|
void | setF (const 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 (const 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 (const std::function< void(const Vec, Mat)> &P, Mat Pmat) |
| Set the function for computing the preconditioner matrix (optional)
|
|
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_.
|
|
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_.
|
|
void | set_form (const 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 (const 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 (const 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 |
| Return number of Krylov iterations elapsed since solve started.
|
|
double | residual () const |
| Return current residual.
|
|
double | residual0 () const |
| Return initial residual.
|
|
MPI_Comm | comm () const |
| Return MPI communicator.
|
|
|
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.
|
|
This class defines a Newton solver for nonlinear systems of equations of the form \(F(x) = 0\).
◆ NewtonSolver()
Create nonlinear solver.
- Parameters
-
[in] | comm | The MPI communicator for the solver |
◆ get_krylov_solver() [1/2]
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]
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()
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 |
Return number of Krylov iterations elapsed since solve started.
- Returns
- Number of iterations.
◆ residual()
double residual |
( |
| ) |
const |
Return current residual.
- Returns
- Current residual
◆ residual0()
double residual0 |
( |
| ) |
const |
Return initial residual.
- Returns
- Initial residual
◆ set_convergence_check()
void set_convergence_check |
( |
const 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] | c | The function that tests for convergence |
◆ set_form()
void set_form |
( |
const 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] | form | The function to call. It takes the latest solution vector x as an argument |
◆ set_update()
void set_update |
( |
const 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] | update | The function that updates the solution |
◆ setF()
void setF |
( |
const 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] | F | Function to compute the residual vector b (x, b) |
[in] | b | The vector to assemble to residual into |
◆ setJ()
void setJ |
( |
const 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] | J | Function to compute the Jacobian matrix b (x, A) |
[in] | Jmat | The matrix to assemble the Jacobian into |
◆ setP()
void setP |
( |
const std::function< void(const Vec, Mat)> & |
P, |
|
|
Mat |
Pmat |
|
) |
| |
Set the function for computing the preconditioner matrix (optional)
- Parameters
-
[in] | P | Function to compute the preconditioner matrix b (x, P) |
[in] | Pmat | The 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
-
- Returns
- (number of Newton iterations, whether iteration converged)
The documentation for this class was generated from the following files:
- /__w/dolfinx/dolfinx/cpp/dolfinx/nls/NewtonSolver.h
- /__w/dolfinx/dolfinx/cpp/dolfinx/nls/NewtonSolver.cpp