ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
Public Types | Public Member Functions | List of all members
ryujin::Euler::HyperbolicSystemView< dim, Number > Class Template Reference

#include <source/euler/hyperbolic_system.h>

Public Types

using ScalarNumber = typename get_value_type< Number >::type
 

Public Member Functions

 HyperbolicSystemView (const HyperbolicSystem &hyperbolic_system)
 
template<int dim2, typename Number2 >
auto view () const
 
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void precomputation_loop (unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
 
template<int component>
DEAL_II_ALWAYS_INLINE auto linearized_eigenvector (const state_type &U, const dealii::Tensor< 1, dim, Number > &normal) const -> std::array< state_type, 2 >
 
template<int component>
DEAL_II_ALWAYS_INLINE auto prescribe_riemann_characteristic (const state_type &U, const state_type &U_bar, const dealii::Tensor< 1, dim, Number > &normal) const -> state_type
 
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto apply_boundary_conditions (dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const -> state_type
 
template<typename ST >
auto expand_state (const ST &state) const -> state_type
 
template<typename ST >
DEAL_II_ALWAYS_INLINE auto from_initial_state (const ST &initial_state) const -> state_type
 
template<typename Lambda >
auto apply_galilei_transform (const state_type &state, const Lambda &lambda) const -> state_type
 
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void precomputation_loop (unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
 
template<int component>
DEAL_II_ALWAYS_INLINE auto prescribe_riemann_characteristic (const state_type &U, const Number &p, const state_type &U_bar, const Number &p_bar, const dealii::Tensor< 1, dim, Number > &normal) const -> state_type
 
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto apply_boundary_conditions (dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const -> state_type
 
template<typename ST >
auto expand_state (const ST &state) const -> state_type
 
template<typename ST >
DEAL_II_ALWAYS_INLINE auto from_initial_state (const ST &initial_state) const -> state_type
 
template<typename Lambda >
auto apply_galilei_transform (const state_type &state, const Lambda &lambda) const -> state_type
 
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void precomputation_loop (unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
 
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto apply_boundary_conditions (dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &, const Lambda &get_dirichlet_data) const -> state_type
 
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void precomputation_loop (unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
 
template<int component>
DEAL_II_ALWAYS_INLINE auto prescribe_riemann_characteristic (const state_type &U, const state_type &U_bar, const dealii::Tensor< 1, dim, Number > &normal) const -> state_type
 
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto apply_boundary_conditions (const dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const -> state_type
 
template<typename ST >
DEAL_II_ALWAYS_INLINE auto expand_state (const ST &state) const -> state_type
 
template<typename ST >
DEAL_II_ALWAYS_INLINE auto from_initial_state (const ST &initial_state) const -> state_type
 
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto apply_galilei_transform (const state_type &state, const Lambda &lambda) const -> state_type
 
Access to runtime parameters
DEAL_II_ALWAYS_INLINE ScalarNumber gamma () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber reference_density () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber vacuum_state_relaxation_small () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber vacuum_state_relaxation_large () const
 
Special functions for boundary states
template<int component>
std::array< state_type, 2 > linearized_eigenvector (const state_type &U, const dealii::Tensor< 1, dim, Number > &normal) const
 
template<int component>
state_type prescribe_riemann_characteristic (const state_type &U, const state_type &U_bar, const dealii::Tensor< 1, dim, Number > &normal) const
 
template<typename Lambda >
state_type apply_boundary_conditions (const dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const
 
State transformations
template<typename ST >
state_type expand_state (const ST &state) const
 
template<typename ST >
state_type from_initial_state (const ST &initial_state) const
 
state_type from_primitive_state (const state_type &primitive_state) const
 
state_type to_primitive_state (const state_type &state) const
 
template<typename Lambda >
state_type apply_galilei_transform (const state_type &state, const Lambda &lambda) const
 

Types and constexpr constants

using state_type = dealii::Tensor< 1, problem_dimension, Number >
 
using flux_type = dealii::Tensor< 1, problem_dimension, dealii::Tensor< 1, dim, Number > >
 
using flux_contribution_type = flux_type
 
using precomputed_type = std::array< Number, n_precomputed_values >
 
using initial_precomputed_type = std::array< Number, n_initial_precomputed_values >
 
using StateVector = Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values >
 
using HyperbolicVector = Vectors::MultiComponentVector< ScalarNumber, problem_dimension >
 
using PrecomputedVector = Vectors::MultiComponentVector< ScalarNumber, n_precomputed_values >
 
using InitialPrecomputedVector = Vectors::MultiComponentVector< ScalarNumber, n_initial_precomputed_values >
 
static constexpr unsigned int problem_dimension = 2 + dim
 
static const auto component_names
 
static const auto primitive_component_names
 
static constexpr unsigned int n_precomputed_values = 2
 
static const auto precomputed_names
 
static constexpr unsigned int n_initial_precomputed_values = 0
 
static const auto initial_precomputed_names
 

Access to cached inverses

A collection of commonly used expressions with gamma that would otherwise need to be recomputed many times putting unnecessary pressure on the div/sqrt ALU unit.

static constexpr bool have_gamma = true
 
static constexpr bool have_eos_interpolation_b = false
 
DEAL_II_ALWAYS_INLINE ScalarNumber gamma_inverse () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber gamma_plus_one_inverse () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber gamma_minus_one_inverse () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber gamma_minus_one_over_gamma_plus_one () const
 

Computing precomputed quantities

static constexpr unsigned int n_precomputation_cycles = 1
 
template<typename DISPATCH , typename SPARSITY >
void precomputation_loop (unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
 

Flux computations

static constexpr bool have_high_order_flux = false
 
flux_type f (const state_type &U) const
 
flux_contribution_type flux_contribution (const PrecomputedVector &pv, const InitialPrecomputedVector &ipv, const unsigned int i, const state_type &U_i) const
 
flux_contribution_type flux_contribution (const PrecomputedVector &pv, const InitialPrecomputedVector &ipv, const unsigned int *js, const state_type &U_j) const
 
state_type flux_divergence (const flux_contribution_type &flux_i, const flux_contribution_type &flux_j, const dealii::Tensor< 1, dim, Number > &c_ij) const
 
state_type high_order_flux_divergence (const flux_contribution_type &flux_i, const flux_contribution_type &flux_j, const dealii::Tensor< 1, dim, Number > &c_ij) const =delete
 

Computing stencil source terms

static constexpr bool have_source_terms = false
 
state_type nodal_source (const PrecomputedVector &pv, const unsigned int i, const state_type &U_i, const ScalarNumber tau) const =delete
 
state_type nodal_source (const PrecomputedVector &pv, const unsigned int *js, const state_type &U_j, const ScalarNumber tau) const =delete
 

Computing derived physical quantities

Number filter_vacuum_density (const Number &rho) const
 
Number pressure (const state_type &U) const
 
Number speed_of_sound (const state_type &U) const
 
Number specific_entropy (const state_type &U) const
 
Number harten_entropy (const state_type &U) const
 
state_type harten_entropy_derivative (const state_type &U) const
 
Number mathematical_entropy (const state_type &U) const
 
state_type mathematical_entropy_derivative (const state_type &U) const
 
bool is_admissible (const state_type &U) const
 
static Number density (const state_type &U)
 
static dealii::Tensor< 1, dim, Number > momentum (const state_type &U)
 
static Number total_energy (const state_type &U)
 
static Number internal_energy (const state_type &U)
 
static state_type internal_energy_derivative (const state_type &U)
 

Detailed Description

template<int dim, typename Number>
class ryujin::Euler::HyperbolicSystemView< dim, Number >

A view of the HyperbolicSystem that makes methods available 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.

Intended usage:

HyperbolicSystem hyperbolic_system;
const auto view = hyperbolic_system.template view<dim, Number>();
const auto flux_i = view.flux_contribution(...);
const auto flux_j = view.flux_contribution(...);
const auto flux_ij = view.flux_divergence(flux_i, flux_j, c_ij);
// etc.

Definition at line 104 of file hyperbolic_system.h.

Member Typedef Documentation

◆ ScalarNumber

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::ScalarNumber = typename get_value_type<Number>::type

The underlying scalar number type.

Definition at line 128 of file hyperbolic_system.h.

◆ state_type

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::state_type = dealii::Tensor<1, problem_dimension, Number>

Storage type for a (conserved) state vector \(\boldsymbol U\).

Definition at line 223 of file hyperbolic_system.h.

◆ flux_type

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::flux_type = dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number> >

Storage type for the flux \(\mathbf{f}\).

Definition at line 228 of file hyperbolic_system.h.

◆ flux_contribution_type

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::flux_contribution_type = flux_type

The storage type used for flux contributions.

Definition at line 234 of file hyperbolic_system.h.

◆ precomputed_type

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputed_type = std::array<Number, n_precomputed_values>

Array type used for precomputed values.

Definition at line 274 of file hyperbolic_system.h.

◆ initial_precomputed_type

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::initial_precomputed_type = std::array<Number, n_initial_precomputed_values>

Array type used for precomputed initial values.

Definition at line 290 of file hyperbolic_system.h.

◆ StateVector

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::StateVector = Vectors:: StateVector<ScalarNumber, problem_dimension, n_precomputed_values>

A compound state vector.

Definition at line 302 of file hyperbolic_system.h.

◆ HyperbolicVector

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::HyperbolicVector = Vectors::MultiComponentVector<ScalarNumber, problem_dimension>

MulticomponentVector for storing the hyperbolic state vector:

Definition at line 308 of file hyperbolic_system.h.

◆ PrecomputedVector

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::PrecomputedVector = Vectors::MultiComponentVector<ScalarNumber, n_precomputed_values>

MulticomponentVector for storing a vector of precomputed states:

Definition at line 314 of file hyperbolic_system.h.

◆ InitialPrecomputedVector

template<int dim, typename Number >
using ryujin::Euler::HyperbolicSystemView< dim, Number >::InitialPrecomputedVector = Vectors::MultiComponentVector<ScalarNumber, n_initial_precomputed_values>

MulticomponentVector for storing a vector of precomputed initial states:

Definition at line 321 of file hyperbolic_system.h.

Constructor & Destructor Documentation

◆ HyperbolicSystemView()

template<int dim, typename Number >
ryujin::Euler::HyperbolicSystemView< dim, Number >::HyperbolicSystemView ( const HyperbolicSystem hyperbolic_system)
inline

Constructor taking a reference to the underlying HyperbolicSystem

Definition at line 111 of file hyperbolic_system.h.

Member Function Documentation

◆ view()

template<int dim, typename Number >
template<int dim2, typename Number2 >
auto ryujin::Euler::HyperbolicSystemView< dim, Number >::view ( ) const
inline

Create a modified view from the current one:

Definition at line 120 of file hyperbolic_system.h.

◆ gamma()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::gamma ( ) const
inline

Definition at line 135 of file hyperbolic_system.h.

◆ reference_density()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::reference_density ( ) const
inline

Definition at line 140 of file hyperbolic_system.h.

◆ vacuum_state_relaxation_small()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::vacuum_state_relaxation_small ( ) const
inline

Definition at line 146 of file hyperbolic_system.h.

◆ vacuum_state_relaxation_large()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::vacuum_state_relaxation_large ( ) const
inline

Definition at line 152 of file hyperbolic_system.h.

◆ gamma_inverse()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::gamma_inverse ( ) const
inline

Definition at line 167 of file hyperbolic_system.h.

◆ gamma_plus_one_inverse()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::gamma_plus_one_inverse ( ) const
inline

Definition at line 172 of file hyperbolic_system.h.

◆ gamma_minus_one_inverse()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::gamma_minus_one_inverse ( ) const
inline

Definition at line 177 of file hyperbolic_system.h.

◆ gamma_minus_one_over_gamma_plus_one()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE ScalarNumber ryujin::Euler::HyperbolicSystemView< dim, Number >::gamma_minus_one_over_gamma_plus_one ( ) const
inline

Definition at line 183 of file hyperbolic_system.h.

◆ precomputation_loop() [1/5]

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
void ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputation_loop ( unsigned int  cycle,
const DISPATCH &  dispatch_check,
const SPARSITY &  sparsity_simd,
StateVector state_vector,
unsigned int  left,
unsigned int  right 
) const

Step 0: precompute values for hyperbolic update. This routine is called within our usual loop() idiom in HyperbolicModule

◆ density()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::density ( const state_type U)
inlinestatic

For a given (2+dim dimensional) state vector U, return the density U[0]

Definition at line 742 of file hyperbolic_system.h.

◆ filter_vacuum_density()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::filter_vacuum_density ( const Number &  rho) const
inline

Given a density rho this function returns 0 if the magniatude of rho is smaller or equal than relaxation_large * rho_cutoff. Otherwise rho is returned unmodified. Here, rho_cutoff is the reference density multiplied by eps.

Definition at line 750 of file hyperbolic_system.h.

◆ momentum()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > ryujin::Euler::HyperbolicSystemView< dim, Number >::momentum ( const state_type U)
inlinestatic

For a given (2+dim dimensional) state vector U, return the momentum vector [U[1], ..., U[1+dim]].

Definition at line 764 of file hyperbolic_system.h.

◆ total_energy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::total_energy ( const state_type U)
inlinestatic

For a given (2+dim dimensional) state vector U, return the total energy U[1+dim]

Definition at line 775 of file hyperbolic_system.h.

◆ internal_energy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::internal_energy ( const state_type U)
inlinestatic

For a given (2+dim dimensional) state vector U, compute and return the internal energy \(\varepsilon = (\rho e)\).

Definition at line 783 of file hyperbolic_system.h.

◆ internal_energy_derivative()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::internal_energy_derivative ( const state_type U)
inlinestatic

For a given (2+dim dimensional) state vector U, compute and return the derivative of the internal energy \(\varepsilon = (\rho e)\).

Definition at line 797 of file hyperbolic_system.h.

◆ pressure()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::pressure ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the pressure \(p\).

We assume that the pressure is given by a polytropic equation of state, i.e.,

\[ p = (\gamma - 1)\;(\rho e) \]

Definition at line 824 of file hyperbolic_system.h.

◆ speed_of_sound()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::speed_of_sound ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute the (physical) speed of sound:

\[ c^2 = \frac{\gamma\,p}{\rho} \]

Definition at line 833 of file hyperbolic_system.h.

◆ specific_entropy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::specific_entropy ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the (scaled) specific entropy

\[ e^{(\gamma-1)s} = \frac{\rho\,e}{\rho^\gamma}. \]

Definition at line 844 of file hyperbolic_system.h.

References ryujin::pow().

◆ harten_entropy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::harten_entropy ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the Harten-type entropy

\[ \eta = (\rho^2 e) ^ {1 / (\gamma + 1)}. \]

Definition at line 855 of file hyperbolic_system.h.

References ryujin::pow().

◆ harten_entropy_derivative()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::harten_entropy_derivative ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the derivative \(\eta'\) of the Harten-type entropy

\[ \eta = (\rho^2 e) ^ {1 / (\gamma + 1)}. \]

Definition at line 870 of file hyperbolic_system.h.

References ryujin::pow().

◆ mathematical_entropy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::Euler::HyperbolicSystemView< dim, Number >::mathematical_entropy ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the entropy \(\eta = p^{1/\gamma}\).

Definition at line 907 of file hyperbolic_system.h.

References ryujin::pow().

◆ mathematical_entropy_derivative()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::mathematical_entropy_derivative ( const state_type U) const
inline

For a given (2+dim dimensional) state vector U, compute and return the derivative \(\eta'\) of the entropy \(\eta = p^{1/\gamma}\).

Definition at line 918 of file hyperbolic_system.h.

References ryujin::pow().

◆ is_admissible()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE bool ryujin::Euler::HyperbolicSystemView< dim, Number >::is_admissible ( const state_type U) const
inline

Returns whether the state U is admissible. If U is a vectorized state then U is admissible if all vectorized values are admissible.

Definition at line 955 of file hyperbolic_system.h.

◆ linearized_eigenvector() [1/2]

template<int dim, typename Number >
template<int component>
std::array< state_type, 2 > ryujin::Euler::HyperbolicSystemView< dim, Number >::linearized_eigenvector ( const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal 
) const

For a given state U and normal direction normal returns the n-th pair of left and right eigenvectors of the linearized normal flux.

◆ prescribe_riemann_characteristic() [1/4]

template<int dim, typename Number >
template<int component>
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::prescribe_riemann_characteristic ( const state_type U,
const state_type U_bar,
const dealii::Tensor< 1, dim, Number > &  normal 
) const

Decomposes a given state U into Riemann invariants and then replaces the first or second Riemann characteristic from the one taken from U_bar state. Note that the U_bar state is just the prescribed dirichlet values.

◆ apply_boundary_conditions() [1/5]

template<int dim, typename Number >
template<typename Lambda >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_boundary_conditions ( const dealii::types::boundary_id  id,
const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal,
const Lambda &  get_dirichlet_data 
) const

Apply boundary conditions.

For the compressible Euler equations we have:

  • Dirichlet boundary conditions by prescribing the return value of get_dirichlet_data() as is.
  • Slip boundary conditions where we remove the normal component of the momentum.
  • No slip boundary conditions where we set the momentum to 0.
  • "Dynamic boundary" conditions that prescribe different Riemann invariants from the return value of get_dirichlet_data() depending on the flow state (supersonic versus subsonic, outflow versus inflow).

◆ f()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::f ( const state_type U) const
inline

Given a state U compute the flux

\[ \begin{pmatrix} \textbf m \\ \textbf v\otimes \textbf m + p\mathbb{I}_d \\ \textbf v(E+p) \end{pmatrix}, \]

Definition at line 1173 of file hyperbolic_system.h.

◆ flux_contribution() [1/2]

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::flux_contribution ( const PrecomputedVector pv,
const InitialPrecomputedVector ipv,
const unsigned int  i,
const state_type U_i 
) const
inline

Given a state U_i and an index i compute flux contributions.

Intended usage:

for (unsigned int i = n_internal; i < n_owned; ++i) {
// ...
const auto flux_i = flux_contribution(precomputed..., i, U_i);
for (unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
// ...
const auto flux_j = flux_contribution(precomputed..., js, U_j);
const auto flux_ij = flux_divergence(flux_i, flux_j, c_ij);
}
}
flux_contribution_type flux_contribution(const PrecomputedVector &pv, const InitialPrecomputedVector &ipv, const unsigned int i, const state_type &U_i) const
state_type flux_divergence(const flux_contribution_type &flux_i, const flux_contribution_type &flux_j, const dealii::Tensor< 1, dim, Number > &c_ij) const

For the Euler equations we simply compute f(U_i).

Definition at line 1195 of file hyperbolic_system.h.

◆ flux_contribution() [2/2]

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::flux_contribution ( const PrecomputedVector pv,
const InitialPrecomputedVector ipv,
const unsigned int *  js,
const state_type U_j 
) const
inline

Definition at line 1207 of file hyperbolic_system.h.

◆ flux_divergence()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::flux_divergence ( const flux_contribution_type flux_i,
const flux_contribution_type flux_j,
const dealii::Tensor< 1, dim, Number > &  c_ij 
) const
inline

Given flux contributions flux_i and flux_j compute the flux (-f(U_i) - f(U_j)

Definition at line 1219 of file hyperbolic_system.h.

References ryujin::add(), and ryujin::contract().

◆ high_order_flux_divergence()

template<int dim, typename Number >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::high_order_flux_divergence ( const flux_contribution_type flux_i,
const flux_contribution_type flux_j,
const dealii::Tensor< 1, dim, Number > &  c_ij 
) const
delete

◆ nodal_source() [1/2]

template<int dim, typename Number >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::nodal_source ( const PrecomputedVector pv,
const unsigned int  i,
const state_type U_i,
const ScalarNumber  tau 
) const
delete

◆ nodal_source() [2/2]

template<int dim, typename Number >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::nodal_source ( const PrecomputedVector pv,
const unsigned int *  js,
const state_type U_j,
const ScalarNumber  tau 
) const
delete

◆ expand_state() [1/4]

template<int dim, typename Number >
template<typename ST >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::expand_state ( const ST &  state) const

Given a state vector associated with a different spatial dimensions than the current one, return an "expanded" version of the state vector associated with dim spatial dimensions where the momentum vector of the conserved state state is expaned with zeros to a total length of dim entries.

Note
dim has to be larger or equal than the dimension of the ST vector.

◆ from_initial_state() [1/4]

template<int dim, typename Number >
template<typename ST >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::from_initial_state ( const ST &  initial_state) const

Given an initial state [rho, u_1, ..., u_?, p] return a conserved state [rho, m_1, ..., m_d, E].

This function simply calls from_primitive_state() and expand_state().

Note
This function is used to conveniently convert (user provided) primitive initial states with pressure values to a conserved state in the EulerInitialStateLibrary. As such, this function is implemented in the Euler::HyperbolicSystem and EulerAEOS::HyperbolicSystem classes.

◆ from_primitive_state()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::from_primitive_state ( const state_type primitive_state) const
inline

Given a primitive state [rho, u_1, ..., u_d, p] return a conserved state

Definition at line 1264 of file hyperbolic_system.h.

◆ to_primitive_state()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::to_primitive_state ( const state_type state) const
inline

Given a conserved state return a primitive state [rho, u_1, ..., u_d, p]

Definition at line 1286 of file hyperbolic_system.h.

◆ apply_galilei_transform() [1/4]

template<int dim, typename Number >
template<typename Lambda >
state_type ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_galilei_transform ( const state_type state,
const Lambda &  lambda 
) const

Transform the current state according to a given operator lambda acting on a dim dimensional momentum (or velocity) vector.

◆ precomputation_loop() [2/5]

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputation_loop ( unsigned int  cycle,
const DISPATCH &  dispatch_check,
const SPARSITY &  sparsity_simd,
StateVector state_vector,
unsigned int  left,
unsigned int  right 
) const
inline

Definition at line 705 of file hyperbolic_system.h.

References RYUJIN_OMP_FOR.

◆ linearized_eigenvector() [2/2]

template<int dim, typename Number >
template<int component>
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::linearized_eigenvector ( const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal 
) const -> std::array<state_type, 2>
inline

Definition at line 985 of file hyperbolic_system.h.

◆ prescribe_riemann_characteristic() [2/4]

template<int dim, typename Number >
template<int component>
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::prescribe_riemann_characteristic ( const state_type U,
const state_type U_bar,
const dealii::Tensor< 1, dim, Number > &  normal 
) const -> state_type
inline

Definition at line 1040 of file hyperbolic_system.h.

References ryujin::pow().

◆ apply_boundary_conditions() [2/5]

template<int dim, typename Number >
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_boundary_conditions ( dealii::types::boundary_id  id,
const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal,
const Lambda &  get_dirichlet_data 
) const -> state_type
inline

◆ expand_state() [2/4]

template<int dim, typename Number >
template<typename ST >
auto ryujin::Euler::HyperbolicSystemView< dim, Number >::expand_state ( const ST &  state) const -> state_type

Definition at line 1230 of file hyperbolic_system.h.

◆ from_initial_state() [2/4]

template<int dim, typename Number >
template<typename ST >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::from_initial_state ( const ST &  initial_state) const -> state_type
inline

Definition at line 1254 of file hyperbolic_system.h.

◆ apply_galilei_transform() [2/4]

template<int dim, typename Number >
template<typename Lambda >
auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_galilei_transform ( const state_type state,
const Lambda &  lambda 
) const -> state_type

Definition at line 1306 of file hyperbolic_system.h.

◆ precomputation_loop() [3/5]

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputation_loop ( unsigned int  cycle,
const DISPATCH &  dispatch_check,
const SPARSITY &  sparsity_simd,
StateVector state_vector,
unsigned int  left,
unsigned int  right 
) const
inline

Definition at line 876 of file hyperbolic_system.h.

References RYUJIN_OMP_FOR, and RYUJIN_OMP_SINGLE.

◆ prescribe_riemann_characteristic() [3/4]

template<int dim, typename Number >
template<int component>
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::prescribe_riemann_characteristic ( const state_type U,
const Number &  p,
const state_type U_bar,
const Number &  p_bar,
const dealii::Tensor< 1, dim, Number > &  normal 
) const -> state_type
inline

Definition at line 1291 of file hyperbolic_system.h.

References ryujin::pow().

◆ apply_boundary_conditions() [3/5]

template<int dim, typename Number >
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_boundary_conditions ( dealii::types::boundary_id  id,
const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal,
const Lambda &  get_dirichlet_data 
) const -> state_type
inline

◆ expand_state() [3/4]

template<int dim, typename Number >
template<typename ST >
auto ryujin::Euler::HyperbolicSystemView< dim, Number >::expand_state ( const ST &  state) const -> state_type

Definition at line 1599 of file hyperbolic_system.h.

◆ from_initial_state() [3/4]

template<int dim, typename Number >
template<typename ST >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::from_initial_state ( const ST &  initial_state) const -> state_type
inline

Definition at line 1623 of file hyperbolic_system.h.

◆ apply_galilei_transform() [3/4]

template<int dim, typename Number >
template<typename Lambda >
auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_galilei_transform ( const state_type state,
const Lambda &  lambda 
) const -> state_type

Definition at line 1683 of file hyperbolic_system.h.

◆ precomputation_loop() [4/5]

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputation_loop ( unsigned int  cycle,
const DISPATCH &  dispatch_check,
const SPARSITY &  sparsity_simd,
StateVector state_vector,
unsigned int  left,
unsigned int  right 
) const
inline

Definition at line 600 of file hyperbolic_system.h.

References RYUJIN_OMP_FOR.

◆ apply_boundary_conditions() [4/5]

template<int dim, typename Number >
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_boundary_conditions ( dealii::types::boundary_id  id,
const state_type U,
const dealii::Tensor< 1, dim, Number > &  ,
const Lambda &  get_dirichlet_data 
) const -> state_type
inline

◆ precomputation_loop() [5/5]

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
DEAL_II_ALWAYS_INLINE void ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputation_loop ( unsigned int  cycle,
const DISPATCH &  dispatch_check,
const SPARSITY &  sparsity_simd,
StateVector state_vector,
unsigned int  left,
unsigned int  right 
) const
inline

Definition at line 667 of file hyperbolic_system.h.

References ryujin::pow(), and RYUJIN_OMP_FOR.

◆ prescribe_riemann_characteristic() [4/4]

template<int dim, typename Number >
template<int component>
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::prescribe_riemann_characteristic ( const state_type U,
const state_type U_bar,
const dealii::Tensor< 1, dim, Number > &  normal 
) const -> state_type
inline

Definition at line 888 of file hyperbolic_system.h.

◆ apply_boundary_conditions() [5/5]

template<int dim, typename Number >
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_boundary_conditions ( const dealii::types::boundary_id  id,
const state_type U,
const dealii::Tensor< 1, dim, Number > &  normal,
const Lambda &  get_dirichlet_data 
) const -> state_type
inline

◆ expand_state() [4/4]

template<int dim, typename Number >
template<typename ST >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::expand_state ( const ST &  state) const -> state_type
inline

Definition at line 1241 of file hyperbolic_system.h.

◆ from_initial_state() [4/4]

template<int dim, typename Number >
template<typename ST >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::from_initial_state ( const ST &  initial_state) const -> state_type
inline

Definition at line 1263 of file hyperbolic_system.h.

◆ apply_galilei_transform() [4/4]

template<int dim, typename Number >
template<typename Lambda >
DEAL_II_ALWAYS_INLINE auto ryujin::Euler::HyperbolicSystemView< dim, Number >::apply_galilei_transform ( const state_type state,
const Lambda &  lambda 
) const -> state_type
inline

Definition at line 1306 of file hyperbolic_system.h.

Member Data Documentation

◆ have_gamma

template<int dim, typename Number >
constexpr bool ryujin::Euler::HyperbolicSystemView< dim, Number >::have_gamma = true
staticconstexpr

constexpr boolean used in the EulerInitialStates namespace

Definition at line 192 of file hyperbolic_system.h.

◆ have_eos_interpolation_b

template<int dim, typename Number >
constexpr bool ryujin::Euler::HyperbolicSystemView< dim, Number >::have_eos_interpolation_b = false
staticconstexpr

constexpr boolean used in the EulerInitialStates namespace

Definition at line 197 of file hyperbolic_system.h.

◆ problem_dimension

template<int dim, typename Number >
constexpr unsigned int ryujin::Euler::HyperbolicSystemView< dim, Number >::problem_dimension = 2 + dim
staticconstexpr

The dimension of the state space.

Definition at line 218 of file hyperbolic_system.h.

◆ component_names

template<int dim, typename Number >
const auto ryujin::Euler::HyperbolicSystemView< dim, Number >::component_names
inlinestatic
Initial value:
=
[]() -> std::array<std::string, problem_dimension> {
if constexpr (dim == 1)
return {"rho", "m", "E"};
else if constexpr (dim == 2)
return {"rho", "m_1", "m_2", "E"};
else if constexpr (dim == 3)
return {"rho", "m_1", "m_2", "m_3", "E"};
__builtin_trap();
}()

An array holding all component names of the conserved state as a string.

Definition at line 240 of file hyperbolic_system.h.

◆ primitive_component_names

template<int dim, typename Number >
const auto ryujin::Euler::HyperbolicSystemView< dim, Number >::primitive_component_names
inlinestatic
Initial value:
=
[]() -> std::array<std::string, problem_dimension> {
if constexpr (dim == 1)
return {"rho", "v", "p"};
else if constexpr (dim == 2)
return {"rho", "v_1", "v_2", "p"};
else if constexpr (dim == 3)
return {"rho", "v_1", "v_2", "v_3", "p"};
__builtin_trap();
}()

An array holding all component names of the primitive state as a string.

Definition at line 255 of file hyperbolic_system.h.

◆ n_precomputed_values

template<int dim, typename Number >
constexpr unsigned int ryujin::Euler::HyperbolicSystemView< dim, Number >::n_precomputed_values = 2
staticconstexpr

The number of precomputed values.

Definition at line 269 of file hyperbolic_system.h.

◆ precomputed_names

template<int dim, typename Number >
const auto ryujin::Euler::HyperbolicSystemView< dim, Number >::precomputed_names
inlinestatic
Initial value:
=
std::array<std::string, n_precomputed_values>{"s", "eta_h"}

An array holding all component names of the precomputed values.

Definition at line 279 of file hyperbolic_system.h.

◆ n_initial_precomputed_values

template<int dim, typename Number >
constexpr unsigned int ryujin::Euler::HyperbolicSystemView< dim, Number >::n_initial_precomputed_values = 0
staticconstexpr

The number of precomputed initial values.

Definition at line 285 of file hyperbolic_system.h.

◆ initial_precomputed_names

template<int dim, typename Number >
const auto ryujin::Euler::HyperbolicSystemView< dim, Number >::initial_precomputed_names
inlinestatic
Initial value:
=
std::array<std::string, n_initial_precomputed_values>{}

An array holding all component names of the precomputed values.

Definition at line 296 of file hyperbolic_system.h.

◆ n_precomputation_cycles

template<int dim, typename Number >
constexpr unsigned int ryujin::Euler::HyperbolicSystemView< dim, Number >::n_precomputation_cycles = 1
staticconstexpr

The number of precomputation cycles.

Definition at line 334 of file hyperbolic_system.h.

◆ have_high_order_flux

template<int dim, typename Number >
constexpr bool ryujin::Euler::HyperbolicSystemView< dim, Number >::have_high_order_flux = false
staticconstexpr

The low-order and high-order fluxes are the same

Definition at line 573 of file hyperbolic_system.h.

◆ have_source_terms

template<int dim, typename Number >
constexpr bool ryujin::Euler::HyperbolicSystemView< dim, Number >::have_source_terms = false
staticconstexpr

We do not have source terms:

Definition at line 587 of file hyperbolic_system.h.


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