ryujin 2.1.1 revision 6dc06e5864abd5d99e5d7ab641dbe621936411d9
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 Member Functions

Constructor and setup
 HyperbolicModule (const MPIEnsemble &mpi_ensemble, 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 ()
 

Typedefs and constexpr constants

using HyperbolicSystem = typename Description::HyperbolicSystem
 
using View = typename Description::template HyperbolicSystemView< dim, Number >
 
using state_type = typename View::state_type
 
using precomputed_type = typename View::precomputed_type
 
using initial_precomputed_type = typename View::initial_precomputed_type
 
using StateVector = typename View::StateVector
 
using InitialPrecomputedVector = typename View::InitialPrecomputedVector
 
static constexpr auto problem_dimension = View::problem_dimension
 
static constexpr auto n_precomputation_cycles
 

Functons for performing explicit time steps

IDViolationStrategy id_violation_strategy_
 
void prepare_state_vector (StateVector &state_vector, Number t) const
 
template<int stages>
Number step (const StateVector &old_state_vector, std::array< std::reference_wrapper< const StateVector >, stages > stage_state_vectors, const std::array< Number, stages > stage_weights, StateVector &new_state_vector, Number tau=Number(0.), std::atomic< Number > tau_max=std::numeric_limits< Number >::max()) const
 
void cfl (Number new_cfl) const
 
const auto & cfl () const
 
const auto & offline_data () const
 
const auto & hyperbolic_system () const
 
const auto & initial_precomputed () const
 
const auto & alpha () const
 
const auto & n_restarts () const
 
const auto & n_corrections () const
 
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 104 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 112 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 114 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 119 of file hyperbolic_module.h.

◆ precomputed_type

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

Definition at line 121 of file hyperbolic_module.h.

◆ initial_precomputed_type

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

Definition at line 123 of file hyperbolic_module.h.

◆ StateVector

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

Definition at line 125 of file hyperbolic_module.h.

◆ InitialPrecomputedVector

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

Definition at line 127 of file hyperbolic_module.h.

Constructor & Destructor Documentation

◆ HyperbolicModule()

template<typename Description , int dim, typename Number >
ryujin::HyperbolicModule< Description, dim, Number >::HyperbolicModule ( const MPIEnsemble mpi_ensemble,
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 29 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 55 of file hyperbolic_module.template.h.

◆ prepare_state_vector()

template<typename Description , int dim, typename Number >
void ryujin::HyperbolicModule< Description, dim, Number >::prepare_state_vector ( StateVector state_vector,
Number  t 
) const

This function preprocesses a given state vector U in preparation for an explicit euler step performed by the step() function. The function performs the following tasks:

  • For a continuous finite element ansatz the method updates the U component of the state vector by enforcing boundary conditions for the supplied time time t. It then updates ghost ranges on U so that the state vector is consistent across MPI ranks.
  • For a discontinuous finite element ansatz it populates a local boundary state vector that is used for computing the boundary jump terms in the step() function when performing a dG update. It then updates ghost ranges on U so that the state vector is consistent across MPI ranks.
  • The function then runs the precomputation loop that populates the "precomputed values" component of the state vector and distributes the result over all MPI ranks by updating ghost ranges of the precomputed values vector.

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

References ryujin::do_nothing, and LIKWID_MARKER_START.

◆ step()

template<typename Description , int dim, typename Number >
template<int stages>
Number ryujin::HyperbolicModule< Description, dim, Number >::step ( const StateVector old_state_vector,
std::array< std::reference_wrapper< const StateVector >, stages >  stage_state_vectors,
const std::array< Number, stages >  stage_weights,
StateVector new_state_vector,
Number  tau = Number(0.),
std::atomic< Number >  tau_max = std::numeric_limits<Number>::max() 
) 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 chosen time step size tau with which the update was performed.

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 minimum of the specified parameter tau_max and the computed maximal time step size according to the CFL condition. 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 only performs an explicit update step on the locally owned dof index range. It neither updates the precomputed block of the state vector, sets boundary conditions (prior) to the update step, nor automatically updates the ghost range of the vector. It is thus necessary to call HyperbolicModule::prepare_state_vector() on old_state_vector prior to calling the step function.

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

◆ 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 263 of file hyperbolic_module.h.

◆ cfl() [2/2]

template<typename Description , int dim, typename Number = double>
const 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 273 of file hyperbolic_module.h.

◆ offline_data()

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

Return a reference to the OfflineData object

Definition at line 278 of file hyperbolic_module.h.

◆ hyperbolic_system()

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

Return a reference to the HyperbolicSystem object

Definition at line 283 of file hyperbolic_module.h.

◆ initial_precomputed()

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

Return a reference to the precomputed initial data vector

Definition at line 288 of file hyperbolic_module.h.

◆ alpha()

template<typename Description , int dim, typename Number = double>
const 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.

Definition at line 295 of file hyperbolic_module.h.

◆ n_restarts()

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

The number of restarts signalled by the step() function. We signal a restart whenever we encounter an ID violation in the low order update.

Definition at line 302 of file hyperbolic_module.h.

◆ n_corrections()

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

The number of corrections performed by the step() function. This function exists to mirror the ParabolicModule interface and will always return 0.

Definition at line 309 of file hyperbolic_module.h.

◆ n_warnings()

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

The number of ID violation warnings encounterd in the step() function. We issue a warning whenever we encounter an ID violation in the low order update (and throwing a restart is disabled).

Definition at line 316 of file hyperbolic_module.h.

Member Data Documentation

◆ problem_dimension

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

Definition at line 117 of file hyperbolic_module.h.

◆ n_precomputation_cycles

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

Definition at line 129 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 319 of file hyperbolic_module.h.


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