DiFfRG
Loading...
Searching...
No Matches
kinsol.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/sundials/kinsol.h>
8
9namespace DiFfRG
10{
11 using namespace dealii;
17 template <typename VectorType_> class KINSOL
18 {
19 public:
20 using VectorType = VectorType_;
21
22 KINSOL(double abstol_ = 1e-10, double reltol_ = 1e-6, double assemble_threshold_ = 0., unsigned int max_steps_ = 21,
23 unsigned int n_stepsize_iterations_ = 21);
24
28 void reinit(const VectorType &u);
29
35 void operator()(VectorType &iterate);
36
42 double threshold(const double thr);
43
44 const VectorType &get_residual() const;
45 const VectorType &get_initial() const;
46 const VectorType &get_iterate() const;
47
48 double average_time_newton() const;
49 unsigned int num_newton_calls() const;
50
57 void set_ignore_nonconv(bool x);
58
59 std::function<void(VectorType &, const VectorType &)> residual;
60 std::function<void(VectorType &, const VectorType &)> lin_solve;
61 std::function<void(const VectorType &)> update_jacobian;
62
66 double get_error();
67
71 unsigned int get_step();
72
76 unsigned int get_jacobians();
77
78 private:
84
87 double resnorm;
88 unsigned int step, jacobians;
89
90 std::vector<double> timings_newton;
91
97 double get_EEst();
98
105 bool check();
106
107 std::shared_ptr<SUNDIALS::KINSOL<VectorType>> kinsol;
108 };
109} // namespace DiFfRG
A newton solver, using local error estimates for each vector component.
Definition kinsol.hh:18
void set_ignore_nonconv(bool x)
Make the newton algorithm return a success even if the evolution did not fulfill the convergence crit...
unsigned int jacobians
Definition kinsol.hh:88
double reltol
Definition kinsol.hh:82
double resnorm
Definition kinsol.hh:87
std::function< void(VectorType &, const VectorType &)> lin_solve
Definition kinsol.hh:60
unsigned int num_newton_calls() const
const VectorType & get_initial() const
const VectorType & get_residual() const
double get_EEst()
Calculate the current error estimate.
void operator()(VectorType &iterate)
Perform a newton iteration on the input vector.
void reinit(const VectorType &u)
Should it be necessary, one can force a recalculation of the jacobian hereby.
double abstol
Definition kinsol.hh:82
std::shared_ptr< SUNDIALS::KINSOL< VectorType > > kinsol
Definition kinsol.hh:107
unsigned int n_stepsize_iterations
Definition kinsol.hh:83
unsigned int max_steps
Definition kinsol.hh:83
VectorType_ VectorType
Definition kinsol.hh:20
bool ignore_nonconv
Definition kinsol.hh:86
std::function< void(VectorType &, const VectorType &)> residual
Definition kinsol.hh:59
double average_time_newton() const
VectorType residual_vector
Definition kinsol.hh:79
VectorType initial_vector
Definition kinsol.hh:80
double assemble_threshold
Definition kinsol.hh:82
std::function< void(const VectorType &)> update_jacobian
Definition kinsol.hh:61
unsigned int get_jacobians()
Get the number of jacobians created in the latest iteration.
unsigned int get_step()
Get the number of steps taken in the latest iteration.
unsigned int step
Definition kinsol.hh:88
double threshold(const double thr)
Set the jacobian recalculation threshold, i.e. below which reduction we force a jacobian update.
std::vector< double > timings_newton
Definition kinsol.hh:90
double get_error()
Get the latest error.
KINSOL(double abstol_=1e-10, double reltol_=1e-6, double assemble_threshold_=0., unsigned int max_steps_=21, unsigned int n_stepsize_iterations_=21)
VectorType iterate_vector
Definition kinsol.hh:81
const VectorType & get_iterate() const
bool converged
Definition kinsol.hh:85
bool check()
Check if the iteration should go on.
Definition complex_math.hh:14