9#include <dolfinx/common/MPI.h>
10#include <dolfinx/la/petsc.h>
57 void setF(
const std::function<
void(
const Vec, Vec)>& F, Vec b);
63 void setJ(
const std::function<
void(
const Vec, Mat)>& J, Mat Jmat);
68 void setP(
const std::function<
void(
const Vec, Mat)>& P, Mat Pmat);
86 void set_form(
const std::function<
void(Vec)>& form);
98 const Vec, Vec)>& update);
105 std::pair<int, bool>
solve(Vec x);
126 MPI_Comm
comm()
const;
154 std::function<void(
const Vec x, Vec b)> _fnF;
159 std::function<void(
const Vec x, Mat J)> _fnJ;
164 std::function<void(
const Vec x, Mat P)> _fnP;
168 std::function<void(
const Vec x)> _system;
174 Mat _matJ =
nullptr, _matP =
nullptr;
177 std::function<std::pair<double, bool>(
const NewtonSolver& solver,
182 std::function<void(
const NewtonSolver& solver,
const Vec dx, Vec x)>
186 int _krylov_iterations;
192 double _residual, _residual0;
A duplicate MPI communicator and manage lifetime of the communicator.
Definition MPI.h:43
This class implements Krylov methods for linear systems of the form Ax = b. It is a wrapper for the K...
Definition petsc.h:451
This class defines a Newton solver for nonlinear systems of equations of the form .
Definition NewtonSolver.h:32
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.
Definition NewtonSolver.cpp:141
double relaxation_parameter
Relaxation parameter.
Definition NewtonSolver.h:148
int krylov_iterations() const
Return number of Krylov iterations elapsed since solve started.
Definition NewtonSolver.cpp:279
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 ...
Definition NewtonSolver.cpp:129
double rtol
Relative tolerance.
Definition NewtonSolver.h:132
~NewtonSolver()
Destructor.
Definition NewtonSolver.cpp:82
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.
Definition NewtonSolver.cpp:94
double atol
Absolute tolerance.
Definition NewtonSolver.h:135
double residual() const
Return current residual.
Definition NewtonSolver.cpp:286
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.
Definition NewtonSolver.cpp:102
double residual0() const
Return initial residual.
Definition NewtonSolver.cpp:288
std::string convergence_criterion
Convergence criterion.
Definition NewtonSolver.h:139
int iteration() const
The number of Newton iterations. It can can called by functions that check for convergence during a s...
Definition NewtonSolver.cpp:284
bool report
Monitor convergence.
Definition NewtonSolver.h:142
bool error_on_nonconvergence
Throw error if solver fails to converge.
Definition NewtonSolver.h:145
MPI_Comm comm() const
Return MPI communicator.
Definition NewtonSolver.cpp:290
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 p...
Definition NewtonSolver.cpp:119
std::pair< int, bool > solve(Vec x)
Solve the nonlinear problem for given and Jacobian .
Definition NewtonSolver.cpp:148
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.
Definition NewtonSolver.cpp:134
int max_it
Maximum number of iterations.
Definition NewtonSolver.h:129
void setP(const std::function< void(const Vec, Mat)> &P, Mat Pmat)
Set the function for computing the preconditioner matrix (optional)
Definition NewtonSolver.cpp:110
Top-level namespace.
Definition defines.h:12