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

#include <source/scalar_conservation/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
 
Access to runtime parameters
DEAL_II_ALWAYS_INLINE const std::string & flux () const
 
DEAL_II_ALWAYS_INLINE ScalarNumber derivative_approximation_delta () const
 
Low-level access to the flux function parser:
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_function (const Number &u) const
 
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_gradient_function (const Number &u) const
 
Special functions for boundary states
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
 
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 &) 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 = 1
 
static const auto component_names
 
static const auto primitive_component_names
 
static constexpr unsigned int n_precomputed_values = 2 * dim
 
static const auto precomputed_names
 
static constexpr unsigned int n_initial_precomputed_values = 0
 
static const auto initial_precomputed_names
 

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
 
dealii::Tensor< 1, dim, Number > construct_flux_tensor (const precomputed_type &precomputed_state) const
 
dealii::Tensor< 1, dim, Number > construct_flux_gradient_tensor (const precomputed_type &precomputed_state) const
 
flux_contribution_type flux_contribution (const PrecomputedVector &pv, const InitialPrecomputedVector &, const unsigned int i, const state_type &) const
 
flux_contribution_type flux_contribution (const PrecomputedVector &pv, const InitialPrecomputedVector &, const unsigned int *js, const state_type &) 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 &, const flux_contribution_type &, 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 square_entropy (const Number &u) const
 
Number square_entropy_derivative (const Number &u) const
 
Number kruzkov_entropy (const Number &k, const Number &u) const
 
Number kruzkov_entropy_derivative (const Number &k, const Number &u) const
 
bool is_admissible (const state_type &) const
 
static Number state (const state_type &U)
 

Detailed Description

template<int dim, typename Number>
class ryujin::ScalarConservation::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 85 of file hyperbolic_system.h.

Member Typedef Documentation

◆ ScalarNumber

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

The underlying scalar number type.

Definition at line 109 of file hyperbolic_system.h.

◆ state_type

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

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

Definition at line 170 of file hyperbolic_system.h.

◆ flux_type

template<int dim, typename Number >
using ryujin::ScalarConservation::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 175 of file hyperbolic_system.h.

◆ flux_contribution_type

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

The storage type used for flux contributions.

Definition at line 181 of file hyperbolic_system.h.

◆ precomputed_type

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

Array type used for precomputed values.

Definition at line 205 of file hyperbolic_system.h.

◆ initial_precomputed_type

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

Array type used for precomputed initial values.

Definition at line 229 of file hyperbolic_system.h.

◆ StateVector

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

A compound state vector.

Definition at line 241 of file hyperbolic_system.h.

◆ HyperbolicVector

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

MulticomponentVector for storing the hyperbolic state vector:

Definition at line 247 of file hyperbolic_system.h.

◆ PrecomputedVector

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

MulticomponentVector for storing a vector of precomputed states:

Definition at line 253 of file hyperbolic_system.h.

◆ InitialPrecomputedVector

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

MulticomponentVector for storing a vector of precomputed initial states:

Definition at line 260 of file hyperbolic_system.h.

Constructor & Destructor Documentation

◆ HyperbolicSystemView()

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

Constructor taking a reference to the underlying HyperbolicSystem

Definition at line 92 of file hyperbolic_system.h.

Member Function Documentation

◆ view()

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

Create a modified view from the current one:

Definition at line 101 of file hyperbolic_system.h.

◆ flux()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE const std::string & ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::flux ( ) const
inline

◆ derivative_approximation_delta()

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

◆ flux_function()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::flux_function ( const Number &  u) const
inline

For a given state \(u\) compute the flux \(f(u)\).

Definition at line 554 of file hyperbolic_system.h.

◆ flux_gradient_function()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::flux_gradient_function ( const Number &  u) const
inline

For a given state \(u\) compute the flux gradient \(f'(u)\).

Definition at line 576 of file hyperbolic_system.h.

◆ precomputation_loop()

template<int dim, typename Number >
template<typename DISPATCH , typename SPARSITY >
void ryujin::ScalarConservation::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

◆ state()

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

Return the scalar value stored in the state "vector". Given the fact that a scalar conservation equation is indeed scalar this simply unwraps the Tensor from the state_type and returns the one and only entry.

Definition at line 647 of file hyperbolic_system.h.

Referenced by ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::apply_galilei_transform(), ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::expand_state(), and ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::to_primitive_state().

◆ square_entropy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::square_entropy ( const Number &  u) const
inline

For a given state U, compute the square entropy

\[ \eta = 1/2 u^2. \]

Definition at line 655 of file hyperbolic_system.h.

◆ square_entropy_derivative()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::square_entropy_derivative ( const Number &  u) const
inline

For a given state U, compute the derivative of the square entropy

\[ \eta' = u. \]

Definition at line 663 of file hyperbolic_system.h.

◆ kruzkov_entropy()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::kruzkov_entropy ( const Number &  k,
const Number &  u 
) const
inline

For a given state U, compute the Krŭzkov entropy

\[ \eta = |u-k|. \]

Definition at line 672 of file hyperbolic_system.h.

◆ kruzkov_entropy_derivative()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE Number ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::kruzkov_entropy_derivative ( const Number &  k,
const Number &  u 
) const
inline

For a given state U, compute the derivative of the Krŭzkov entropy:

\[ \eta' = \text{sgn}(u-k). \]

Definition at line 681 of file hyperbolic_system.h.

◆ is_admissible()

template<int dim, typename Number >
bool ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::is_admissible ( const state_type ) 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 341 of file hyperbolic_system.h.

◆ apply_boundary_conditions()

template<int dim, typename Number >
template<typename Lambda >
state_type ryujin::ScalarConservation::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.

◆ construct_flux_tensor()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::construct_flux_tensor ( const precomputed_type precomputed_state) const
inline

Helper function that given a precomputed_state constructs a dealii::Tensor with the flux \(f(u)\).

Definition at line 745 of file hyperbolic_system.h.

◆ construct_flux_gradient_tensor()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::construct_flux_gradient_tensor ( const precomputed_type precomputed_state) const
inline

Helper function that given a precomputed_state constructs a dealii::Tensor with the gradient of the flux \(f(u)\) with respect to the state \(u\).

Definition at line 772 of file hyperbolic_system.h.

◆ flux_contribution() [1/2]

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::flux_contribution ( const PrecomputedVector pv,
const InitialPrecomputedVector ,
const unsigned int  i,
const state_type  
) 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 &, const unsigned int i, const state_type &) 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 799 of file hyperbolic_system.h.

◆ flux_contribution() [2/2]

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

Definition at line 815 of file hyperbolic_system.h.

◆ flux_divergence()

template<int dim, typename Number >
DEAL_II_ALWAYS_INLINE auto ryujin::ScalarConservation::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 831 of file hyperbolic_system.h.

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

◆ high_order_flux_divergence()

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

◆ nodal_source() [1/2]

template<int dim, typename Number >
state_type ryujin::ScalarConservation::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::ScalarConservation::HyperbolicSystemView< dim, Number >::nodal_source ( const PrecomputedVector pv,
const unsigned int *  js,
const state_type U_j,
const ScalarNumber  tau 
) const
delete

◆ expand_state()

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

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.

Definition at line 467 of file hyperbolic_system.h.

References ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::state().

◆ from_primitive_state()

template<int dim, typename Number >
state_type ryujin::ScalarConservation::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 476 of file hyperbolic_system.h.

◆ to_primitive_state()

template<int dim, typename Number >
state_type ryujin::ScalarConservation::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 485 of file hyperbolic_system.h.

References ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::state().

◆ apply_galilei_transform()

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

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

Definition at line 496 of file hyperbolic_system.h.

References ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::state().

Member Data Documentation

◆ problem_dimension

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

The dimension of the state space.

Definition at line 165 of file hyperbolic_system.h.

◆ component_names

template<int dim, typename Number >
const auto ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::component_names
inlinestatic
Initial value:
=
std::array<std::string, problem_dimension>{"u"}

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

Definition at line 187 of file hyperbolic_system.h.

◆ primitive_component_names

template<int dim, typename Number >
const auto ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::primitive_component_names
inlinestatic
Initial value:
=
std::array<std::string, problem_dimension>{"u"}

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

Definition at line 194 of file hyperbolic_system.h.

◆ n_precomputed_values

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

The number of precomputed values.

Definition at line 200 of file hyperbolic_system.h.

◆ precomputed_names

template<int dim, typename Number >
const auto ryujin::ScalarConservation::HyperbolicSystemView< dim, Number >::precomputed_names
inlinestatic
Initial value:
=
[]() -> std::array<std::string, n_precomputed_values> {
if constexpr (dim == 1)
return {"f", "df"};
else if constexpr (dim == 2)
return {"f_1", "f_2", "df_1", "df_2"};
else if constexpr (dim == 3)
return {"f_1", "f_2", "f_3", "df_1", "df_2", "df_3"};
__builtin_trap();
}()

An array holding all component names of the precomputed values.

Definition at line 210 of file hyperbolic_system.h.

◆ n_initial_precomputed_values

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

The number of precomputed initial values.

Definition at line 224 of file hyperbolic_system.h.

◆ initial_precomputed_names

template<int dim, typename Number >
const auto ryujin::ScalarConservation::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 235 of file hyperbolic_system.h.

◆ n_precomputation_cycles

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

The number of precomputation cycles.

Definition at line 274 of file hyperbolic_system.h.

◆ have_high_order_flux

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

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

Definition at line 424 of file hyperbolic_system.h.

◆ have_source_terms

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

We do not have source terms

Definition at line 438 of file hyperbolic_system.h.


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