DiFfRG
Loading...
Searching...
No Matches
model.hh
Go to the documentation of this file.
1#pragma once
2
3// external libraries
5#include <deal.II/base/point.h>
6#include <deal.II/base/tensor.h>
7
8// DiFfRG
9#include <DiFfRG/model/ad.hh>
12
13namespace DiFfRG
14{
15 using namespace dealii;
16
21 namespace def
22 {
31 template <typename Model, typename Components_> class AbstractModel
32 {
33 Model &asImp() { return static_cast<Model &>(*this); }
34 const Model &asImp() const { return static_cast<const Model &>(*this); }
35
36 protected:
37 Components_ m_components;
38 auto &components() { return m_components; }
39
40 public:
41 const auto &get_components() const { return m_components; }
42 using Components = Components_;
47
57 template <int dim, typename Vector> void initial_condition(const Point<dim> &x, Vector &u_i) const = delete;
58
76 template <int dim, typename NumberType, typename Vector, typename Vector_dot, size_t n_fe_functions>
77 void mass(std::array<NumberType, n_fe_functions> &m_i, const Point<dim> &x, const Vector &u_i,
78 const Vector_dot &dt_u_i) const
79 {
80 // Just to avoid warnings
81 (void)x;
82 (void)u_i;
83
84 for (uint i = 0; i < n_fe_functions; ++i)
85 m_i[i] = dt_u_i[i];
86 }
87
102 template <int dim, typename NumberType, size_t n_fe_functions>
103 void mass(std::array<std::array<NumberType, n_fe_functions>, n_fe_functions> &m_ij, const Point<dim> &x) const
104 {
105 // Just to avoid warnings
106 (void)x;
107
108 for (uint i = 0; i < n_fe_functions; ++i)
109 for (uint j = 0; j < n_fe_functions; ++j)
110 m_ij[i][j] = 0.;
111 for (uint i = 0; i < n_fe_functions; ++i)
112 m_ij[i][i] = 1.;
113 }
114
135 template <int dim, typename NumberType, typename Solutions, size_t n_fe_functions>
136 void flux(std::array<Tensor<1, dim, NumberType>, n_fe_functions> &F_i, const Point<dim> &x,
137 const Solutions &sol) const
138 {
139 // Just to avoid warnings
140 (void)F_i;
141 (void)x;
142 (void)sol;
143 }
144
165 template <int dim, typename NumberType, typename Solutions, size_t n_fe_functions>
166 void source(std::array<NumberType, n_fe_functions> &s_i, const Point<dim> &x, const Solutions &sol) const
167 {
168 // Just to avoid warnings
169 (void)s_i;
170 (void)x;
171 (void)sol;
172 }
173
183 template <uint dim> std::vector<bool> differential_components() const
184 {
185 std::vector<bool> differential_components(Model::Components::count_fe_functions(), false);
186
187 // First we need two reference solutions u_i and dt_u_i, which we then both fill with 1.s
188 std::array<double, Model::Components::count_fe_functions()> u_i{{}};
189 std::array<double, Model::Components::count_fe_functions()> dt_u_i{{}};
190 for (uint i = 0; i < Model::Components::count_fe_functions(); ++i) {
191 u_i[i] = 1.;
192 dt_u_i[i] = 1.;
193 }
194 // Set the point to be at 1. in all directions
195 Point<dim> x;
196 for (uint i = 0; i < dim; ++i)
197 x[i] = 1.;
198 // Get the mass function m_i
199 std::array<double, Model::Components::count_fe_functions()> m_i{{}};
200 asImp().mass(m_i, x, u_i, dt_u_i);
201
202 // Now we check which components are differential by changing dt_u_i slightly and checking whether the mass
203 // function changes.
204 for (uint i = 0; i < Model::Components::count_fe_functions(); ++i) {
205 dt_u_i[i] = 1. + 1e-1;
206 std::array<double, Model::Components::count_fe_functions()> m_i_new{{}};
207 asImp().mass(m_i_new, x, u_i, dt_u_i);
208 dt_u_i[i] = 1.;
209 for (uint j = 0; j < Model::Components::count_fe_functions(); ++j)
210 if (!is_close(m_i[j], m_i_new[j])) differential_components[j] = true;
211 }
212
214 }
215
217
221
222 template <typename Vector> void initial_condition_variables(Vector &v_a) const
223 {
224 // Just to avoid warnings
225 (void)v_a;
226 }
227
228 template <typename Vector, typename Solution> void dt_variables(Vector &r_a, const Solution &sol) const
229 {
230 // Just to avoid warnings
231 (void)r_a;
232 (void)sol;
233 }
234
236
240
241 template <int dim, typename Vector, typename Solutions>
242 void extract(Vector &, const Point<dim> &, const Solutions &) const
243 {
244 }
245
247
251
274 template <uint dependent, int dim, typename NumberType, typename Vector, size_t n_fe_functions_dep>
275 void ldg_flux(std::array<Tensor<1, dim, NumberType>, n_fe_functions_dep> &F, const Point<dim> &x,
276 const Vector &u) const
277 {
278 (void)F;
279 (void)x;
280 (void)u;
281 }
282
305 template <uint dependent, int dim, typename NumberType, typename Vector, size_t n_fe_functions_dep>
306 void ldg_source(std::array<NumberType, n_fe_functions_dep> &s, const Point<dim> &x, const Vector &u) const
307 {
308 (void)s;
309 (void)x;
310 (void)u;
311 }
312
314
318
319 template <int dim, typename NumberType, typename Solutions_s, typename Solutions_n>
320 void face_indicator(std::array<NumberType, 2> & /*indicator*/, const Tensor<1, dim> & /*normal*/,
321 const Point<dim> & /*p*/, const Solutions_s & /*sol_s*/, const Solutions_n & /*sol_n*/) const
322 {
323 }
324
325 template <int dim, typename NumberType, typename Solution>
326 void cell_indicator(NumberType & /*indicator*/, const Point<dim> & /*p*/, const Solution & /*sol*/) const
327 {
328 }
330
334
335 template <int dim, typename Vector> std::array<double, dim> EoM(const Point<dim> &x, const Vector &u) const
336 {
337 // Just to avoid warnings
338 (void)x;
339 return std::array<double, dim>{{u[0]}};
340 }
341
342 template <int dim, typename Vector> Point<dim> EoM_postprocess(const Point<dim> &EoM, const Vector &) const
343 {
344 return EoM;
345 }
346
347 template <typename FUN, typename DataOut> void readouts_multiple(FUN &helper, DataOut &) const
348 {
349 helper([&](const auto &x, const auto &u_i) { return asImp().EoM(x, u_i); }, // chiral EoM
350 [&](auto &output, const auto &x, const auto &sol) { asImp().readouts(output, x, sol); });
351 }
352
353 template <int dim, typename DataOut, typename Solutions>
354 void readouts(DataOut &output, const Point<dim> &x, const Solutions &sol) const
355 {
356 // Just to avoid warnings
357 (void)output;
358 (void)x;
359 (void)sol;
360 }
361
362 template <int dim, typename Constraints>
363 void affine_constraints(Constraints &constraints, const std::vector<IndexSet> &component_boundary_dofs,
364 const std::vector<std::vector<Point<dim>>> &component_boundary_points)
365 {
366 // Just to avoid warnings
367 (void)constraints;
368 (void)component_boundary_dofs;
369 (void)component_boundary_points;
370 }
371
373 };
374
375 class Time
376 {
377 public:
378 void set_time(double t);
379 const double &get_time() const;
380
381 protected:
382 double t;
383 };
384
388 class fRG
389 {
390 public:
396 fRG(double Lambda);
397
403 fRG(const JSONValue &json);
404
410 void set_time(double t);
411
417 const double &get_time() const;
418
419 protected:
420 const double Lambda;
421 double t, k, k2, k3, k4, k5, k6;
422 };
423 } // namespace def
424} // namespace DiFfRG
A wrapper around the boost json value class.
Definition json.hh:19
The abstract interface for any numerical model. Most methods have a standard implementation,...
Definition model.hh:32
std::vector< bool > differential_components() const
A method to find out which components of the mass function are differential when using a DAE.
Definition model.hh:183
void dt_variables(Vector &r_a, const Solution &sol) const
Definition model.hh:228
void source(std::array< NumberType, n_fe_functions > &s_i, const Point< dim > &x, const Solutions &sol) const
The source function is implemented by this method.
Definition model.hh:166
auto & components()
Definition model.hh:38
const Model & asImp() const
Definition model.hh:34
Components_ m_components
Definition model.hh:37
void cell_indicator(NumberType &, const Point< dim > &, const Solution &) const
Definition model.hh:326
void mass(std::array< std::array< NumberType, n_fe_functions >, n_fe_functions > &m_ij, const Point< dim > &x) const
If not using a DAE, the mass matrix is implemented in this method.
Definition model.hh:103
Model & asImp()
Definition model.hh:33
const auto & get_components() const
Definition model.hh:41
void flux(std::array< Tensor< 1, dim, NumberType >, n_fe_functions > &F_i, const Point< dim > &x, const Solutions &sol) const
The flux function is implemented by this method.
Definition model.hh:136
void ldg_flux(std::array< Tensor< 1, dim, NumberType >, n_fe_functions_dep > &F, const Point< dim > &x, const Vector &u) const
The LDG flux function is implemented by this method.
Definition model.hh:275
void face_indicator(std::array< NumberType, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Solutions_s &, const Solutions_n &) const
Definition model.hh:320
void mass(std::array< NumberType, n_fe_functions > &m_i, const Point< dim > &x, const Vector &u_i, const Vector_dot &dt_u_i) const
The mass function is implemented in this method.
Definition model.hh:77
void extract(Vector &, const Point< dim > &, const Solutions &) const
Definition model.hh:242
void readouts_multiple(FUN &helper, DataOut &) const
Definition model.hh:347
void initial_condition(const Point< dim > &x, Vector &u_i) const =delete
This method implements the initial condition for the FE functions.
void affine_constraints(Constraints &constraints, const std::vector< IndexSet > &component_boundary_dofs, const std::vector< std::vector< Point< dim > > > &component_boundary_points)
Definition model.hh:363
void ldg_source(std::array< NumberType, n_fe_functions_dep > &s, const Point< dim > &x, const Vector &u) const
The LDG source function is implemented by this method.
Definition model.hh:306
std::array< double, dim > EoM(const Point< dim > &x, const Vector &u) const
Definition model.hh:335
void initial_condition_variables(Vector &v_a) const
Definition model.hh:222
void readouts(DataOut &output, const Point< dim > &x, const Solutions &sol) const
Definition model.hh:354
Point< dim > EoM_postprocess(const Point< dim > &EoM, const Vector &) const
Definition model.hh:342
Components_ Components
Definition model.hh:42
Definition model.hh:376
const double & get_time() const
void set_time(double t)
double t
Definition model.hh:382
The fRG class is used to keep track of the RG time and the cutoff scale.
Definition model.hh:389
double k3
Definition model.hh:421
double k
Definition model.hh:421
double k4
Definition model.hh:421
double k6
Definition model.hh:421
fRG(const JSONValue &json)
Construct a new fRG object from a given JSONValue object.
double t
Definition model.hh:421
fRG(double Lambda)
Construct a new fRG object from a given initial cutoff scale.
double k5
Definition model.hh:421
const double Lambda
Definition model.hh:420
void set_time(double t)
Set the time of the fRG object, updating the cutoff scale and its powers.
const double & get_time() const
Get the time of the fRG object.
double k2
Definition model.hh:421
Definition complex_math.hh:14
bool __forceinline__ __host__ __device__ is_close(T1 a, T2 b, T3 eps_)
Function to evaluate whether two floats are equal to numerical precision. Tests for both relative and...
Definition math.hh:160
unsigned int uint
Definition utils.hh:22