DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
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#ifdef HAS_PETSC
10
11#include <dolfinx/common/MPI.h>
12#include <dolfinx/la/petsc.h>
13#include <functional>
14#include <memory>
15#include <petscmat.h>
16#include <petscvec.h>
17#include <utility>
18
19namespace dolfinx
20{
21
22namespace la::petsc
23{
24class KrylovSolver;
25} // namespace la::petsc
26
27namespace nls::petsc
28{
29
32
34{
35public:
38 explicit NewtonSolver(MPI_Comm comm);
39
41 NewtonSolver(NewtonSolver&& solver) = default;
42
43 // Copy constructor (deleted)
44 NewtonSolver(const NewtonSolver& solver) = delete;
45
47 NewtonSolver& operator=(NewtonSolver&& solver) = default;
48
49 // Assignment operator (deleted)
50 NewtonSolver& operator=(const NewtonSolver& solver) = delete;
51
54
59 void setF(std::function<void(const Vec, Vec)> F, Vec b);
60
65 void setJ(std::function<void(const Vec, Mat)> J, Mat Jmat);
66
71 void setP(std::function<void(const Vec, Mat)> P, Mat Pmat);
72
80
88
93 void set_form(std::function<void(Vec)> form);
94
99 std::function<std::pair<double, bool>(const NewtonSolver&, const Vec)> c);
100
104 void set_update(
105 std::function<void(const NewtonSolver& solver, const Vec, Vec)> update);
106
112 std::pair<int, bool> solve(Vec x);
113
117 int iteration() const;
118
122 int krylov_iterations() const;
123
126 double residual() const;
127
130 double residual0() const;
131
133 MPI_Comm comm() const;
134
136 int max_it = 50;
137
139 double rtol = 1e-9;
140
142 double atol = 1e-10;
143
144 // FIXME: change to string to enum
146 std::string convergence_criterion = "residual";
147
149 bool report = true;
150
153
156
157private:
158 // Function for computing the residual vector. The first argument is
159 // the latest solution vector x and the second argument is the
160 // residual vector.
161 std::function<void(const Vec x, Vec b)> _fnF;
162
163 // Function for computing the Jacobian matrix operator. The first
164 // argument is the latest solution vector x and the second argument is
165 // the matrix operator.
166 std::function<void(const Vec x, Mat J)> _fnJ;
167
168 // Function for computing the preconditioner matrix operator. The
169 // first argument is the latest solution vector x and the second
170 // argument is the matrix operator.
171 std::function<void(const Vec x, Mat P)> _fnP;
172
173 // Function called before the residual and Jacobian function at each
174 // iteration.
175 std::function<void(const Vec x)> _system;
176
177 // Residual vector
178 Vec _b = nullptr;
179
180 // Jacobian matrix and preconditioner matrix
181 Mat _matJ = nullptr, _matP = nullptr;
182
183 // Function to check for convergence
184 std::function<std::pair<double, bool>(const NewtonSolver& solver,
185 const Vec r)>
186 _converged;
187
188 // Function to update the solution once convergence is reached
189 std::function<void(const NewtonSolver& solver, const Vec dx, Vec x)>
190 _update_solution;
191
192 // Accumulated number of Krylov iterations since solve began
193 int _krylov_iterations;
194
195 // Number of iterations
196 int _iteration;
197
198 // Most recent residual and initial residual
199 double _residual, _residual0;
200
201 // Linear solver
203
204 // Solution vector
205 Vec _dx = nullptr;
206
207 // MPI communicator
208 dolfinx::MPI::Comm _comm;
209};
210} // namespace nls::petsc
211} // namespace dolfinx
212
213#endif
A duplicate MPI communicator and manage lifetime of the communicator.
Definition MPI.h:43
Definition petsc.h:453
Definition NewtonSolver.h:34
double relaxation_parameter
Relaxation parameter.
Definition NewtonSolver.h:155
int krylov_iterations() const
Get number of Krylov iterations elapsed since solve started.
Definition NewtonSolver.cpp:276
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 ...
Definition NewtonSolver.cpp:128
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.
Definition NewtonSolver.cpp:93
NewtonSolver(NewtonSolver &&solver)=default
Move constructor.
double rtol
Relative tolerance.
Definition NewtonSolver.h:139
~NewtonSolver()
Destructor.
Definition NewtonSolver.cpp:81
NewtonSolver & operator=(NewtonSolver &&solver)=default
Move assignment constructor.
double atol
Absolute tolerance.
Definition NewtonSolver.h:142
NewtonSolver(MPI_Comm comm)
Definition NewtonSolver.cpp:69
double residual() const
Get current residual.
Definition NewtonSolver.cpp:283
void setP(std::function< void(const Vec, Mat)> P, Mat Pmat)
Set the function for computing the preconditioner matrix (optional).
Definition NewtonSolver.cpp:109
double residual0() const
Return initial residual.
Definition NewtonSolver.cpp:285
std::string convergence_criterion
Convergence criterion.
Definition NewtonSolver.h:146
int iteration() const
The number of Newton iterations. It can can called by functions that check for convergence during a s...
Definition NewtonSolver.cpp:281
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.
Definition NewtonSolver.cpp:139
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.
Definition NewtonSolver.cpp:101
bool report
Monitor convergence.
Definition NewtonSolver.h:149
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.
Definition NewtonSolver.cpp:133
bool error_on_nonconvergence
Throw error if solver fails to converge.
Definition NewtonSolver.h:152
MPI_Comm comm() const
Get MPI communicator.
Definition NewtonSolver.cpp:287
const la::petsc::KrylovSolver & get_krylov_solver() const
Get the internal Krylov solver used to solve for the Newton updates (const version).
Definition NewtonSolver.cpp:118
std::pair< int, bool > solve(Vec x)
Solve the nonlinear problem for given and Jacobian .
Definition NewtonSolver.cpp:145
int max_it
Maximum number of iterations.
Definition NewtonSolver.h:136
Top-level namespace.
Definition defines.h:12