ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
riemann_solver.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0
3// [LANL Copyright Statement]
4// Copyright (C) 2022 - 2024 by the ryujin authors
5// Copyright (C) 2023 - 2024 by Triad National Security, LLC
6//
7
8#pragma once
9
10#include <compile_time_options.h>
11
12#include "riemann_solver.h"
13
14#include <newton.h>
15#include <simd.h>
16
17// #define DEBUG_RIEMANN_SOLVER
18
19namespace ryujin
20{
21 namespace ShallowWater
22 {
23 using namespace dealii;
24
25 template <int dim, typename Number>
26 DEAL_II_ALWAYS_INLINE inline Number
28 const Number &h) const
29 {
30 const auto view = hyperbolic_system.view<dim, Number>();
31 const ScalarNumber gravity = view.gravity();
32
33 const auto &[h_Z, u_Z, a_Z] = riemann_data_Z;
34
35 const auto left_value = ScalarNumber(2.) * (std::sqrt(gravity * h) - a_Z);
36
37 const Number radicand =
38 ScalarNumber(0.5) * gravity * (h + h_Z) / (h * h_Z);
39 const Number right_value = (h - h_Z) * std::sqrt(radicand);
40
41 return dealii::compare_and_apply_mask<
42 dealii::SIMDComparison::less_than_or_equal>(
43 h, h_Z, left_value, right_value);
44 }
45
46
47 template <int dim, typename Number>
48 DEAL_II_ALWAYS_INLINE inline Number
50 const primitive_type &riemann_data_j,
51 const Number &h) const
52 {
53 const Number &u_i = riemann_data_i[1];
54 const Number &u_j = riemann_data_j[1];
55
56#ifdef DEBUG_RIEMANN_SOLVER
57 std::cout << "f_L --> " << f(riemann_data_i, h) << std::endl;
58 std::cout << "f_R --> " << f(riemann_data_j, h) << std::endl;
59#endif
60 return f(riemann_data_i, h) + f(riemann_data_j, h) + u_j - u_i;
61 }
62
63
64 template <int dim, typename Number>
65 DEAL_II_ALWAYS_INLINE inline Number
67 const primitive_type &riemann_data, const Number h_star) const
68 {
69 const auto &[h, u, a] = riemann_data;
70
71 const Number factor = positive_part((h_star - h) / h);
72 const Number half_factor = ScalarNumber(0.5) * factor;
73
74 return u - a * std::sqrt((ScalarNumber(1.) + half_factor) *
75 (ScalarNumber(1.) + factor));
76 }
77
78
79 template <int dim, typename Number>
80 DEAL_II_ALWAYS_INLINE inline Number
82 const Number h_star) const
83 {
84 const auto &[h, u, a] = riemann_data;
85
86 const Number factor = positive_part((h_star - h) / h);
87 const Number half_factor = ScalarNumber(0.5) * factor;
88
89 return u + a * std::sqrt((ScalarNumber(1.) + half_factor) *
90 (ScalarNumber(1.) + factor));
91 }
92
93
94 template <int dim, typename Number>
95 DEAL_II_ALWAYS_INLINE inline Number
97 const primitive_type &riemann_data_i,
98 const primitive_type &riemann_data_j,
99 const Number h_star) const
100 {
101 const Number lambda1 = lambda1_minus(riemann_data_i, h_star);
102 const Number lambda3 = lambda3_plus(riemann_data_j, h_star);
103
104 return std::max(negative_part(lambda1), positive_part(lambda3));
105 }
106
107
108 template <int dim, typename Number>
109 DEAL_II_ALWAYS_INLINE inline Number
111 const primitive_type &riemann_data_i,
112 const primitive_type &riemann_data_j) const
113 {
114 const auto view = hyperbolic_system.view<dim, Number>();
115 const ScalarNumber gravity = view.gravity();
116 const auto gravity_inverse = ScalarNumber(1.) / gravity;
117
118 const auto &[h_i, u_i, a_i] = riemann_data_i;
119 const auto &[h_j, u_j, a_j] = riemann_data_j;
120
121 const Number h_min = std::min(h_i, h_j);
122 const Number h_max = std::max(h_i, h_j);
123
124#ifdef DEBUG_RIEMANN_SOLVER
125 std::cout << h_min << " <- h_min/max -> " << h_max << std::endl;
126#endif
127
128 const Number a_min = std::sqrt(gravity * h_min);
129 const Number a_max = std::sqrt(gravity * h_max);
130
131#ifdef DEBUG_RIEMANN_SOLVER
132 std::cout << a_min << " <- a_min/max -> " << a_max << std::endl;
133#endif
134
135 const Number sqrt_two = std::sqrt(ScalarNumber(2.));
136
137 /* x0 = (2 sqrt(2) - 1)^2 */
138 const Number x0 = Number(9.) - ScalarNumber(4.) * sqrt_two;
139
140 const Number phi_value_min =
141 phi(riemann_data_i, riemann_data_j, x0 * h_min);
142#ifdef DEBUG_RIEMANN_SOLVER
143 std::cout << "phi_value_min ->" << phi_value_min << std::endl;
144#endif
145
146 const Number phi_value_max =
147 phi(riemann_data_i, riemann_data_j, x0 * h_max);
148#ifdef DEBUG_RIEMANN_SOLVER
149 std::cout << "phi_value_max ->" << phi_value_max << std::endl;
150#endif
151
152 /* We compute the three h_star quantities */
153
154 Number tmp;
155
156 /* Double rarefaction case (h_star left): */
157
158 tmp = positive_part(u_i - u_j + ScalarNumber(2.) * (a_i + a_j));
159 const Number h_star_left =
160 ScalarNumber(0.0625) * gravity_inverse * tmp * tmp;
161
162#ifdef DEBUG_RIEMANN_SOLVER
163 std::cout << "left: " << h_star_left << std::endl;
164#endif
165
166 /* Double modified shock (h_star middle): */
167
168 tmp = Number(1.) + sqrt_two * (u_i - u_j) / (a_min + a_max);
169 const Number h_star_middle = std::sqrt(h_min * h_max) * tmp;
170
171#ifdef DEBUG_RIEMANN_SOLVER
172 std::cout << "middle: " << h_star_middle << std::endl;
173#endif
174
175 /* Expansion and modified shock (h_star right): */
176
177 const auto left_radicand =
178 ScalarNumber(3.) * h_min +
179 ScalarNumber(2.) * sqrt_two * std::sqrt(h_min * h_max);
180
181 const auto right_radicand =
182 sqrt_two * std::sqrt(gravity_inverse * h_min) * (u_i - u_j);
183
184 tmp = std::sqrt(positive_part(left_radicand + right_radicand));
185 tmp -= sqrt_two * std::sqrt(h_min);
186
187 const Number h_star_right = tmp * tmp;
188
189#ifdef DEBUG_RIEMANN_SOLVER
190 std::cout << "right: " << h_star_right << std::endl;
191#endif
192
193 /* Finally define h_star */
194
195 Number h_star = dealii::compare_and_apply_mask<
196 dealii::SIMDComparison::less_than_or_equal>(
197 Number(0.), phi_value_min, h_star_left, h_star_right);
198
199 h_star =
200 dealii::compare_and_apply_mask<dealii::SIMDComparison::less_than>(
201 phi_value_max, Number(0.), h_star_middle, h_star_right);
202
203 return h_star;
204 }
205
206
207 template <int dim, typename Number>
208 DEAL_II_ALWAYS_INLINE inline auto
210 const state_type &U, const dealii::Tensor<1, dim, Number> &n_ij) const
212 {
213 const auto view = hyperbolic_system.view<dim, Number>();
214
215 const Number h = view.water_depth_sharp(U);
216 const Number gravity = view.gravity();
217
218 const auto velocity = view.momentum(U) / h;
219 const auto projected_velocity = n_ij * velocity;
220 const auto a = std::sqrt(h * gravity);
221
222 return {{h, projected_velocity, a}};
223 }
224
225
226 template <int dim, typename Number>
228 const primitive_type &riemann_data_i,
229 const primitive_type &riemann_data_j) const
230 {
231 const Number h_star = compute_h_star(riemann_data_i, riemann_data_j);
232
233 const Number lambda_max =
234 compute_lambda(riemann_data_i, riemann_data_j, h_star);
235
236 return lambda_max;
237 }
238
239
240 template <int dim, typename Number>
242 const state_type &U_i,
243 const state_type &U_j,
244 const unsigned int /*i*/,
245 const unsigned int * /*js*/,
246 const dealii::Tensor<1, dim, Number> &n_ij) const
247 {
248 const auto riemann_data_i = riemann_data_from_state(U_i, n_ij);
249 const auto riemann_data_j = riemann_data_from_state(U_j, n_ij);
250 return compute(riemann_data_i, riemann_data_j);
251 }
252
253 } // namespace ShallowWater
254} // namespace ryujin
Number phi(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j, const Number &h) const
primitive_type riemann_data_from_state(const state_type &U, const dealii::Tensor< 1, dim, Number > &n_ij) const
typename std::array< Number, riemann_data_size > primitive_type
typename View::ScalarNumber ScalarNumber
Number lambda3_plus(const primitive_type &riemann_data, const Number h_star) const
Number lambda1_minus(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
Number compute_h_star(const primitive_type &riemann_data_i, const primitive_type &riemann_data_j) const
typename View::state_type state_type
Number compute(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
DEAL_II_ALWAYS_INLINE Number negative_part(const Number number)
Definition: simd.h:124
DEAL_II_ALWAYS_INLINE Number positive_part(const Number number)
Definition: simd.h:112