8#include <compile_time_options.h>
21 template <
int dim,
typename Number>
22 DEAL_II_ALWAYS_INLINE
inline Number
24 const Number p_star)
const
26 const auto view = hyperbolic_system.view<dim, Number>();
27 const auto &gamma = view.gamma();
29 const auto &[rho, u, p, a] = riemann_data;
31 const Number Az =
ScalarNumber(2.) / (rho * (gamma + Number(1.)));
34 const Number radicand = Az / (p_star + Bz);
35 const Number true_value = (p_star - p) * std::sqrt(radicand);
39 const Number factor =
ryujin::pow(p_star / p, exponent) - Number(1.);
40 const auto false_value =
43 return dealii::compare_and_apply_mask<
44 dealii::SIMDComparison::greater_than_or_equal>(
45 p_star, p, true_value, false_value);
49 template <
int dim,
typename Number>
50 DEAL_II_ALWAYS_INLINE
inline Number
52 const Number &p_star)
const
54 const auto view = hyperbolic_system.view<dim, Number>();
57 const auto &gamma = view.gamma();
58 const auto &gamma_inverse = view.gamma_inverse();
59 const auto &gamma_minus_one_inverse = view.gamma_minus_one_inverse();
60 const auto &gamma_plus_one_inverse = view.gamma_plus_one_inverse();
62 const auto &[rho, u, p, a] = riemann_data;
64 const Number radicand_inverse =
ScalarNumber(0.5) * rho *
67 const Number denominator =
68 (p_star + (gamma -
ScalarNumber(1.)) * gamma_plus_one_inverse * p);
69 const Number true_value =
71 (denominator * std::sqrt(radicand_inverse));
78 const auto false_value =
81 return dealii::compare_and_apply_mask<
82 dealii::SIMDComparison::greater_than_or_equal>(
83 p_star, p, true_value, false_value);
87 template <
int dim,
typename Number>
88 DEAL_II_ALWAYS_INLINE
inline Number
91 const Number p_in)
const
93 const Number &u_i = riemann_data_i[1];
94 const Number &u_j = riemann_data_j[1];
96 return f(riemann_data_i, p_in) + f(riemann_data_j, p_in) + u_j - u_i;
100 template <
int dim,
typename Number>
101 DEAL_II_ALWAYS_INLINE
inline Number
104 const Number &p)
const
106 return df(riemann_data_i, p) + df(riemann_data_j, p);
122 template <
int dim,
typename Number>
123 DEAL_II_ALWAYS_INLINE
inline Number
128 const auto view = hyperbolic_system.view<dim, Number>();
129 const auto &gamma = view.gamma();
131 const auto &[rho_i, u_i, p_i, a_i] = riemann_data_i;
132 const auto &[rho_j, u_j, p_j, a_j] = riemann_data_j;
134 const Number p_max = std::max(p_i, p_j);
136 const Number radicand_inverse_i =
ScalarNumber(0.5) * rho_i *
140 const Number value_i = (p_max - p_i) / std::sqrt(radicand_inverse_i);
142 const Number radicand_inverse_j =
ScalarNumber(0.5) * rho_j *
146 const Number value_j = (p_max - p_j) / std::sqrt(radicand_inverse_j);
148 return value_i + value_j + u_j - u_i;
164 template <
int dim,
typename Number>
165 DEAL_II_ALWAYS_INLINE
inline Number
169 const auto view = hyperbolic_system.view<dim, Number>();
170 const auto &gamma = view.gamma();
171 const auto &gamma_inverse = view.gamma_inverse();
175 const auto &[rho, u, p, a] = riemann_data;
180 return u - a * std::sqrt(
ScalarNumber(1.0) + factor * tmp);
189 template <
int dim,
typename Number>
190 DEAL_II_ALWAYS_INLINE
inline Number
194 const auto view = hyperbolic_system.view<dim, Number>();
195 const auto &gamma = view.gamma();
196 const auto &gamma_inverse = view.gamma_inverse();
197 const Number factor =
200 const auto &[rho, u, p, a] = primitive_state;
204 return u + a * std::sqrt(Number(1.0) + factor * tmp);
217 template <
int dim,
typename Number>
218 DEAL_II_ALWAYS_INLINE
inline std::array<Number, 2>
220 const std::array<Number, 4> &riemann_data_i,
221 const std::array<Number, 4> &riemann_data_j,
223 const Number p_2)
const
225 const Number nu_11 = lambda1_minus(riemann_data_i, p_2 );
226 const Number nu_12 = lambda1_minus(riemann_data_i, p_1 );
228 const Number nu_31 = lambda3_plus(riemann_data_j, p_1);
229 const Number nu_32 = lambda3_plus(riemann_data_j, p_2);
231 const Number lambda_max =
235 std::max(std::abs(nu_32 - nu_31), std::abs(nu_12 - nu_11));
237 return {{gap, lambda_max}};
252 template <
int dim,
typename Number>
253 DEAL_II_ALWAYS_INLINE
inline Number
257 const Number p_star)
const
259 const Number nu_11 = lambda1_minus(riemann_data_i, p_star);
260 const Number nu_32 = lambda3_plus(riemann_data_j, p_star);
274 template <
int dim,
typename Number>
275 DEAL_II_ALWAYS_INLINE
inline Number
280 const auto view = hyperbolic_system.view<dim, Number>();
281 const auto &gamma = view.gamma();
282 const auto &gamma_inverse = view.gamma_inverse();
283 const auto &gamma_minus_one_inverse = view.gamma_minus_one_inverse();
285 const auto &[rho_i, u_i, p_i, a_i] = riemann_data_i;
286 const auto &[rho_j, u_j, p_j, a_j] = riemann_data_j;
306 const Number numerator =
positive_part(a_i + a_j - factor * (u_j - u_i));
307 const Number denominator =
308 a_i *
ryujin::pow(p_i * inv_p_j, -factor * gamma_inverse) + a_j;
310 const auto exponent =
ScalarNumber(2.0) * gamma * gamma_minus_one_inverse;
312 const auto p_1_tilde =
313 p_j *
ryujin::pow(numerator / denominator, exponent);
315#ifdef DEBUG_RIEMANN_SOLVER
316 std::cout <<
"p_star_two_rarefaction = " << p_1_tilde << std::endl;
330 template <
int dim,
typename Number>
331 DEAL_II_ALWAYS_INLINE
inline Number
336 const auto view = hyperbolic_system.view<dim, Number>();
337 const auto &gamma = view.gamma();
339 const auto &[rho_i, u_i, p_i, a_i] = riemann_data_i;
340 const auto &[rho_j, u_j, p_j, a_j] = riemann_data_j;
348 const Number p_max = std::max(p_i, p_j);
352 rho_i * ((gamma + Number(1.)) * p_max + (gamma - Number(1.)) * p_i);
354 const Number x_i = std::sqrt(radicand_i);
358 rho_j * ((gamma + Number(1.)) * p_max + (gamma - Number(1.)) * p_j);
360 const Number x_j = std::sqrt(radicand_j);
362 const Number a = x_i + x_j;
363 const Number b = u_j - u_i;
364 const Number c = -p_i * x_i - p_j * x_j;
366 const Number base = (-b + std::sqrt(b * b -
ScalarNumber(4.) * a * c)) /
368 const Number p_2_tilde = base * base;
370#ifdef DEBUG_RIEMANN_SOLVER
371 std::cout <<
"p_star_failsafe = " << p_2_tilde << std::endl;
377 template <
int dim,
typename Number>
378 DEAL_II_ALWAYS_INLINE
inline auto
380 const state_type &U,
const dealii::Tensor<1, dim, Number> &n_ij)
const
383 const auto view = hyperbolic_system.view<dim, Number>();
385 const auto rho = view.density(U);
386 const auto rho_inverse = Number(1.0) / rho;
388 const auto m = view.momentum(U);
389 const auto proj_m = n_ij * m;
390 const auto perp = m - proj_m * n_ij;
393 view.total_energy(U) - Number(0.5) * perp.norm_square() * rho_inverse;
395 using state_type_1d =
397 const auto view_1d = hyperbolic_system.view<1, Number>();
399 const auto state = state_type_1d{{rho, proj_m, E}};
400 const auto p = view_1d.pressure(state);
401 const auto a = view_1d.speed_of_sound(state);
402 return {{rho, proj_m * rho_inverse, p, a}};
406 template <
int dim,
typename Number>
458 const auto &[rho_i, u_i, p_i, a_i] = riemann_data_i;
459 const auto &[rho_j, u_j, p_j, a_j] = riemann_data_j;
461#ifdef DEBUG_RIEMANN_SOLVER
462 std::cout <<
"rho_left: " << rho_i << std::endl;
463 std::cout <<
"u_left: " << u_i << std::endl;
464 std::cout <<
"p_left: " << p_i << std::endl;
465 std::cout <<
"a_left: " << a_i << std::endl;
466 std::cout <<
"rho_right: " << rho_j << std::endl;
467 std::cout <<
"u_right: " << u_j << std::endl;
468 std::cout <<
"p_right: " << p_j << std::endl;
469 std::cout <<
"a_right: " << a_j << std::endl;
472 const Number p_max = std::max(p_i, p_j);
474 const Number rarefaction =
475 p_star_two_rarefaction(riemann_data_i, riemann_data_j);
476 const Number failsafe = p_star_failsafe(riemann_data_i, riemann_data_j);
477 const Number p_star_tilde = std::min(rarefaction, failsafe);
479 const Number phi_p_max = phi_of_p_max(riemann_data_i, riemann_data_j);
482 dealii::compare_and_apply_mask<dealii::SIMDComparison::less_than>(
486 std::min(p_max, p_star_tilde));
488#ifdef DEBUG_RIEMANN_SOLVER
489 std::cout <<
" p^*_tilde = " << p_2 <<
"\n";
490 std::cout <<
" phi(p_*_t) = "
491 << phi(riemann_data_i, riemann_data_j, p_2) << std::endl;
498 if (parameters.newton_max_iterations() == 0) {
499 const auto lambda_max =
500 compute_lambda(riemann_data_i, riemann_data_j, p_2);
502#ifdef DEBUG_RIEMANN_SOLVER
503 std::cout <<
"-> lambda_max = " << lambda_max << std::endl;
514 const Number p_min = std::min(riemann_data_i[2], riemann_data_j[2]);
517 dealii::compare_and_apply_mask<dealii::SIMDComparison::less_than>(
518 phi_p_max, Number(0.), p_max, p_min);
520 p_1 = dealii::compare_and_apply_mask<
521 dealii::SIMDComparison::less_than_or_equal>(p_1, p_2, p_1, p_2);
529 auto [gap, lambda_max] =
530 compute_gap(riemann_data_i, riemann_data_j, p_1, p_2);
532#ifdef DEBUG_RIEMANN_SOLVER
533 std::cout << std::fixed << std::setprecision(16);
534 std::cout <<
"p_1: (start) " << p_1 << std::endl;
535 std::cout <<
"p_2: (start) " << p_2 << std::endl;
536 std::cout <<
"gap: (start) " << gap << std::endl;
537 std::cout <<
"l_m: (start) " << lambda_max << std::endl;
540 for (
unsigned int i = 0; i < parameters.newton_max_iterations(); ++i) {
543 const Number tolerance(parameters.newton_tolerance());
544 if (std::max(Number(0.), gap - tolerance) == Number(0.)) {
545#ifdef DEBUG_RIEMANN_SOLVER
546 std::cout <<
"converged after " << i <<
" iterations." << std::endl;
552 const Number phi_p_1 = phi(riemann_data_i, riemann_data_j, p_1);
553 const Number phi_p_2 = phi(riemann_data_i, riemann_data_j, p_2);
554 const Number dphi_p_1 = dphi(riemann_data_i, riemann_data_j, p_1);
555 const Number dphi_p_2 = dphi(riemann_data_i, riemann_data_j, p_2);
560 auto [gap_new, lambda_max_new] =
561 compute_gap(riemann_data_i, riemann_data_j, p_1, p_2);
563 lambda_max = lambda_max_new;
565#ifdef DEBUG_RIEMANN_SOLVER
566 std::cout <<
"phi_p_1: " << phi_p_1 << std::endl;
567 std::cout <<
"phi_p_2: " << phi_p_2 << std::endl;
568 std::cout <<
"dphi_p_1: " << dphi_p_1 << std::endl;
569 std::cout <<
"dphi_p_2: " << dphi_p_2 << std::endl;
570 std::cout <<
"p_1: ( " << i <<
" ) " << p_1 << std::endl;
571 std::cout <<
"p_2: ( " << i <<
" ) " << p_2 << std::endl;
572 std::cout <<
"gap: " << gap << std::endl;
573 std::cout <<
"l_m: " << lambda_max << std::endl;
577#ifdef DEBUG_RIEMANN_SOLVER
578 std::cout <<
"-> lambda_max = " << lambda_max << std::endl;
585 template <
int dim,
typename Number>
590 const unsigned int * ,
591 const dealii::Tensor<1, dim, Number> &n_ij)
const
593 const auto riemann_data_i = riemann_data_from_state(U_i, n_ij);
594 const auto riemann_data_j = riemann_data_from_state(U_j, n_ij);
596 return compute(riemann_data_i, riemann_data_j);
dealii::Tensor< 1, problem_dimension, Number > state_type
Number phi_of_p_max(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const
Number compute_lambda(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number p_star) const
typename View::ScalarNumber ScalarNumber
Number df(const primitive_type &riemann_data, const Number &p_star) const
Number lambda1_minus(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 p_star_two_rarefaction(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
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
typename View::state_type state_type
Number compute(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const
std::array< Number, riemann_data_size > primitive_type
Number p_star_failsafe(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const
Number dphi(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number &p) const
Number f(const primitive_type &riemann_data, const Number p_star) const
Number lambda3_plus(const primitive_type &primitive_state, const Number p_star) const
T pow(const T x, const T b)
DEAL_II_ALWAYS_INLINE Number negative_part(const Number number)
DEAL_II_ALWAYS_INLINE Number positive_part(const Number number)
DEAL_II_ALWAYS_INLINE void quadratic_newton_step(Number &p_1, Number &p_2, const Number phi_p_1, const Number phi_p_2, const Number dphi_p_1, const Number dphi_p_2, const Number sign=Number(1.0))