ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
|
#include <source/hyperbolic_module.h>
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 |
auto & | cfl () const |
auto & | offline_data () const |
auto & | hyperbolic_system () const |
auto & | initial_precomputed () const |
auto & | alpha () const |
auto & | n_restarts () const |
auto & | n_warnings () const |
Explicit forward Euler time-stepping for hyperbolic systems with convex limiting.
This module is described in detail in [13], Alg. 1.
Definition at line 79 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 87 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number> |
Definition at line 89 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::state_type = typename View::state_type |
Definition at line 94 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::precomputed_type = typename View::precomputed_type |
Definition at line 96 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::initial_precomputed_type = typename View::initial_precomputed_type |
Definition at line 98 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::StateVector = typename View::StateVector |
Definition at line 100 of file hyperbolic_module.h.
using ryujin::HyperbolicModule< Description, dim, Number >::InitialPrecomputedVector = typename View::InitialPrecomputedVector |
Definition at line 102 of file hyperbolic_module.h.
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 >" |
||
) |
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 54 of file hyperbolic_module.template.h.
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:
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.U
so that the state vector is consistent across MPI ranks.Definition at line 98 of file hyperbolic_module.template.h.
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}\).
old_state_vector
prior to calling the step function. Definition at line 237 of file hyperbolic_module.template.h.
References ryujin::ScalarConservation::Indicator< dim, Number >::accumulate(), ryujin::ScalarConservation::Limiter< dim, Number >::accumulate(), ryujin::ScalarConservation::Indicator< dim, Number >::alpha(), ryujin::ScalarConservation::Limiter< dim, Number >::bounds(), CALLGRIND_START_INSTRUMENTATION, CALLGRIND_STOP_INSTRUMENTATION, ryujin::ScalarConservation::Limiter< dim, Number >::combine_bounds(), ryujin::ScalarConservation::RiemannSolver< dim, Number >::compute(), LIKWID_MARKER_START, LIKWID_MARKER_STOP, ryujin::ScalarConservation::Limiter< dim, Number >::limit(), ryujin::raise_exception, ryujin::ScalarConservation::Indicator< dim, Number >::reset(), ryujin::ScalarConservation::Limiter< dim, Number >::reset(), RYUJIN_OMP_FOR, RYUJIN_PARALLEL_REGION_BEGIN, RYUJIN_PARALLEL_REGION_END, and ryujin::warn.
|
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 238 of file hyperbolic_module.h.
|
inline |
Returns the relative CFL number used for computing an appropriate time-step size.
Definition at line 248 of file hyperbolic_module.h.
|
inline |
Return a reference to the OfflineData object
Definition at line 253 of file hyperbolic_module.h.
|
inline |
Return a reference to the HyperbolicSystem object
Definition at line 258 of file hyperbolic_module.h.
|
inline |
Return a reference to the precomputed initial data vector
Definition at line 263 of file hyperbolic_module.h.
|
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 270 of file hyperbolic_module.h.
|
inline |
The number of restarts issued by the step() function.
Definition at line 275 of file hyperbolic_module.h.
|
inline |
The number of ID violation warnings encounterd in the step() function.
Definition at line 281 of file hyperbolic_module.h.
|
staticconstexpr |
Definition at line 92 of file hyperbolic_module.h.
|
staticconstexpr |
Definition at line 104 of file hyperbolic_module.h.
|
mutable |
Definition at line 284 of file hyperbolic_module.h.