5#include <deal.II/base/point.h>
6#include <deal.II/base/tensor.h>
15 using namespace dealii;
33 Model &
asImp() {
return static_cast<Model &
>(*this); }
34 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
57 template <
int dim,
typename Vector>
void initial_condition(
const Point<dim> &x, Vector &u_i)
const =
delete;
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
84 for (
uint i = 0; i < n_fe_functions; ++i)
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
108 for (
uint i = 0; i < n_fe_functions; ++i)
109 for (
uint j = 0; j < n_fe_functions; ++j)
111 for (
uint i = 0; i < n_fe_functions; ++i)
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
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
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) {
196 for (
uint i = 0; i < dim; ++i)
199 std::array<double, Model::Components::count_fe_functions()> m_i{{}};
200 asImp().mass(m_i, x, u_i, dt_u_i);
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);
209 for (
uint j = 0; j < Model::Components::count_fe_functions(); ++j)
228 template <
typename Vector,
typename Solution>
void dt_variables(Vector &r_a,
const Solution &sol)
const
241 template <
int dim,
typename Vector,
typename Solutions>
242 void extract(Vector &,
const Point<dim> &,
const Solutions &)
const
274 template <u
int 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
305 template <u
int 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
319 template <
int dim,
typename NumberType,
typename Solutions_s,
typename Solutions_n>
321 const Point<dim> & ,
const Solutions_s & ,
const Solutions_n & )
const
325 template <
int dim,
typename NumberType,
typename Solution>
335 template <
int dim,
typename Vector> std::array<double, dim>
EoM(
const Point<dim> &x,
const Vector &u)
const
339 return std::array<double, dim>{{u[0]}};
342 template <
int dim,
typename Vector> Point<dim>
EoM_postprocess(
const Point<dim> &
EoM,
const Vector &)
const
349 helper([&](
const auto &x,
const auto &u_i) {
return asImp().EoM(x, u_i); },
350 [&](
auto &output,
const auto &x,
const auto &sol) {
asImp().readouts(output, x, sol); });
353 template <
int dim,
typename DataOut,
typename Solutions>
354 void readouts(DataOut &output,
const Point<dim> &x,
const Solutions &sol)
const
362 template <
int dim,
typename Constra
ints>
363 void affine_constraints(Constraints &constraints,
const std::vector<IndexSet> &component_boundary_dofs,
364 const std::vector<std::vector<Point<dim>>> &component_boundary_points)
368 (void)component_boundary_dofs;
369 (void)component_boundary_points;
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
const double & get_time() const
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