9#include <compile_time_options.h>
13#include <deal.II/base/tensor.h>
21 namespace EulerInitialStates
35 template <
typename Description,
int dim,
typename Number>
50 const std::string subsection)
53 , hyperbolic_system_(hyperbolic_system)
56 if constexpr (!View::have_gamma) {
57 this->add_parameter(
"gamma", gamma_,
"The ratio of specific heats");
60 primitive_left_[0] = 1.4;
61 primitive_left_[1] = 0.0;
62 primitive_left_[2] = 1.0;
64 "primitive state left",
66 "Initial 1d primitive state (rho, u, p) on the left");
68 primitive_right_[0] = 1.4;
69 primitive_right_[1] = 0.0;
70 primitive_right_[2] = 1.0;
72 "primitive state right",
74 "Initial 1d primitive state (rho, u, p) on the right");
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();
83 const Number p_L = primitive_left_[2];
84 const Number p_R = primitive_right_[2];
86 p_star_ = compute_pstar(p_L, p_R, primitive_left_, primitive_right_);
88 const Number u_L = primitive_left_[1];
89 u_star_ = u_L - fZofP(p_star_, primitive_left_);
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;
101 lambda_left_minus_ = lambda(p_star_, primitive_left_, -1.);
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.);
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_
118 this->parse_parameters_call_back.connect(prepare_riemann_data);
119 prepare_riemann_data();
125 const auto view = hyperbolic_system_.template view<dim, Number>();
127 const double &x = point[0];
129 const Number xi = x / t;
131 dealii::Tensor<1, 3, Number> primitive_state;
133 if (t < 1.e-14 && x < 0.) {
134 primitive_state = primitive_left_;
136 std::cout <<
"Left primitive state: " << primitive_state << std::endl;
139 }
else if (t < 1.e-14 && x > 0.) {
140 primitive_state = primitive_right_;
142 std::cout <<
"Right primitive state: " << primitive_state
146 }
else if (xi < lambda_left_minus_) {
148 primitive_state = primitive_left_;
150 std::cout <<
"Left primitive state: " << primitive_state << std::endl;
153 }
else if (xi < lambda_left_plus_) {
155 expansion_solution(p_star_, xi, primitive_left_, -1.);
156 primitive_state = c_LL;
158 std::cout <<
"Left expansion state: " << primitive_state << std::endl;
161 }
else if (xi < u_star_) {
162 primitive_state = cstar_solution(p_star_, u_star_, primitive_left_);
164 const Number p_L = primitive_left_[2];
166 primitive_state = expansion_solution(
167 p_star_, lambda_left_plus_, primitive_left_, -1.);
169 std::cout <<
"Left cstar state: " << primitive_state << std::endl;
172 }
else if (xi < lambda_right_minus_) {
173 primitive_state = cstar_solution(p_star_, u_star_, primitive_right_);
175 const Number p_R = primitive_right_[2];
177 primitive_state = expansion_solution(
178 p_star_, lambda_right_minus_, primitive_right_, 1.);
180 std::cout <<
"Right cstar state: " << primitive_state << std::endl;
183 }
else if (xi < lambda_right_plus_) {
185 expansion_solution(p_star_, xi, primitive_right_, 1.);
187 std::cout <<
"Right expansion state: " << primitive_state
193 primitive_state = primitive_right_;
195 std::cout <<
"Right primitive state: " << primitive_state
200 return view.from_initial_state(primitive_state);
212 dealii::Tensor<1, 3, Number> primitive_left_;
213 dealii::Tensor<1, 3, Number> primitive_right_;
225 Number lambda_left_minus_;
226 Number lambda_left_plus_;
227 Number lambda_right_minus_;
228 Number lambda_right_plus_;
236 Number fZofP(
const Number &p_in,
237 const dealii::Tensor<1, 3, Number> &data_in)
const
240 const Number rho_Z = data_in[0];
241 const Number p_Z = data_in[2];
243 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
245 const Number A_Z = 2. / (gamma_ + 1.) / rho_Z;
246 const Number B_Z = (gamma_ - 1.) / (gamma_ + 1.) * p_Z;
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.);
252 Number f_of_p = (p_in - p_Z) * std::sqrt(A_Z / (p_in + B_Z));
261 Number dfZofP(
const Number &p_in,
262 const dealii::Tensor<1, 3, Number> &data_in)
const
265 const Number rho_Z = data_in[0];
266 const Number p_Z = data_in[2];
268 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
270 const Number A_Z = 2. / (gamma_ + 1.) / rho_Z;
271 const Number B_Z = (gamma_ - 1.) / (gamma_ + 1.) * p_Z;
273 Number exp = 0.5 * (gamma_ - 1.) / gamma_;
274 Number left_brach = 2. * c_Z / (gamma_ - 1.) * exp;
277 left_brach *=
std::pow(p_in / p_Z, exp - 1.) / p_Z;
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);
282 Number df_of_p = right_branch;
285 df_of_p = left_brach;
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
295 return dfZofP(p_in, data_left) + dfZofP(p_in, data_right);
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
303 const Number u_L = data_left[1];
304 const Number u_R = data_right[1];
306 return fZofP(p_in, data_right) + fZofP(p_in, data_left) + u_R - u_L;
310 Number lambda(
const Number &p_in,
311 const dealii::Tensor<1, 3, Number> &data_in,
312 const Number &sign)
const
315 const Number rho_Z = data_in[0];
316 const Number u_Z = data_in[1];
317 const Number p_Z = data_in[2];
319 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
321 const Number radicand =
322 1. + 0.5 * (gamma_ + 1.) / gamma_ * std::max(p_in / p_Z - 1., 0.);
324 return u_Z + sign * c_Z * std::sqrt(radicand);
328 Number lambda_intermediate(
const Number &p_in,
329 const dealii::Tensor<1, 3, Number> &data_in,
330 const Number &sign)
const
332 const Number rho_Z = data_in[0];
333 const Number u_Z = data_in[1];
334 const Number p_Z = data_in[2];
336 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
338 const auto lambda_value = lambda(p_in, data_in, sign);
340 const Number f_of_p = fZofP(p_in, data_in);
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));
346 Number result = lambda_value;
348 result = expansion_speed;
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
359 const Number rho_Z = data_in[0];
360 const Number p_Z = data_in[2];
363 const Number p_ratio = p_star / p_Z;
364 const Number gamma_ratio = (gamma_ - 1.) / (gamma_ + 1.);
366 const Number numerator = rho_Z * (p_ratio + gamma_ratio);
367 const Number denominator = gamma_ratio * p_ratio + 1.;
369 Number rho_star = numerator / denominator;
371 auto result = data_in;
372 result[0] = rho_star;
380 dealii::Tensor<1, 3, Number>
381 expansion_solution(
const Number & ,
383 const dealii::Tensor<1, 3, Number> &data_in,
384 const Number &sign)
const
386 const Number rho_Z = data_in[0];
387 const Number u_Z = data_in[1];
388 const Number p_Z = data_in[2];
390 const Number c_Z = std::sqrt(gamma_ * p_Z / rho_Z);
393 const Number gamma_ratio = (gamma_ - 1.) / (gamma_ + 1.);
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.);
399 Number rho_expansion = rho_Z *
std::pow(first - sign * second, exp);
402 const Number factor = p_Z /
std::pow(rho_Z, gamma_);
403 const Number p_expansion = factor *
std::pow(rho_expansion, gamma_);
406 const Number u_expansion = u_Z + sign * fZofP(p_expansion, data_in);
408 auto result = data_in;
409 result[0] = rho_expansion;
410 result[1] = u_expansion;
411 result[2] = p_expansion;
420 double compute_pstar(
double p_1,
422 dealii::Tensor<1, 3, Number> data_1,
423 dealii::Tensor<1, 3, Number> data_2)
425 constexpr Number eps = std::numeric_limits<Number>::epsilon();
431 std::swap(data_1, data_2);
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.,
440 "Euler::ExactRiemannSolver: failed to compute p_star."));
451 std::cout <<
"Computing p_star with a bisection method." << std::endl;
454 unsigned int iter = 0;
455 for (; iter < 200; ++iter) {
458 if (std::abs(p_2 - p_1) < 10. * eps * std::max(p_1, p_2)) {
462 const double phi_2 = phi(p_2, data_1, data_2);
465 const double phi_1 = phi(p_1, data_1, data_2);
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";
474 const auto p_m = 0.5 * (p_2 + p_1);
475 const double phi_m = phi(p_m, data_1, data_2);
477 if (phi_m * phi_2 >= 0.) {
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;
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
typename Description::HyperbolicSystem HyperbolicSystem
ExactRiemannSolution(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
typename View::ScalarNumber ScalarNumber
state_type compute(const dealii::Point< dim > &point, Number t) final
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem