ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
|
#include <source/euler/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 = 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 = 4 |
Internal functions used in the Riemann solver | |
Number | f (const primitive_type &riemann_data, const Number p_star) const |
Number | df (const primitive_type &riemann_data, const Number &p_star) const |
Number | phi (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number p_in) const |
Number | dphi (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number &p) const |
Number | phi_of_p_max (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const |
Number | lambda1_minus (const primitive_type &riemann_data, const Number p_star) const |
Number | lambda3_plus (const primitive_type &primitive_state, const Number p_star) const |
std::array< Number, 2 > | compute_gap (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number p_1, const Number p_2) const |
Number | compute_lambda (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number p_star) const |
Number | p_star_two_rarefaction (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const |
Number | p_star_failsafe (const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) 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 1D Riemann problem. The solver ensures that the estimate \(\lambda_{\text{max}}\) that is returned for the maximal wavespeed is a strict upper bound.
The solver is based on [5].
Definition at line 60 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::View = HyperbolicSystemView<dim, Number> |
Definition at line 68 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::ScalarNumber = typename View::ScalarNumber |
Definition at line 70 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::state_type = typename View::state_type |
Definition at line 74 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::primitive_type = 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 86 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::precomputed_type = typename View::precomputed_type |
Definition at line 88 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::PrecomputedVector = typename View::PrecomputedVector |
Definition at line 90 of file riemann_solver.h.
using ryujin::Euler::RiemannSolver< dim, Number >::Parameters = RiemannSolverParameters<ScalarNumber> |
Definition at line 92 of file riemann_solver.h.
|
inline |
Constructor taking a HyperbolicSystem instance as argument
Definition at line 103 of file riemann_solver.h.
Number ryujin::Euler::RiemannSolver< dim, Number >::compute | ( | const primitive_type & | riemann_data_i, |
const primitive_type & | riemann_data_j | ||
) | const |
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 407 of file riemann_solver.template.h.
References ryujin::quadratic_newton_step().
|
inline |
For two given states U_i a U_j and a (normalized) "direction" n_ij compute an estimation of an upper bound for lambda.
Returns a tuple consisting of lambda max and the number of Newton iterations used in the solver to find it.
Definition at line 586 of file riemann_solver.template.h.
|
inlineprotected |
See [5], page 912, (3.4).
Cost: 1x pow, 1x division, 2x sqrt
Definition at line 23 of file riemann_solver.template.h.
References ryujin::pow().
|
inlineprotected |
See [5], page 912, (3.4).
Cost: 1x pow, 3x division, 1x sqrt
Definition at line 51 of file riemann_solver.template.h.
References ryujin::pow().
|
inlineprotected |
See [5], page 912, (3.3).
Cost: 2x pow, 6x division, 2x sqrt
Definition at line 89 of file riemann_solver.template.h.
|
inlineprotected |
See [5], page 912, (3.3).
Cost: 2x pow, 6x division, 2x sqrt
Definition at line 102 of file riemann_solver.template.h.
|
inlineprotected |
See [5], page 912, (3.3).
The approximate Riemann solver is based on a function phi(p) that is montone increasing in p, concave down and whose (weak) third derivative is non-negative and locally bounded [1, p. 912]. Because we actually do not perform any iteration for computing our wavespeed estimate we can get away by only implementing a specialized variant of the phi function that computes phi(p_max). It inlines the implementation of the "f" function and eliminates all unnecessary branches in "f".
Cost: 0x pow, 2x division, 2x sqrt
Definition at line 124 of file riemann_solver.template.h.
|
inlineprotected |
see [5], page 912, (3.7)
Cost: 0x pow, 1x division, 1x sqrt
Definition at line 166 of file riemann_solver.template.h.
References ryujin::positive_part().
|
inlineprotected |
see [5], page 912, (3.8)
Cost: 0x pow, 1x division, 1x sqrt
Definition at line 191 of file riemann_solver.template.h.
References ryujin::positive_part().
|
inlineprotected |
For two given primitive states riemann_data_i
and riemann_data_j
, and two guesses p_1 <= p* <= p_2, compute the gap in lambda between both guesses.
See [5], page 914, (4.4a), (4.4b), (4.5), and (4.6)
Cost: 0x pow, 4x division, 4x sqrt
For two given primitive states riemann_data_i
and riemann_data_j
, and two guesses p_1 <= p* <= p_2, compute the gap in lambda between both guesses.
See [1], page 914, (4.4a), (4.4b), (4.5), and (4.6)
Cost: 0x pow, 4x division, 4x sqrt
Definition at line 219 of file riemann_solver.template.h.
References ryujin::negative_part(), and ryujin::positive_part().
|
inlineprotected |
see [5], page 912, (3.9)
For two given primitive states riemann_data_i
and riemann_data_j
, and a guess p_2, compute an upper bound for lambda.
Cost: 0x pow, 2x division, 2x sqrt (inclusive)
Definition at line 254 of file riemann_solver.template.h.
References ryujin::negative_part(), and ryujin::positive_part().
|
inlineprotected |
Two-rarefaction approximation to p_star computed for two primitive states riemann_data_i
and riemann_data_j
.
See [5], page 914, (4.3)
Cost: 2x pow, 2x division, 0x sqrt
Definition at line 276 of file riemann_solver.template.h.
References ryujin::positive_part(), and ryujin::pow().
|
inlineprotected |
Failsafe approximation to p_star computed for two primitive states riemann_data_i
and riemann_data_j
.
See [2], (5.11):
Cost: 0x pow, 3x division, 3x sqrt
Definition at line 332 of file riemann_solver.template.h.
|
inlineprotected |
For a given (2+dim dimensional) state vector U
, and a (normalized) "direction" n_ij, first compute the corresponding projected state in the corresponding 1D Riemann problem, and then compute and return the Riemann data [rho, u, p, a] (used in the approximative Riemann solver).
Definition at line 379 of file riemann_solver.template.h.
|
staticconstexpr |
Definition at line 72 of file riemann_solver.h.
|
staticconstexpr |
Number of components in a primitive state, we store \([\rho, v, p, a]\), thus, 4.
Definition at line 80 of file riemann_solver.h.