ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
initial_state_becker_solution.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2022 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
11
12namespace ryujin
13{
14 namespace EulerInitialStates
15 {
46 template <typename Description, int dim, typename Number>
47 class BeckerSolution : public InitialState<Description, dim, Number>
48 {
49 public:
51 using View =
52 typename Description::template HyperbolicSystemView<dim, Number>;
53 using state_type = typename View::state_type;
54
55 BeckerSolution(const HyperbolicSystem &hyperbolic_system,
56 const std::string &subsection)
57 : InitialState<Description, dim, Number>("becker solution",
58 subsection)
59 , hyperbolic_system_(hyperbolic_system)
60 {
61 gamma_ = 1.4;
62 if constexpr (!View::have_gamma) {
63 this->add_parameter("gamma", gamma_, "The ratio of specific heats");
64 }
65
66 velocity_ = 0.2;
67 this->add_parameter("velocity galilean frame",
68 velocity_,
69 "Velocity used to apply a Galilean transformation "
70 "to the otherwise stationary solution");
71
72 velocity_left_ = 1.0;
73 this->add_parameter(
74 "velocity left", velocity_left_, "Left limit velocity");
75
76 velocity_right_ = 7. / 27.;
77 this->add_parameter(
78 "velocity right", velocity_right_, "Right limit velocity");
79
80 density_left_ = 1.0;
81 this->add_parameter(
82 "density left", density_left_, "Left limit density");
83
84 mu_ = 0.01;
85 this->add_parameter("mu", mu_, "Shear viscosity");
86
87 /* Callback: */
88
89 dealii::ParameterAcceptor::parse_parameters_call_back.connect([this]() {
90 const auto view = hyperbolic_system_.template view<dim, Number>();
91
92 if constexpr (View::have_gamma) {
93 gamma_ = view.gamma();
94 }
95
96 AssertThrow(
97 velocity_left_ > velocity_right_,
98 dealii::ExcMessage("The left limiting velocity must be greater "
99 "than the right limiting velocity"));
100
101 AssertThrow(velocity_left_ > 0.,
102 dealii::ExcMessage(
103 "The left limiting velocity must be positive"));
104
105 /*
106 * Set up all helper functions and quantities:
107 */
108
109 const double velocity_origin =
110 std::sqrt(velocity_left_ * velocity_right_);
111
112 /* Prefactor as given in: (7.1) */
113
114 const double Pr = 0.75;
115 const double factor = 2. * gamma_ / (gamma_ + 1.) //
116 * mu_ / (density_left_ * velocity_left_ * Pr);
117
118 psi = [=, this](double x, double v) {
119 const double c_l =
120 velocity_left_ / (velocity_left_ - velocity_right_);
121 const double c_r =
122 velocity_right_ / (velocity_left_ - velocity_right_);
123 const double log_l = std::log(velocity_left_ - v) -
124 std::log(velocity_left_ - velocity_origin);
125 const double log_r = std::log(v - velocity_right_) -
126 std::log(velocity_origin - velocity_right_);
127
128 const double value = factor * (c_l * log_l - c_r * log_r) - x;
129
130 const double derivative = factor * (-c_l / (velocity_left_ - v) -
131 c_r / (v - velocity_right_));
132
133 return std::make_tuple(value, derivative);
134 };
135
136 /* Determine cut-off points: */
137
138 constexpr double tol = 1.e-12;
139
140 const double x_left = std::get<0>(
141 psi(0., (1. - tol) * velocity_left_ + tol * velocity_right_));
142
143 const double x_right = std::get<0>(
144 psi(0., tol * velocity_left_ + (1. - tol) * velocity_right_));
145
146 const double norm = (x_right - x_left) * tol;
147
148 /* Root finding algorithm: */
149
150 find_velocity = [=, this](double x) {
151 /* Return extremal cases: */
152 if (x <= x_left)
153 return double(velocity_left_);
154 if (x >= x_right)
155 return double(velocity_right_);
156
157 /* Interpolate initial guess: */
158 const auto nu =
159 0.5 * std::tanh(10. * (x - 0.5 * (x_right + x_left)) /
160 (x_right - x_left));
161 double v =
162 velocity_left_ * (0.5 - nu) + velocity_right_ * (nu + 0.5);
163
164 auto [f, df] = psi(x, v);
165
166 while (std::abs(f) > norm) {
167 const double v_next = v - f / df;
168
169 /* Also break if we made no progress: */
170 if (std::abs(v_next - v) <
171 tol * 0.5 * (velocity_right_ + velocity_left_)) {
172 v = v_next;
173 break;
174 }
175
176 if (v_next < velocity_right_)
177 v = 0.5 * (velocity_right_ + v);
178 else if (v_next > velocity_left_)
179 v = 0.5 * (velocity_left_ + v);
180 else
181 v = v_next;
182
183 const auto [new_f, new_df] = psi(x, v);
184 f = new_f;
185 df = new_df;
186 }
187
188 return v;
189 }; /* find_velocity */
190 });
191 }
192
193 state_type compute(const dealii::Point<dim> &point, Number t) final
194 {
195 const auto view = hyperbolic_system_.template view<dim, Number>();
196
197 /* (7.2) */
198 const double R_infty = (gamma_ + 1) / (gamma_ - 1);
199
200 /* (7.3) */
201 const double x = point[0] - velocity_ * t;
202 const double v = find_velocity(x);
203 Assert(v >= velocity_right_, dealii::ExcInternalError());
204 Assert(v <= velocity_left_, dealii::ExcInternalError());
205 const double rho = density_left_ * velocity_left_ / v;
206 Assert(rho > 0., dealii::ExcInternalError());
207 const double e = 1. / (2. * gamma_) *
208 (R_infty * velocity_left_ * velocity_right_ - v * v);
209 Assert(e > 0., dealii::ExcInternalError());
210
211 using state_type_1d = typename Description::
213 const auto state_1d = state_type_1d{
214 {Number(rho),
215 Number(rho * (velocity_ + v)),
216 Number(rho * (e + 0.5 * (velocity_ + v) * (velocity_ + v)))}};
217
218 return view.expand_state(state_1d);
219 }
220
221 private:
222 const HyperbolicSystem &hyperbolic_system_;
223 Number gamma_;
224
225 Number velocity_;
226 Number velocity_left_;
227 Number velocity_right_;
228 Number density_left_;
229 Number mu_;
230 std::function<std::tuple<double, double>(double, double)> psi;
231 std::function<double(double)> find_velocity;
232 };
233
234 } // namespace EulerInitialStates
235} // namespace ryujin
typename Description::template HyperbolicSystemView< dim, Number > View
state_type compute(const dealii::Point< dim > &point, Number t) final
BeckerSolution(const HyperbolicSystem &hyperbolic_system, const std::string &subsection)
typename Description::HyperbolicSystem HyperbolicSystem
dealii::Tensor< 1, problem_dimension, Number > state_type
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:34