ryujin 2.1.1 revision 7b33fa74945ccba789651fca714670471a34aef7
|
#include <source/shallow_water/riemann_solver.h>
Public Member Functions | |
Compute wavespeed estimates | |
RiemannSolver (const HyperbolicSystem &hyperbolic_system, const Parameters ¶meters, const PrecomputedVector &precomputed_values) | |
Number | compute (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) 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 |
Typedefs and constexpr constants | |
using | View = HyperbolicSystemView< dim, Number > |
using | ScalarNumber = typename View::ScalarNumber |
using | state_type = typename View::state_type |
using | primitive_type = typename std::array< Number, riemann_data_size > |
using | precomputed_type = typename View::precomputed_type |
using | PrecomputedVector = typename View::PrecomputedVector |
using | Parameters = RiemannSolverParameters< ScalarNumber > |
static constexpr auto | problem_dimension = View::problem_dimension |
static constexpr unsigned int | riemann_data_size = 3 |
Internal functions used in the Riemann solver | |
Number | compute_h_star (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const |
Number | f (const primitive_type &primitive_state, const Number &h_star) const |
Number | phi (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number &h) const |
Number | lambda1_minus (const primitive_type &riemann_data, const Number h_star) const |
Number | lambda3_plus (const primitive_type &riemann_data, const Number h_star) const |
Number | compute_lambda (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number h_star) const |
primitive_type | riemann_data_from_state (const state_type &U, const dealii::Tensor< 1, dim, Number > &n_ij) const |
A fast approximative solver for the associated 1D Riemann problem. The solver has to ensure that the estimate \(\lambda_{\text{max}}\) that is returned for the maximal wavespeed is a strict upper bound.
Definition at line 41 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::View = HyperbolicSystemView<dim, Number> |
Definition at line 49 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::ScalarNumber = typename View::ScalarNumber |
Definition at line 51 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::state_type = typename View::state_type |
Definition at line 55 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::primitive_type = typename std::array<Number, riemann_data_size> |
The array type to store the expanded primitive state for the Riemann solver \([\rho, v, p, a]\)
Definition at line 67 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::precomputed_type = typename View::precomputed_type |
Definition at line 69 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::PrecomputedVector = typename View::PrecomputedVector |
Definition at line 71 of file riemann_solver.h.
using ryujin::ShallowWater::RiemannSolver< dim, Number >::Parameters = RiemannSolverParameters<ScalarNumber> |
Definition at line 73 of file riemann_solver.h.
|
inline |
Constructor taking a HyperbolicSystem instance as argument
Definition at line 84 of file riemann_solver.h.
|
inline |
For two given 1D primitive states riemann_data_i and riemann_data_j, compute an estimation of an upper bound for the maximum wavespeed lambda.
Definition at line 227 of file riemann_solver.template.h.
Number ryujin::ShallowWater::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 |
For two given states U_i a U_j and a (normalized) "direction" n_ij compute an estimation of an upper bound for lambda.
Definition at line 241 of file riemann_solver.template.h.
|
inlineprotected |
Definition at line 27 of file riemann_solver.template.h.
|
inlineprotected |
Definition at line 49 of file riemann_solver.template.h.
|
inlineprotected |
Definition at line 66 of file riemann_solver.template.h.
References ryujin::positive_part().
|
inlineprotected |
Definition at line 81 of file riemann_solver.template.h.
References ryujin::positive_part().
|
inlineprotected |
Definition at line 96 of file riemann_solver.template.h.
References ryujin::negative_part(), and ryujin::positive_part().
|
inline |
Definition at line 110 of file riemann_solver.template.h.
References ryujin::positive_part().
|
inlineprotected |
Definition at line 209 of file riemann_solver.template.h.
|
staticconstexpr |
Definition at line 53 of file riemann_solver.h.
|
staticconstexpr |
Number of components in a primitive state, we store \([\rho, v, p, a]\), thus, 4.
Definition at line 61 of file riemann_solver.h.