DiFfRG
Loading...
Searching...
No Matches
newton.hh
Go to the documentation of this file.
1#pragma once
2
3// standard library
4#include <functional>
5
6// external libraries
7#include <deal.II/base/config.h>
8#include <deal.II/base/timer.h>
9#include <deal.II/lac/vector_memory.h>
10
11namespace DiFfRG
12{
13 using namespace dealii;
14
20 template <typename VectorType_> class Newton
21 {
22 public:
23 using VectorType = VectorType_;
24
25 Newton(double abstol_ = 1e-10, double reltol_ = 1e-6, double assemble_threshold_ = 0., uint max_steps_ = 21,
26 uint n_stepsize_iterations_ = 21)
27 : abstol(abstol_), reltol(reltol_), assemble_threshold(assemble_threshold_), max_steps(max_steps_),
28 n_stepsize_iterations(n_stepsize_iterations_), ignore_nonconv(false)
29 {
30 }
31
35 void reinit(const VectorType &u) { update_jacobian(u); }
36
42 void operator()(VectorType &iterate)
43 {
44 Timer timer;
45
46 converged = false;
47 step = 0;
48 resnorm = 1.;
49 jacobians = 0;
50
51 initial_vector = iterate;
52 iterate_vector = iterate;
54
55 GrowingVectorMemory<VectorType> mem;
56 typename VectorMemory<VectorType>::Pointer Du(mem);
57
58 Du->reinit(u);
59 residual_vector.reinit(u);
60
61 // fill res with (f(u), v)
63 if (!std::isfinite(residual_vector.l2_norm())) throw std::runtime_error("Initial residual non-finite!");
64 resnorm = get_EEst();
65 double old_residual = 0.;
66 double assembly_time = 0.;
67
68 while (check()) {
69 step++;
70
71 // assemble (Df(u), v)
72 if ((step > 1) && (resnorm / old_residual >= assemble_threshold)) {
73 Timer timer_assembly;
75 assembly_time += timer_assembly.wall_time();
76 jacobians++;
77 }
78
79 Du->reinit(u);
80 try {
81 lin_solve(*Du.get(), residual_vector);
82 } catch (std::exception &) {
83 throw std::runtime_error("Linear solve failed!");
84 }
85
86 u.add(-1., *Du);
87 old_residual = resnorm;
88
90 resnorm = get_EEst();
91 if (!std::isfinite(residual_vector.l2_norm())) {
92 std::cerr << "Residual non-finite!" << std::endl;
93 throw std::runtime_error("Residual non-finite!");
94 }
95
96 // Simple line search: decrease the step size until the residual is smaller than the one in the step before
97 uint step_size = 0;
98 while (resnorm >= old_residual) {
99 ++step_size;
100 if (step_size > n_stepsize_iterations) {
101 break;
102 }
103 u.add(1. / double(1 << step_size), *Du);
105 resnorm = get_EEst();
106 }
107 }
108
109 timings_newton.push_back(timer.wall_time() - assembly_time);
110
111 // in case of failure: throw exception
112 if (!(converged || ignore_nonconv) || !std::isfinite(resnorm)) {
113 throw std::runtime_error("Newton did not converge!");
114 }
115 // otherwise exit as normal
116 else
117 iterate = iterate_vector;
118 }
119
125 double threshold(const double thr)
126 {
127 const double t = assemble_threshold;
128 assemble_threshold = thr;
129 return t;
130 }
131
132 const VectorType &get_residual() const { return residual_vector; }
133 const VectorType &get_initial() const { return initial_vector; }
134 const VectorType &get_iterate() const { return iterate_vector; }
135
136 double average_time_newton() const
137 {
138 double t = 0.;
139 double n = timings_newton.size();
140 for (const auto &t_ : timings_newton)
141 t += t_ / n;
142 return t;
143 }
144 uint num_newton_calls() const { return timings_newton.size(); }
145
152 void set_ignore_nonconv(bool x) { ignore_nonconv = x; }
153
154 std::function<void(VectorType &, const VectorType &)> residual;
155 std::function<void(VectorType &, const VectorType &)> lin_solve;
156 std::function<void(const VectorType &)> update_jacobian;
157
161 double get_error() { return resnorm; }
162
166 uint get_step() { return step; }
167
172
173 private:
179
182 double resnorm;
184
185 std::vector<double> timings_newton;
186
192 double get_EEst()
193 {
194 const VectorType &err = residual_vector;
195 const VectorType &u_prev = initial_vector;
196 const VectorType &u = iterate_vector;
197
198 GrowingVectorMemory<VectorType> mem;
199 typename VectorMemory<VectorType>::Pointer tmp(mem);
200 tmp->reinit(u, true);
201
202 using std::abs, std::max, std::sqrt;
203 tbb::parallel_for(tbb::blocked_range<uint>(0, tmp->size()), [&](tbb::blocked_range<uint> r) {
204 for (uint n = r.begin(); n < r.end(); ++n) {
205 (*tmp)[n] = abs(err[n]) / (abstol + max(abs(u_prev[n]), abs(u[n])) * reltol);
206 }
207 });
208 return tmp->l2_norm() / sqrt(tmp->size());
209 return tmp->linfty_norm(); // tmp->l2_norm() / sqrt(tmp->size());
210 }
211
218 bool check()
219 {
220 if (step > max_steps) return false;
221 if (step == 0) return true;
222 double err = get_EEst();
223 converged = err <= 1.;
224 return !converged;
225 }
226 };
227} // namespace DiFfRG
A newton solver, using local error estimates for each vector component.
Definition newton.hh:21
Newton(double abstol_=1e-10, double reltol_=1e-6, double assemble_threshold_=0., uint max_steps_=21, uint n_stepsize_iterations_=21)
Definition newton.hh:25
uint num_newton_calls() const
Definition newton.hh:144
bool ignore_nonconv
Definition newton.hh:181
double average_time_newton() const
Definition newton.hh:136
double get_EEst()
Calculate the current error estimate.
Definition newton.hh:192
uint n_stepsize_iterations
Definition newton.hh:178
bool check()
Check if the iteration should go on.
Definition newton.hh:218
std::function< void(const VectorType &)> update_jacobian
Definition newton.hh:156
const VectorType & get_iterate() const
Definition newton.hh:134
void set_ignore_nonconv(bool x)
Make the newton algorithm return a success even if the evolution did not fulfill the convergence crit...
Definition newton.hh:152
std::function< void(VectorType &, const VectorType &)> residual
Definition newton.hh:154
double get_error()
Get the latest error.
Definition newton.hh:161
const VectorType & get_initial() const
Definition newton.hh:133
uint max_steps
Definition newton.hh:178
double assemble_threshold
Definition newton.hh:177
double resnorm
Definition newton.hh:182
VectorType iterate_vector
Definition newton.hh:176
uint get_step()
Get the number of steps taken in the latest iteration.
Definition newton.hh:166
const VectorType & get_residual() const
Definition newton.hh:132
uint get_jacobians()
Get the number of jacobians created in the latest iteration.
Definition newton.hh:171
void reinit(const VectorType &u)
Should it be necessary, one can force a recalculation of the jacobian hereby.
Definition newton.hh:35
VectorType initial_vector
Definition newton.hh:175
double threshold(const double thr)
Set the jacobian recalculation threshold, i.e. below which reduction we force a jacobian update.
Definition newton.hh:125
VectorType residual_vector
Definition newton.hh:174
bool converged
Definition newton.hh:180
std::function< void(VectorType &, const VectorType &)> lin_solve
Definition newton.hh:155
VectorType_ VectorType
Definition newton.hh:23
void operator()(VectorType &iterate)
Perform a newton iteration on the input vector.
Definition newton.hh:42
double abstol
Definition newton.hh:177
std::vector< double > timings_newton
Definition newton.hh:185
uint step
Definition newton.hh:183
uint jacobians
Definition newton.hh:183
double reltol
Definition newton.hh:177
Definition complex_math.hh:14
unsigned int uint
Definition utils.hh:22