DiFfRG
Loading...
Searching...
No Matches
DiFfRG::def::AbstractModel< Model, Components_ > Class Template Reference

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

template<typename Model, typename Components_>
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
ModelThe model which implements this interface. (CRTP)
Components_The components of the model, this must be a DiFfRG::ComponentDescriptor.

Member Typedef Documentation

◆ Components

template<typename Model , typename Components_ >
using DiFfRG::def::AbstractModel< Model, Components_ >::Components = Components_

Member Function Documentation

◆ affine_constraints()

template<typename Model , typename Components_ >
template<int dim, typename Constraints >
void DiFfRG::def::AbstractModel< Model, Components_ >::affine_constraints ( Constraints & constraints,
const std::vector< IndexSet > & component_boundary_dofs,
const std::vector< std::vector< Point< dim > > > & component_boundary_points )
inline

◆ asImp() [1/2]

template<typename Model , typename Components_ >
Model & DiFfRG::def::AbstractModel< Model, Components_ >::asImp ( )
inlineprivate

◆ asImp() [2/2]

template<typename Model , typename Components_ >
const Model & DiFfRG::def::AbstractModel< Model, Components_ >::asImp ( ) const
inlineprivate

◆ cell_indicator()

template<typename Model , typename Components_ >
template<int dim, typename NumberType , typename Solution >
void DiFfRG::def::AbstractModel< Model, Components_ >::cell_indicator ( NumberType & ,
const Point< dim > & ,
const Solution &  ) const
inline

◆ components()

template<typename Model , typename Components_ >
auto & DiFfRG::def::AbstractModel< Model, Components_ >::components ( )
inlineprotected

◆ differential_components()

template<typename Model , typename Components_ >
template<uint dim>
std::vector< bool > DiFfRG::def::AbstractModel< Model, Components_ >::differential_components ( ) const
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 true for differential components and false for algebraic components.

◆ dt_variables()

template<typename Model , typename Components_ >
template<typename Vector , typename Solution >
void DiFfRG::def::AbstractModel< Model, Components_ >::dt_variables ( Vector & r_a,
const Solution & sol ) const
inline

◆ EoM()

template<typename Model , typename Components_ >
template<int dim, typename Vector >
std::array< double, dim > DiFfRG::def::AbstractModel< Model, Components_ >::EoM ( const Point< dim > & x,
const Vector & u ) const
inline

◆ EoM_postprocess()

template<typename Model , typename Components_ >
template<int dim, typename Vector >
Point< dim > DiFfRG::def::AbstractModel< Model, Components_ >::EoM_postprocess ( const Point< dim > & EoM,
const Vector &  ) const
inline

◆ extract()

template<typename Model , typename Components_ >
template<int dim, typename Vector , typename Solutions >
void DiFfRG::def::AbstractModel< Model, Components_ >::extract ( Vector & ,
const Point< dim > & ,
const Solutions &  ) const
inline

◆ face_indicator()

template<typename Model , typename Components_ >
template<int dim, typename NumberType , typename Solutions_s , typename Solutions_n >
void DiFfRG::def::AbstractModel< Model, Components_ >::face_indicator ( std::array< NumberType, 2 > & ,
const Tensor< 1, dim > & ,
const Point< dim > & ,
const Solutions_s & ,
const Solutions_n &  ) const
inline

◆ flux()

template<typename Model , typename Components_ >
template<int dim, typename NumberType , typename Solutions , size_t n_fe_functions>
void DiFfRG::def::AbstractModel< Model, Components_ >::flux ( std::array< Tensor< 1, dim, NumberType >, n_fe_functions > & F_i,
const Point< dim > & x,
const Solutions & sol ) const
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_ithe resulting flux function \(F_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation.
xa d-dimensional dealii::Point<dim> representing field coordinates.
sola std::tuple<...> which contains
  1. the array u_j
  2. the array of arrays \(\partial_x u_j\)
  3. the array of arrays of arrays \(\partial_x^2 u_j\)
  4. the array of extractors \(e_b\)

◆ get_components()

template<typename Model , typename Components_ >
const auto & DiFfRG::def::AbstractModel< Model, Components_ >::get_components ( ) const
inline

◆ initial_condition()

template<typename Model , typename Components_ >
template<int dim, typename Vector >
void DiFfRG::def::AbstractModel< Model, Components_ >::initial_condition ( const Point< dim > & x,
Vector & u_i ) const
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
xa d-dimensional dealii::Point<dim> representing field coordinates.
u_ithe field values \(u_i(x)\) at the point x. This method should fill this argument with the desired initial condition.

◆ initial_condition_variables()

template<typename Model , typename Components_ >
template<typename Vector >
void DiFfRG::def::AbstractModel< Model, Components_ >::initial_condition_variables ( Vector & v_a) const
inline

◆ ldg_flux()

template<typename Model , typename Components_ >
template<uint dependent, int dim, typename NumberType , typename Vector , size_t n_fe_functions_dep>
void DiFfRG::def::AbstractModel< Model, Components_ >::ldg_flux ( std::array< Tensor< 1, dim, NumberType >, n_fe_functions_dep > & F,
const Point< dim > & x,
const Vector & u ) const
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
dependentthe index \(i\) of the dependent variable \(l_i\) which is constructed from the previous level \(l_{i-1}\).
Parameters
Fthe 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.
xa d-dimensional dealii::Point<dim> representing field coordinates.
uthe field values of \(l_j(x)\) at the point x.

◆ ldg_source()

template<typename Model , typename Components_ >
template<uint dependent, int dim, typename NumberType , typename Vector , size_t n_fe_functions_dep>
void DiFfRG::def::AbstractModel< Model, Components_ >::ldg_source ( std::array< NumberType, n_fe_functions_dep > & s,
const Point< dim > & x,
const Vector & u ) const
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
dependentthe index \(i\) of the dependent variable \(l_i\) which is constructed from the previous level \(l_{i-1}\).
Parameters
sthe 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.
xa d-dimensional dealii::Point<dim> representing field coordinates.
uthe field values of \(l_j(x)\) at the point x.

◆ mass() [1/2]

template<typename Model , typename Components_ >
template<int dim, typename NumberType , typename Vector , typename Vector_dot , size_t n_fe_functions>
void DiFfRG::def::AbstractModel< Model, Components_ >::mass ( std::array< NumberType, n_fe_functions > & m_i,
const Point< dim > & x,
const Vector & u_i,
const Vector_dot & dt_u_i ) const
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_ithe resulting mass function \(m_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation.
xa d-dimensional dealii::Point<dim> representing field coordinates.
u_ithe field values \(u_i(x)\) at the point x.
dt_u_ithe time derivative of the field values \(\partial_t u_i(x)\) at the point x.

◆ mass() [2/2]

template<typename Model , typename Components_ >
template<int dim, typename NumberType , size_t n_fe_functions>
void DiFfRG::def::AbstractModel< Model, Components_ >::mass ( std::array< std::array< NumberType, n_fe_functions >, n_fe_functions > & m_ij,
const Point< dim > & x ) const
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_ijthe 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.
xa d-dimensional dealii::Point<dim> representing field coordinates.

◆ readouts()

template<typename Model , typename Components_ >
template<int dim, typename DataOut , typename Solutions >
void DiFfRG::def::AbstractModel< Model, Components_ >::readouts ( DataOut & output,
const Point< dim > & x,
const Solutions & sol ) const
inline

◆ readouts_multiple()

template<typename Model , typename Components_ >
template<typename FUN , typename DataOut >
void DiFfRG::def::AbstractModel< Model, Components_ >::readouts_multiple ( FUN & helper,
DataOut &  ) const
inline

◆ source()

template<typename Model , typename Components_ >
template<int dim, typename NumberType , typename Solutions , size_t n_fe_functions>
void DiFfRG::def::AbstractModel< Model, Components_ >::source ( std::array< NumberType, n_fe_functions > & s_i,
const Point< dim > & x,
const Solutions & sol ) const
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_ithe resulting source function \(s_i\), with \(N_f\) components. This method should fill this argument with the desired structure of the flow equation.
xa d-dimensional dealii::Point<dim> representing field coordinates.
sola std::tuple<...> which contains
  1. the array u_j
  2. the array of arrays \(\partial_x u_j\)
  3. the array of arrays of arrays \(\partial_x^2 u_j\)
  4. the array of extractors \(e_b\)

Member Data Documentation

◆ m_components

template<typename Model , typename Components_ >
Components_ DiFfRG::def::AbstractModel< Model, Components_ >::m_components
protected

The documentation for this class was generated from the following file: