ryujin 2.1.1 revision 0f269b6e84acefa08132307953bd597a510add24
Public Types | Public Member Functions | Static Public Attributes | List of all members
ryujin::HyperbolicModule< Description, dim, Number > Class Template Referencefinal

#include <source/hyperbolic_module.h>

Inheritance diagram for ryujin::HyperbolicModule< Description, dim, Number >:
Inheritance graph
[legend]
Collaboration diagram for ryujin::HyperbolicModule< Description, dim, Number >:
Collaboration graph
[legend]

Public Types

using HyperbolicSystem = typename Description::HyperbolicSystem
 
using View = typename Description::template HyperbolicSystemView< dim, Number >
 
using state_type = typename View::state_type
 
using scalar_type = typename OfflineData< dim, Number >::scalar_type
 
using vector_type = typename View::vector_type
 
using precomputed_vector_type = typename View::precomputed_vector_type
 
using precomputed_initial_vector_type = typename View::precomputed_initial_vector_type
 

Public Member Functions

 HyperbolicModule (const MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="/HyperbolicModule")
 
void prepare ()
 

Static Public Attributes

static constexpr unsigned int problem_dimension = View::problem_dimension
 
static constexpr unsigned int n_precomputed_values
 
static constexpr unsigned int n_precomputation_cycles
 
static constexpr unsigned int n_precomputed_initial_values
 

Functons for performing explicit time steps

bool precompute_only_
 
IDViolationStrategy id_violation_strategy_
 
template<int stages>
Number step (const vector_type &old_U, std::array< std::reference_wrapper< const vector_type >, stages > stage_U, std::array< std::reference_wrapper< const precomputed_vector_type >, stages > stage_precomputed, const std::array< Number, stages > stage_weights, vector_type &new_U, precomputed_vector_type &new_precomputed, Number tau=Number(0.)) const
 
void apply_boundary_conditions (vector_type &U, Number t) const
 
void cfl (Number new_cfl) const
 
auto & cfl () const
 
auto & offline_data () const
 
auto & hyperbolic_system () const
 
auto & precomputed_initial () const
 
auto & alpha () const
 
auto & n_restarts () const
 
auto & n_warnings () const
 

Detailed Description

template<typename Description, int dim, typename Number = double>
class ryujin::HyperbolicModule< Description, dim, Number >

Explicit forward Euler time-stepping for hyperbolic systems with convex limiting.

This module is described in detail in [13], Alg. 1.

Definition at line 75 of file hyperbolic_module.h.

Member Typedef Documentation

◆ HyperbolicSystem

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem

Definition at line 81 of file hyperbolic_module.h.

◆ View

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number>

Definition at line 86 of file hyperbolic_module.h.

◆ state_type

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::state_type = typename View::state_type

Definition at line 97 of file hyperbolic_module.h.

◆ scalar_type

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::scalar_type = typename OfflineData<dim, Number>::scalar_type

Shorthand typedef for dealii::LinearAlgebra::distributed::Vector<Number>.

Definition at line 102 of file hyperbolic_module.h.

◆ vector_type

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::vector_type = typename View::vector_type

MulticomponentVector for storing a vector of conserved states:

Definition at line 107 of file hyperbolic_module.h.

◆ precomputed_vector_type

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::precomputed_vector_type = typename View::precomputed_vector_type

Typedef for a MultiComponentVector storing precomputed values.

Definition at line 124 of file hyperbolic_module.h.

◆ precomputed_initial_vector_type

template<typename Description , int dim, typename Number = double>
using ryujin::HyperbolicModule< Description, dim, Number >::precomputed_initial_vector_type = typename View::precomputed_initial_vector_type

Typedef for a MultiComponentVector storing precomputed initial_values.

Definition at line 135 of file hyperbolic_module.h.

Constructor & Destructor Documentation

◆ HyperbolicModule()

template<typename Description , int dim, typename Number >
ryujin::HyperbolicModule< Description, dim, Number >::HyperbolicModule ( const MPI_Comm &  mpi_communicator,
std::map< std::string, dealii::Timer > &  computing_timer,
const OfflineData< dim, Number > &  offline_data,
const HyperbolicSystem hyperbolic_system,
const InitialValues< Description, dim, Number > &  initial_values,
const std::string &  subsection = "/HyperbolicModule< Description, dim, Number >" 
)

Constructor.

Definition at line 28 of file hyperbolic_module.template.h.

References ryujin::warn.

Member Function Documentation

◆ prepare()

template<typename Description , int dim, typename Number >
void ryujin::HyperbolicModule< Description, dim, Number >::prepare

Prepare time stepping. A call to prepare() allocates temporary storage and is necessary before any of the following time-stepping functions can be called.

Definition at line 59 of file hyperbolic_module.template.h.

◆ step()

template<typename Description , int dim, typename Number >
template<int stages>
Number ryujin::HyperbolicModule< Description, dim, Number >::step ( const vector_type old_U,
std::array< std::reference_wrapper< const vector_type >, stages >  stage_U,
std::array< std::reference_wrapper< const precomputed_vector_type >, stages >  stage_precomputed,
const std::array< Number, stages >  stage_weights,
vector_type new_U,
precomputed_vector_type new_precomputed,
Number  tau = Number(0.) 
) const

Given a reference to a previous state vector old_U perform an explicit euler step (and store the result in new_U). The function returns the computed maximal time step size tau_max according to the CFL condition.

The time step is performed with either tau_max (if tau is set to 0), or tau (if tau is nonzero). Here, tau_max is the computed maximal time step size and tau is the last parameter of the function.

The function takes an optional array of states stage_U together with a an array of weights stage_weights to construct a modified high-order flux. The standard high-order flux reads (cf [13], Eq. 12):

\begin{align} \newcommand{\bF}{{\boldsymbol F}} \newcommand{\bU}{{\boldsymbol U}} \newcommand\bUni{\bU^n_i} \newcommand\bUnj{\bU^n_j} \newcommand{\polf}{{\mathbb f}} \newcommand\Ii{\mathcal{I}(i)} \newcommand{\bc}{{\boldsymbol c}} \sum_{j\in\Ii} \frac{m_{ij}}{m_{j}} \; \frac{m_{j}}{\tau_n}\big( \tilde\bU_j^{H,n+1} - \bU_j^{n}\big) \;=\; \bF^n_i + \sum_{j\in\Ii}d_{ij}^{H,n}\big(\bUnj-\bUni\big), \qquad\text{with}\quad \bF^n_i\;:=\; \sum_{j\in\Ii}\Big(-(\polf(\bUni)+\polf(\bUnj)) \cdot\bc_{ij}\Big). \end{align}

Instead, the function assembles the modified high-order flux:

\begin{align} \newcommand{\bF}{{\boldsymbol F}} \newcommand{\bU}{{\boldsymbol U}} \newcommand\bUnis{\bU^{s,n}_i} \newcommand\bUnjs{\bU^{s,n}_j} \newcommand{\polf}{{\mathbb f}} \newcommand\Ii{\mathcal{I}(i)} \newcommand{\bc}{{\boldsymbol c}} \tilde{\bF}^n_i\;:=\; \big(1-\sum_{s=\{1:\text{stages}\}}\omega_s\big)\bF^n_i \;+\; \sum_{s=\{1:stages\}}\omega_s \bF^{s,n}_i \qquad\text{with}\quad \bF^{s,n}_i\;:=\; \sum_{j\in\Ii}\Big(-(\polf(\bUnis)+\polf(\bUnjs)) \cdot\bc_{ij}\Big). \end{align}

where \(\omega_s\) denotes the weigths for the given stages \(\bU^{s,n}\).

Note
The routine does not automatically update ghost vectors of the distributed vector new_U. It is best to simply call HyperbolicModule::apply_boundary_conditions() on the appropriate vector immediately after performing a time step.

Definition at line 128 of file hyperbolic_module.template.h.

◆ apply_boundary_conditions()

template<typename Description , int dim, typename Number >
void ryujin::HyperbolicModule< Description, dim, Number >::apply_boundary_conditions ( vector_type U,
Number  t 
) const

This function postprocesses a given state U to conform with all prescribed boundary conditions at time t. This implies that on slip (and no-slip) boundaries the normal momentum is set to zero; on Dirichlet boundaries the appropriate state at time t is substituted; and on "flexible" boundaries depending on the fact whether we have supersonic or subsonic inflow/outflow the appropriate Riemann invariant is prescribed. See [12] for details.

Note
The routine does update ghost vectors of the distributed vector U

Definition at line 1033 of file hyperbolic_module.template.h.

References ryujin::do_nothing, LIKWID_MARKER_START, and LIKWID_MARKER_STOP.

◆ cfl() [1/2]

template<typename Description , int dim, typename Number = double>
void ryujin::HyperbolicModule< Description, dim, Number >::cfl ( Number  new_cfl) const
inline

Sets the relative CFL number used for computing an appropriate time-step size to the given value. The CFL number must be a positive value. If chosen to be within the interval \((0,1)\) then the low-order update and limiting stages guarantee invariant domain preservation.

Definition at line 253 of file hyperbolic_module.h.

◆ cfl() [2/2]

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::cfl ( ) const
inline

Returns the relative CFL number used for computing an appropriate time-step size.

Definition at line 263 of file hyperbolic_module.h.

◆ offline_data()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::offline_data ( ) const
inline

Return a reference to the OfflineData object

Definition at line 268 of file hyperbolic_module.h.

◆ hyperbolic_system()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::hyperbolic_system ( ) const
inline

Return a reference to the HyperbolicSystem object

Definition at line 273 of file hyperbolic_module.h.

◆ precomputed_initial()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::precomputed_initial ( ) const
inline

Return a reference to the precomputed initial data vector

Definition at line 278 of file hyperbolic_module.h.

◆ alpha()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::alpha ( ) const
inline

Return a reference to alpha vector storing indicator values. Note that the values stored in alpha correspond to the last step executed by this class. This value can be recomputed for a given state vector by setting precompute_only_ to true and calling the step() function.

Definition at line 286 of file hyperbolic_module.h.

◆ n_restarts()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::n_restarts ( ) const
inline

The number of restarts issued by the step() function.

Definition at line 291 of file hyperbolic_module.h.

◆ n_warnings()

template<typename Description , int dim, typename Number = double>
auto & ryujin::HyperbolicModule< Description, dim, Number >::n_warnings ( ) const
inline

The number of ID violation warnings encounterd in the step() function.

Definition at line 297 of file hyperbolic_module.h.

Member Data Documentation

◆ problem_dimension

template<typename Description , int dim, typename Number = double>
constexpr unsigned int ryujin::HyperbolicModule< Description, dim, Number >::problem_dimension = View::problem_dimension
staticconstexpr

Definition at line 92 of file hyperbolic_module.h.

◆ n_precomputed_values

template<typename Description , int dim, typename Number = double>
constexpr unsigned int ryujin::HyperbolicModule< Description, dim, Number >::n_precomputed_values
staticconstexpr
Initial value:
=
View::n_precomputed_values

Definition at line 112 of file hyperbolic_module.h.

◆ n_precomputation_cycles

template<typename Description , int dim, typename Number = double>
constexpr unsigned int ryujin::HyperbolicModule< Description, dim, Number >::n_precomputation_cycles
staticconstexpr
Initial value:
=
View::n_precomputation_cycles

The number of precomputation cycles.

Definition at line 118 of file hyperbolic_module.h.

◆ n_precomputed_initial_values

template<typename Description , int dim, typename Number = double>
constexpr unsigned int ryujin::HyperbolicModule< Description, dim, Number >::n_precomputed_initial_values
staticconstexpr
Initial value:
=
View::n_precomputed_initial_values

Definition at line 129 of file hyperbolic_module.h.

◆ precompute_only_

template<typename Description , int dim, typename Number = double>
bool ryujin::HyperbolicModule< Description, dim, Number >::precompute_only_
mutable

Definition at line 300 of file hyperbolic_module.h.

◆ id_violation_strategy_

template<typename Description , int dim, typename Number = double>
IDViolationStrategy ryujin::HyperbolicModule< Description, dim, Number >::id_violation_strategy_
mutable

Definition at line 303 of file hyperbolic_module.h.


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