AbstractModel< Model, Components_ > Class Template Reference#
|
DiFfRG
|
The abstract interface for any numerical model. Most methods have a standard implementation, which can be overwritten if needed. To see how the models are used, refer to the DiFfRG::AbstractAssembler class and the guide. More...
#include <model.hh>
Public Types | |
| using | Components = Components_ |
Public Member Functions | |
| const auto & | get_components () const |
Spatial discretization | |
| template<int dim, typename Vector > | |
| void | initial_condition (const Point< dim > &x, Vector &u_i) const =delete |
| This method implements the initial condition for the FE functions. | |
| template<int dim, typename NumberType , typename Vector , typename Vector_dot , size_t n_fe_functions> | |
| 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 \(m_i(\partial_t u_j, u_j, x)\) is implemented in this method. | |
| template<int dim, typename NumberType , size_t n_fe_functions> | |
| 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 \(m_{ij}(x)\) is implemented in this method. | |
| template<int dim, typename NumberType , typename Solutions , size_t n_fe_functions> | |
| void | flux (std::array< Tensor< 1, dim, NumberType >, n_fe_functions > &F_i, const Point< dim > &x, const Solutions &sol) const |
| The flux function \(F_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x)\) is implemented by this method. | |
| template<int dim, typename NumberType , typename Solutions , size_t n_fe_functions> | |
| void | source (std::array< NumberType, n_fe_functions > &s_i, const Point< dim > &x, const Solutions &sol) const |
| The source function \(s_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x)\) is implemented by this method. | |
| template<uint dim> | |
| std::vector< bool > | differential_components () const |
| A method to find out which components of the mass function are differential when using a DAE. | |
Other variables | |
| template<typename Vector > | |
| void | initial_condition_variables (Vector &v_a) const |
| template<typename Vector , typename Solution > | |
| void | dt_variables (Vector &r_a, const Solution &sol) const |
Extractors | |
| template<int dim, typename Vector , typename Solutions > | |
| void | extract (Vector &, const Point< dim > &, const Solutions &) const |
LDG equations | |
| template<uint dependent, int dim, typename NumberType , typename Vector , size_t n_fe_functions_dep> | |
| 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 \(F^{LDG}_i(u_j, x),\,i>0\) is implemented by this method. | |
| template<uint dependent, int dim, typename NumberType , typename Vector , size_t n_fe_functions_dep> | |
| void | ldg_source (std::array< NumberType, n_fe_functions_dep > &s, const Point< dim > &x, const Vector &u) const |
| The LDG source function \(s^{LDG}_i(u_j, x),\,i>0\) is implemented by this method. | |
Adaptive mesh refinement | |
| template<int dim, typename NumberType , typename Solutions_s , typename Solutions_n > | |
| void | face_indicator (std::array< NumberType, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Solutions_s &, const Solutions_n &) const |
| template<int dim, typename NumberType , typename Solution > | |
| void | cell_indicator (NumberType &, const Point< dim > &, const Solution &) const |
Other | |
| template<int dim, typename Vector > | |
| std::array< double, dim > | EoM (const Point< dim > &x, const Vector &u) const |
| template<int dim, typename Vector > | |
| Point< dim > | EoM_postprocess (const Point< dim > &EoM, const Vector &) const |
| template<typename FUN , typename DataOut > | |
| void | readouts_multiple (FUN &helper, DataOut &) const |
| template<int dim, typename DataOut , typename Solutions > | |
| void | readouts (DataOut &output, const Point< dim > &x, const Solutions &sol) const |
| template<int dim, typename Constraints > | |
| void | affine_constraints (Constraints &constraints, const std::vector< IndexSet > &component_boundary_dofs, const std::vector< std::vector< Point< dim > > > &component_boundary_points) |
Protected Member Functions | |
| auto & | components () |
Protected Attributes | |
| Components_ | m_components |
Private Member Functions | |
| Model & | asImp () |
| const Model & | asImp () const |
Detailed Description
class DiFfRG::def::AbstractModel< Model, Components_ >
The abstract interface for any numerical model. Most methods have a standard implementation, which can be overwritten if needed. To see how the models are used, refer to the DiFfRG::AbstractAssembler class and the guide.
- Template Parameters
-
Model The model which implements this interface. (CRTP) Components_ The components of the model, this must be a DiFfRG::ComponentDescriptor.
Member Typedef Documentation
◆ Components
| using DiFfRG::def::AbstractModel< Model, Components_ >::Components = Components_ |
Member Function Documentation
◆ affine_constraints()
|
inline |
◆ asImp() [1/2]
|
inlineprivate |
◆ asImp() [2/2]
|
inlineprivate |
◆ cell_indicator()
|
inline |
◆ components()
|
inlineprotected |
◆ differential_components()
|
inline |
A method to find out which components of the mass function are differential when using a DAE.
- Note
- The standard implementation of this method tests whether the mass function changes when changing the time derivative of one component slightly. For highly complicated models, this method might not be able to set all differential components correctly.
- Returns
- std::vector<bool> with
truefor differential components andfalsefor algebraic components.
◆ dt_variables()
|
inline |
◆ EoM()
|
inline |
◆ EoM_postprocess()
|
inline |
◆ extract()
|
inline |
◆ face_indicator()
|
inline |
◆ flux()
|
inline |
The flux function \(F_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x)\) is implemented by this method.
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model. It is however necessary to leave at least the NumberType, Vector, and Vector_dot template parameters, as these can differ between calls (e.g. when doing automatic differentiation).
- Note
- The standard implementation of this method simply sets \(F_i = 0\).
- Parameters
-
F_i the resulting flux function \(F_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates. sol a std::tuple<...>which contains- the array u_j
- the array of arrays \(\partial_x u_j\)
- the array of arrays of arrays \(\partial_x^2 u_j\)
- the array of extractors \(e_b\)
◆ get_components()
|
inline |
◆ initial_condition()
|
delete |
This method implements the initial condition for the FE functions.
- Note
- No standard implementation is given, this method has to be reimplemented whenever one uses FE functions.
- Parameters
-
x a d-dimensional dealii::Point<dim> representing field coordinates. u_i the field values \(u_i(x)\) at the point x. This method should fill this argument with the desired initial condition.
◆ initial_condition_variables()
|
inline |
◆ ldg_flux()
|
inline |
The LDG flux function \(F^{LDG}_i(u_j, x),\,i>0\) is implemented by this method.
The assembler constructs the i-th LDG function l_i from the i-1-th level as
\[l_i = \partial_x F^{LDG}_i(l_{i-1}, x) + s^{LDG}_i(l_{i-1}, x)\]
Here, \(l_0\) is the solution itself (with all its components).
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model.
- Note
- The standard implementation of this method simply sets \(F^{LDG}_i = 0\).
- Template Parameters
-
dependent the index \(i\) of the dependent variable \(l_i\) which is constructed from the previous level \(l_{i-1}\).
- Parameters
-
F the resulting LDG flux function \(F^{LDG}_i\), with n_fe_functions_dep components. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates. u the field values of \(l_j(x)\) at the point x.
◆ ldg_source()
|
inline |
The LDG source function \(s^{LDG}_i(u_j, x),\,i>0\) is implemented by this method.
The assembler constructs the i-th LDG function l_i from the i-1-th level as
\[l_i = \partial_x F^{LDG}_i(l_{i-1}, x) + s^{LDG}_i(l_{i-1}, x)\]
Here, \(l_0\) is the solution itself (with all its components).
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model.
- Note
- The standard implementation of this method simply sets \(s^{LDG}_i = 0\).
- Template Parameters
-
dependent the index \(i\) of the dependent variable \(l_i\) which is constructed from the previous level \(l_{i-1}\).
- Parameters
-
s the resulting LDG source function \(s^{LDG}_i\), with n_fe_functions_dep components. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates. u the field values of \(l_j(x)\) at the point x.
◆ mass() [1/2]
|
inline |
The mass function \(m_i(\partial_t u_j, u_j, x)\) is implemented in this method.
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model. It is however necessary to leave at least the NumberType, Vector, and Vector_dot template parameters, as these can differ between calls (e.g. when doing automatic differentiation).
- Note
- The standard implementation of this method simply sets \(m_i = \partial_t u_i\).
- Parameters
-
m_i the resulting mass function \(m_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates. u_i the field values \(u_i(x)\) at the point x.dt_u_i the time derivative of the field values \(\partial_t u_i(x)\) at the point x.
◆ mass() [2/2]
|
inline |
If not using a DAE, the mass matrix \(m_{ij}(x)\) is implemented in this method.
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model. It is however necessary to leave at least the NumberType, Vector, and Vector_dot template parameters, as these can differ between calls (e.g. when doing automatic differentiation).
- Note
- The standard implementation of this method simply sets \(m_{ij} = \delta_{ij}\).
- Parameters
-
m_ij the resulting mass matrix \(m_{ij}\), with \(N_f\) components in each dimension. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates.
◆ readouts()
|
inline |
◆ readouts_multiple()
|
inline |
◆ source()
|
inline |
The source function \(s_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x)\) is implemented by this method.
- Remarks
- Note, that the precise template structure is not important, the only important thing is that the types are consistent with the rest of the model. It is however necessary to leave at least the NumberType, Vector, and Vector_dot template parameters, as these can differ between calls (e.g. when doing automatic differentiation).
- Note
- The standard implementation of this method simply sets \(s_i = 0\).
- Parameters
-
s_i the resulting source function \(s_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation. x a d-dimensional dealii::Point<dim> representing field coordinates. sol a std::tuple<...>which contains- the array u_j
- the array of arrays \(\partial_x u_j\)
- the array of arrays of arrays \(\partial_x^2 u_j\)
- the array of extractors \(e_b\)
Member Data Documentation
◆ m_components
|
protected |
The documentation for this class was generated from the following file:
- /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/model/model.hh
Generated by