DOLFINx 0.9.0
DOLFINx C++ interface
|
#include <NewtonSolver.h>
Public Member Functions | |
NewtonSolver (MPI_Comm comm) | |
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 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::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) |
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. | |
This class defines a Newton solver for nonlinear systems of equations of the form \(F(x) = 0\).
|
explicit |
Create nonlinear solver
[in] | comm | The MPI communicator for the solver |
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 |
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.
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 | The 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 | The 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 | ) |
Set function that is called after each Newton iteration to update the solution.
[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 and the vector to the assemble the residual into.
[in] | F | Function to compute the residual vector b (x, b) |
[in] | b | The vector to assemble to 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.
[in] | J | Function to compute the Jacobian matrix b (x, A) |
[in] | Jmat | 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 (optional).
[in] | P | Function to compute the preconditioner matrix b (x, P) |
[in] | Pmat | The matrix to assemble the preconditioner into |
std::pair< int, bool > solve | ( | Vec | x | ) |
Solve the nonlinear problem \(`F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\).
[in,out] | x | The vector |