ryujin 2.1.1 revision ef0fcd4010d109b860652ece4a7b8963fb7d46b1
initial_state_exact_riemann_solution.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 - 2024 by the ryujin authors
4// Copyright (C) 2025 by Triad National Security, LLC
5//
6
7#pragma once
8
9#include <compile_time_options.h>
10
12
13#include <deal.II/base/tensor.h>
14
15#include <cmath>
16
17// #define DEBUG_SOLUTION
18
19namespace ryujin
20{
21 namespace EulerInitialStates
22 {
35 template <typename Description, int dim, typename Number>
36 class ExactRiemannSolution : public InitialState<Description, dim, Number>
37 {
38 public:
40
42 using View =
43 typename Description::template HyperbolicSystemView<dim, Number>;
44 using state_type = typename View::state_type;
45
46 using ScalarNumber = typename View::ScalarNumber;
47
48
49 ExactRiemannSolution(const HyperbolicSystem &hyperbolic_system,
50 const std::string subsection)
51 : InitialState<Description, dim, Number>("exact riemann solution",
52 subsection)
53 , hyperbolic_system_(hyperbolic_system)
54 {
55 gamma_ = 1.4;
56 if constexpr (!View::have_gamma) {
57 this->add_parameter("gamma", gamma_, "The ratio of specific heats");
58 }
59
60 primitive_left_[0] = 1.4;
61 primitive_left_[1] = 0.0;
62 primitive_left_[2] = 1.0;
63 this->add_parameter(
64 "primitive state left",
65 primitive_left_,
66 "Initial 1d primitive state (rho, u, p) on the left");
67
68 primitive_right_[0] = 1.4;
69 primitive_right_[1] = 0.0;
70 primitive_right_[2] = 1.0;
71 this->add_parameter(
72 "primitive state right",
73 primitive_right_,
74 "Initial 1d primitive state (rho, u, p) on the right");
75
76 // Convert the primitive states to conserved states
77 const auto prepare_riemann_data = [&]() {
78 const auto view = hyperbolic_system_.template view<dim, Number>();
79 if constexpr (View::have_gamma) {
80 gamma_ = view.gamma();
81 }
82
83 const Number p_L = primitive_left_[2];
84 const Number p_R = primitive_right_[2];
85
86 p_star_ = compute_pstar(p_L, p_R, primitive_left_, primitive_right_);
87
88 const Number u_L = primitive_left_[1];
89 u_star_ = u_L - fZofP(p_star_, primitive_left_);
90
91#ifdef DEBUG_SOLUTION
92 const Number u_R = primitive_right_[1];
93 std::cout << "left data = " << primitive_left_
94 << "\nright data = " << primitive_right_
95 << "\np_star = " << p_star_
96 << "\nu_star = " << u_star_
97 << "\nVerifying u_star = "
98 << u_R + fZofP(p_star_, primitive_right_) << std::endl;
99#endif
100
101 lambda_left_minus_ = lambda(p_star_, primitive_left_, -1.);
102 lambda_left_plus_ =
103 lambda_intermediate(p_star_, primitive_left_, -1.);
104 lambda_right_minus_ =
105 lambda_intermediate(p_star_, primitive_right_, 1.);
106 lambda_right_plus_ = lambda(p_star_, primitive_right_, 1.);
107
108
109#ifdef DEBUG_SOLUTION
110 std::cout << "lambda_left_minus = " << lambda_left_minus_
111 << "\nlambda_left_plus = " << lambda_left_plus_
112 << "\nlambda_right_minus = " << lambda_right_minus_
113 << "\nlambda_right_plus = " << lambda_right_plus_
114 << std::endl;
115#endif
116 };
117
118 this->parse_parameters_call_back.connect(prepare_riemann_data);
119 prepare_riemann_data();
120 }
121
122
123 state_type compute(const dealii::Point<dim> &point, Number t) final
124 {
125 const auto view = hyperbolic_system_.template view<dim, Number>();
126
127 const double &x = point[0];
128
129 const Number xi = x / t;
130
131 dealii::Tensor<1, 3, Number> primitive_state;
132
133 if (t < 1.e-14 && x < 0.) {
134 primitive_state = primitive_left_;
135#ifdef DEBUG_SOLUTION
136 std::cout << "Left primitive state: " << primitive_state << std::endl;
137#endif
138
139 } else if (t < 1.e-14 && x > 0.) {
140 primitive_state = primitive_right_;
141#ifdef DEBUG_SOLUTION
142 std::cout << "Right primitive state: " << primitive_state
143 << std::endl;
144#endif
145
146 } else if (xi < lambda_left_minus_) {
147 /* Left state: */
148 primitive_state = primitive_left_;
149#ifdef DEBUG_SOLUTION
150 std::cout << "Left primitive state: " << primitive_state << std::endl;
151#endif
152
153 } else if (xi < lambda_left_plus_) {
154 const auto c_LL =
155 expansion_solution(p_star_, xi, primitive_left_, -1.);
156 primitive_state = c_LL;
157#ifdef DEBUG_SOLUTION
158 std::cout << "Left expansion state: " << primitive_state << std::endl;
159#endif
160
161 } else if (xi < u_star_) {
162 primitive_state = cstar_solution(p_star_, u_star_, primitive_left_);
163
164 const Number p_L = primitive_left_[2];
165 if (p_star_ < p_L)
166 primitive_state = expansion_solution(
167 p_star_, lambda_left_plus_, primitive_left_, -1.);
168#ifdef DEBUG_SOLUTION
169 std::cout << "Left cstar state: " << primitive_state << std::endl;
170#endif
171
172 } else if (xi < lambda_right_minus_) {
173 primitive_state = cstar_solution(p_star_, u_star_, primitive_right_);
174
175 const Number p_R = primitive_right_[2];
176 if (p_star_ < p_R)
177 primitive_state = expansion_solution(
178 p_star_, lambda_right_minus_, primitive_right_, 1.);
179#ifdef DEBUG_SOLUTION
180 std::cout << "Right cstar state: " << primitive_state << std::endl;
181#endif
182
183 } else if (xi < lambda_right_plus_) {
184 primitive_state =
185 expansion_solution(p_star_, xi, primitive_right_, 1.);
186#ifdef DEBUG_SOLUTION
187 std::cout << "Right expansion state: " << primitive_state
188 << std::endl;
189#endif
190
191 } else {
192 /* Right state: */
193 primitive_state = primitive_right_;
194#ifdef DEBUG_SOLUTION
195 std::cout << "Right primitive state: " << primitive_state
196 << std::endl;
197#endif
198 }
199
200 return view.from_initial_state(primitive_state);
201 }
202
203 private:
205
209
210 Number gamma_;
211
212 dealii::Tensor<1, 3, Number> primitive_left_;
213 dealii::Tensor<1, 3, Number> primitive_right_;
214
216
220
221 const HyperbolicSystem &hyperbolic_system_;
222
223 Number p_star_;
224 Number u_star_;
225 Number lambda_left_minus_;
226 Number lambda_left_plus_;
227 Number lambda_right_minus_;
228 Number lambda_right_plus_;
229
231
235
236 Number fZofP(const Number &p_in,
237 const dealii::Tensor<1, 3, Number> &data_in) const
238 {
239 // Get left/right data
240 const Number rho_Z = data_in[0];
241 const Number p_Z = data_in[2];
242
243 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
244
245 const Number A_Z = 2. / (gamma_ + 1.) / rho_Z;
246 const Number B_Z = (gamma_ - 1.) / (gamma_ + 1.) * p_Z;
247
248 const Number exp = 0.5 * (gamma_ - 1.) / gamma_;
249 Number left_brach = 2. * c_Z / (gamma_ - 1.);
250 left_brach *= (std::pow(p_in / p_Z, exp) - 1.);
251
252 Number f_of_p = (p_in - p_Z) * std::sqrt(A_Z / (p_in + B_Z));
253
254 if (p_in <= p_Z)
255 f_of_p = left_brach;
256
257 return f_of_p;
258 }
259
260
261 Number dfZofP(const Number &p_in,
262 const dealii::Tensor<1, 3, Number> &data_in) const
263 {
264 // Get left/right data
265 const Number rho_Z = data_in[0];
266 const Number p_Z = data_in[2];
267
268 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
269
270 const Number A_Z = 2. / (gamma_ + 1.) / rho_Z;
271 const Number B_Z = (gamma_ - 1.) / (gamma_ + 1.) * p_Z;
272
273 Number exp = 0.5 * (gamma_ - 1.) / gamma_;
274 Number left_brach = 2. * c_Z / (gamma_ - 1.) * exp;
275 exp -= 1.;
276
277 left_brach *= std::pow(p_in / p_Z, exp - 1.) / p_Z;
278
279 Number right_branch = std::pow(A_Z / (p_in + B_Z), 1.5);
280 right_branch *= (2. * B_Z + p_in + p_Z) / (2. * A_Z);
281
282 Number df_of_p = right_branch;
283
284 if (p_in <= p_Z)
285 df_of_p = left_brach;
286
287 return df_of_p;
288 }
289
290
291 Number dphi(const Number &p_in,
292 const dealii::Tensor<1, 3, Number> &data_left,
293 const dealii::Tensor<1, 3, Number> &data_right) const
294 {
295 return dfZofP(p_in, data_left) + dfZofP(p_in, data_right);
296 }
297
298
299 Number phi(const Number &p_in,
300 const dealii::Tensor<1, 3, Number> &data_left,
301 const dealii::Tensor<1, 3, Number> &data_right) const
302 {
303 const Number u_L = data_left[1];
304 const Number u_R = data_right[1];
305
306 return fZofP(p_in, data_right) + fZofP(p_in, data_left) + u_R - u_L;
307 }
308
309
310 Number lambda(const Number &p_in,
311 const dealii::Tensor<1, 3, Number> &data_in,
312 const Number &sign) const
313 {
314 // Get left/right data
315 const Number rho_Z = data_in[0];
316 const Number u_Z = data_in[1];
317 const Number p_Z = data_in[2];
318
319 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
320
321 const Number radicand =
322 1. + 0.5 * (gamma_ + 1.) / gamma_ * std::max(p_in / p_Z - 1., 0.);
323
324 return u_Z + sign * c_Z * std::sqrt(radicand);
325 }
326
327
328 Number lambda_intermediate(const Number &p_in,
329 const dealii::Tensor<1, 3, Number> &data_in,
330 const Number &sign) const
331 {
332 const Number rho_Z = data_in[0];
333 const Number u_Z = data_in[1];
334 const Number p_Z = data_in[2];
335
336 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
337
338 const auto lambda_value = lambda(p_in, data_in, sign);
339
340 const Number f_of_p = fZofP(p_in, data_in);
341
342 const Number exp = 0.5 * (gamma_ - 1.) / gamma_;
343 const Number expansion_speed =
344 u_Z + sign * (f_of_p + c_Z * std::pow(p_in / p_Z, exp));
345
346 Number result = lambda_value;
347 if (p_in < p_Z)
348 result = expansion_speed;
349
350 return result;
351 }
352
353
354 dealii::Tensor<1, 3, Number>
355 cstar_solution(const Number &p_star,
356 const Number &u_star,
357 const dealii::Tensor<1, 3, Number> &data_in) const
358 {
359 const Number rho_Z = data_in[0];
360 const Number p_Z = data_in[2];
361
362 // Define rho_star
363 const Number p_ratio = p_star / p_Z;
364 const Number gamma_ratio = (gamma_ - 1.) / (gamma_ + 1.);
365
366 const Number numerator = rho_Z * (p_ratio + gamma_ratio);
367 const Number denominator = gamma_ratio * p_ratio + 1.;
368
369 Number rho_star = numerator / denominator;
370
371 auto result = data_in;
372 result[0] = rho_star;
373 result[1] = u_star;
374 result[2] = p_star;
375
376 return result;
377 }
378
379
380 dealii::Tensor<1, 3, Number>
381 expansion_solution(const Number & /*p_star*/,
382 const Number &xi,
383 const dealii::Tensor<1, 3, Number> &data_in,
384 const Number &sign) const
385 {
386 const Number rho_Z = data_in[0];
387 const Number u_Z = data_in[1];
388 const Number p_Z = data_in[2];
389
390 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
391
392 // Define rho_expansion
393 const Number gamma_ratio = (gamma_ - 1.) / (gamma_ + 1.);
394
395 const Number first = 2. / (gamma_ + 1.);
396 const Number second = gamma_ratio / c_Z * (u_Z - xi);
397 const Number exp = 2. / (gamma_ - 1.);
398
399 Number rho_expansion = rho_Z * std::pow(first - sign * second, exp);
400
401 // Define p_expansion
402 const Number factor = p_Z / std::pow(rho_Z, gamma_);
403 const Number p_expansion = factor * std::pow(rho_expansion, gamma_);
404
405 // Define u_expansion
406 const Number u_expansion = u_Z + sign * fZofP(p_expansion, data_in);
407
408 auto result = data_in;
409 result[0] = rho_expansion;
410 result[1] = u_expansion;
411 result[2] = p_expansion;
412
413 return result;
414 }
415
416
420 double compute_pstar(double p_1,
421 double p_2,
422 dealii::Tensor<1, 3, Number> data_1,
423 dealii::Tensor<1, 3, Number> data_2)
424 {
425 constexpr Number eps = std::numeric_limits<Number>::epsilon();
426
427 // Ensure that p_1 <= p_2
428
429 if (p_1 > p_2) {
430 std::swap(p_1, p_2);
431 std::swap(data_1, data_2);
432 }
433
434#ifdef DEBUG
435 {
436 const double phi_1 = phi(p_1, data_1, data_2);
437 const double phi_2 = phi(p_2, data_1, data_2);
438 Assert(phi_1 * phi_2 <= 0.,
439 dealii::ExcMessage(
440 "Euler::ExactRiemannSolver: failed to compute p_star."));
441 }
442#endif
443
444 //
445 // We simply compute the root of phi with a bisection method down
446 // to machine precision. This is not terribly efficient but luckily
447 // happens only once during initialization.
448 //
449
450#ifdef DEBUG_SOLUTION
451 std::cout << "Computing p_star with a bisection method." << std::endl;
452#endif
453
454 unsigned int iter = 0;
455 for (; iter < 200; ++iter) {
456
457 // Check for convergence:
458 if (std::abs(p_2 - p_1) < 10. * eps * std::max(p_1, p_2)) {
459 break;
460 }
461
462 const double phi_2 = phi(p_2, data_1, data_2);
463
464#ifdef DEBUG_SOLUTION
465 const double phi_1 = phi(p_1, data_1, data_2);
466
467 std::cout << "\niter: " << iter << "\n";
468 std::cout << "p_1: " << p_1 << "\n";
469 std::cout << "p_2: " << p_2 << "\n";
470 std::cout << "phi_1: " << phi_1 << "\n";
471 std::cout << "phi_2: " << phi_2 << "\n";
472#endif
473
474 const auto p_m = 0.5 * (p_2 + p_1);
475 const double phi_m = phi(p_m, data_1, data_2);
476
477 if (phi_m * phi_2 >= 0.) {
478 p_2 = p_m;
479 } else {
480 p_1 = p_m;
481 }
482 }
483
484#ifdef DEBUG_SOLUTION
485 const double phi_2 = phi(p_2, data_1, data_2);
486 std::cout << "After " << iter << " iterations:"
487 << "\np_star = " << p_2 << "\nphi(p_star) = " << phi_2
488 << "\n|p_2 - p_1| = " << std::abs(p_2 - p_1) << std::endl;
489#endif
490
491 return p_2;
492 }
493
495 };
496 } // namespace EulerInitialStates
497} // namespace ryujin
typename Description::template HyperbolicSystemView< dim, Number > View
ExactRiemannSolution(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
state_type compute(const dealii::Point< dim > &point, Number t) final
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:34