Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d1/d88/NewtonSolver_8h_source.html
DOLFINx  0.5.1
DOLFINx C++ interface
NewtonSolver.h
1 // Copyright (C) 2005-2021 Garth N. Wells
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <dolfinx/common/MPI.h>
10 #include <dolfinx/la/petsc.h>
11 #include <functional>
12 #include <memory>
13 #include <petscmat.h>
14 #include <petscvec.h>
15 #include <utility>
16 
17 namespace dolfinx
18 {
19 
20 namespace la::petsc
21 {
22 class KrylovSolver;
23 } // namespace la::petsc
24 
25 namespace nls::petsc
26 {
27 
30 
32 {
33 public:
36  explicit NewtonSolver(MPI_Comm comm);
37 
38  // Move constructor (deleted)
39  NewtonSolver(NewtonSolver&& solver) = delete;
40 
41  // Copy constructor (deleted)
42  NewtonSolver(const NewtonSolver& solver) = delete;
43 
44  // Assignment operator (deleted)
45  NewtonSolver& operator=(const NewtonSolver& solver) = delete;
46 
47  // Move assignment constructor (deleted)
48  NewtonSolver& operator=(const NewtonSolver&& solver) = delete;
49 
51  ~NewtonSolver();
52 
57  void setF(const std::function<void(const Vec, Vec)>& F, Vec b);
58 
63  void setJ(const std::function<void(const Vec, Mat)>& J, Mat Jmat);
64 
68  void setP(const std::function<void(const Vec, Mat)>& P, Mat Pmat);
69 
75 
81 
86  void set_form(const std::function<void(Vec)>& form);
87 
91  void set_convergence_check(const std::function<std::pair<double, bool>(
92  const NewtonSolver&, const Vec)>& c);
93 
97  void set_update(const std::function<void(const NewtonSolver& solver,
98  const Vec, Vec)>& update);
99 
105  std::pair<int, bool> solve(Vec x);
106 
110  int iteration() const;
111 
115  int krylov_iterations() const;
116 
119  double residual() const;
120 
123  double residual0() const;
124 
126  MPI_Comm comm() const;
127 
129  int max_it = 50;
130 
132  double rtol = 1e-9;
133 
135  double atol = 1e-10;
136 
137  // FIXME: change to string to enum
139  std::string convergence_criterion = "residual";
140 
142  bool report = true;
143 
146 
148  double relaxation_parameter = 1.0;
149 
150 private:
151  // Function for computing the residual vector. The first argument is
152  // the latest solution vector x and the second argument is the
153  // residual vector.
154  std::function<void(const Vec x, Vec b)> _fnF;
155 
156  // Function for computing the Jacobian matrix operator. The first
157  // argument is the latest solution vector x and the second argument is
158  // the matrix operator.
159  std::function<void(const Vec x, Mat J)> _fnJ;
160 
161  // Function for computing the preconditioner matrix operator. The
162  // first argument is the latest solution vector x and the second
163  // argument is the matrix operator.
164  std::function<void(const Vec x, Mat P)> _fnP;
165 
166  // Function called before the residual and Jacobian function at each
167  // iteration.
168  std::function<void(const Vec x)> _system;
169 
170  // Residual vector
171  Vec _b = nullptr;
172 
173  // Jacobian matrix and preconditioner matrix
174  Mat _matJ = nullptr, _matP = nullptr;
175 
176  // Function to check for convergence
177  std::function<std::pair<double, bool>(const NewtonSolver& solver,
178  const Vec r)>
179  _converged;
180 
181  // Function to update the solution once convergence is reached
182  std::function<void(const NewtonSolver& solver, const Vec dx, Vec x)>
183  _update_solution;
184 
185  // Accumulated number of Krylov iterations since solve began
186  int _krylov_iterations;
187 
188  // Number of iterations
189  int _iteration;
190 
191  // Most recent residual and initial residual
192  double _residual, _residual0;
193 
194  // Linear solver
195  la::petsc::KrylovSolver _solver;
196 
197  // Solution vector
198  Vec _dx = nullptr;
199 
200  // MPI communicator
201  dolfinx::MPI::Comm _comm;
202 };
203 } // namespace nls::petsc
204 } // namespace dolfinx
A duplicate MPI communicator and manage lifetime of the communicator.
Definition: MPI.h:41
This class implements Krylov methods for linear systems of the form Ax = b. It is a wrapper for the K...
Definition: petsc.h:450
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
NewtonSolver(MPI_Comm comm)
Create nonlinear solver.
Definition: NewtonSolver.cpp:67
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 interations. It can can called by functions that check for convergence during a ...
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