ryujin 2.1.1 revision d0a94ad2ccc0c4c2e8c2485c52b06b90e2fc9853
Public Types | Static Public Attributes | List of all members
ryujin::ScalarConservation::RiemannSolver< dim, Number > Class Template Reference

#include <source/scalar_conservation/riemann_solver.h>

Public Types

using View = HyperbolicSystemView< dim, Number >
 
using state_type = typename View::state_type
 
using precomputed_state_type = typename View::precomputed_state_type
 
using ScalarNumber = typename get_value_type< Number >::type
 
using Parameters = RiemannSolverParameters< ScalarNumber >
 

Static Public Attributes

static constexpr unsigned int n_precomputed_values
 

Compute wavespeed estimates

 RiemannSolver (const HyperbolicSystem &hyperbolic_system, const Parameters &parameters, const MultiComponentVector< ScalarNumber, n_precomputed_values > &precomputed_values)
 
Number compute (const Number &u_i, const Number &u_j, const precomputed_state_type &prec_i, const precomputed_state_type &prec_j, const dealii::Tensor< 1, dim, Number > &n_ij) const
 
Number compute (const state_type &U_i, const state_type &U_j, const unsigned int i, const unsigned int *js, const dealii::Tensor< 1, dim, Number > &n_ij) const
 

Detailed Description

template<int dim, typename Number = double>
class ryujin::ScalarConservation::RiemannSolver< dim, Number >

A fast estimate for a sufficient maximal wavespeed of the 1D Riemann problem. The wavespeed estimate is based on a guaranteed upper bound on the maximal wavespeed for convex fluxes, see Example 79.17 on page 333 of [GuermondErn2021]. As well as an augmented "Roe average" based on an entropy inequality of a suitable Krŭzkov entropy, see [7] Section 4.

Definition at line 75 of file riemann_solver.h.

Member Typedef Documentation

◆ View

template<int dim, typename Number = double>
using ryujin::ScalarConservation::RiemannSolver< dim, Number >::View = HyperbolicSystemView<dim, Number>

A view on the HyperbolicSystem for a given dimension dim and choice of number type Number (which can be a scalar float, or double, as well as a VectorizedArray holding packed scalars.

Definition at line 81 of file riemann_solver.h.

◆ state_type

template<int dim, typename Number = double>
using ryujin::ScalarConservation::RiemannSolver< dim, Number >::state_type = typename View::state_type

The storage type used for a (conserved) state vector \(\boldsymbol U\).

Definition at line 86 of file riemann_solver.h.

◆ precomputed_state_type

template<int dim, typename Number = double>
using ryujin::ScalarConservation::RiemannSolver< dim, Number >::precomputed_state_type = typename View::precomputed_state_type

Array type used for precomputed values.

Definition at line 97 of file riemann_solver.h.

◆ ScalarNumber

template<int dim, typename Number = double>
using ryujin::ScalarConservation::RiemannSolver< dim, Number >::ScalarNumber = typename get_value_type<Number>::type

The underlying scalar number type.

Definition at line 102 of file riemann_solver.h.

◆ Parameters

template<int dim, typename Number = double>
using ryujin::ScalarConservation::RiemannSolver< dim, Number >::Parameters = RiemannSolverParameters<ScalarNumber>

Definition at line 107 of file riemann_solver.h.

Constructor & Destructor Documentation

◆ RiemannSolver()

template<int dim, typename Number = double>
ryujin::ScalarConservation::RiemannSolver< dim, Number >::RiemannSolver ( const HyperbolicSystem hyperbolic_system,
const Parameters parameters,
const MultiComponentVector< ScalarNumber, n_precomputed_values > &  precomputed_values 
)
inline

Constructor taking a HyperbolicSystem instance as argument

Definition at line 117 of file riemann_solver.h.

Member Function Documentation

◆ compute() [1/2]

template<int dim, typename Number >
Number ryujin::ScalarConservation::RiemannSolver< dim, Number >::compute ( const Number &  u_i,
const Number &  u_j,
const precomputed_state_type prec_i,
const precomputed_state_type prec_j,
const dealii::Tensor< 1, dim, Number > &  n_ij 
) const

For two states u_i, u_j, precomputed values prec_i, prec_j, and a (normalized) "direction" n_ij compute an upper bound estimate for the wavespeed.

Definition at line 23 of file riemann_solver.template.h.

◆ compute() [2/2]

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::ScalarConservation::RiemannSolver< dim, Number >::compute ( const state_type U_i,
const state_type U_j,
const unsigned int  i,
const unsigned int *  js,
const dealii::Tensor< 1, dim, Number > &  n_ij 
) const
inline

For two given states U_i a U_j and a (normalized) "direction" n_ij compute an estimate for an upper bound of lambda.

Definition at line 196 of file riemann_solver.template.h.

Member Data Documentation

◆ n_precomputed_values

template<int dim, typename Number = double>
constexpr unsigned int ryujin::ScalarConservation::RiemannSolver< dim, Number >::n_precomputed_values
staticconstexpr
Initial value:

The number of precomputed values.

Definition at line 91 of file riemann_solver.h.


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