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