8#include <compile_time_options.h>
14 namespace EulerInitialStates
46 template <
typename Description,
int dim,
typename Number>
56 const std::string &subsection)
59 , hyperbolic_system_(hyperbolic_system)
62 if constexpr (!View::have_gamma) {
63 this->add_parameter(
"gamma", gamma_,
"The ratio of specific heats");
67 this->add_parameter(
"velocity galilean frame",
69 "Velocity used to apply a Galilean transformation "
70 "to the otherwise stationary solution");
74 "velocity left", velocity_left_,
"Left limit velocity");
76 velocity_right_ = 7. / 27.;
78 "velocity right", velocity_right_,
"Right limit velocity");
82 "density left", density_left_,
"Left limit density");
85 this->add_parameter(
"mu", mu_,
"Shear viscosity");
89 dealii::ParameterAcceptor::parse_parameters_call_back.connect([
this]() {
90 const auto view = hyperbolic_system_.template view<dim, Number>();
92 if constexpr (View::have_gamma) {
93 gamma_ = view.gamma();
97 velocity_left_ > velocity_right_,
98 dealii::ExcMessage(
"The left limiting velocity must be greater "
99 "than the right limiting velocity"));
101 AssertThrow(velocity_left_ > 0.,
103 "The left limiting velocity must be positive"));
109 const double velocity_origin =
110 std::sqrt(velocity_left_ * velocity_right_);
114 const double Pr = 0.75;
115 const double factor = 2. * gamma_ / (gamma_ + 1.)
116 * mu_ / (density_left_ * velocity_left_ * Pr);
118 psi = [=,
this](
double x,
double v) {
120 velocity_left_ / (velocity_left_ - velocity_right_);
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_);
128 const double value = factor * (c_l * log_l - c_r * log_r) - x;
130 const double derivative = factor * (-c_l / (velocity_left_ - v) -
131 c_r / (v - velocity_right_));
133 return std::make_tuple(value, derivative);
138 constexpr double tol = 1.e-12;
140 const double x_left = std::get<0>(
141 psi(0., (1. - tol) * velocity_left_ + tol * velocity_right_));
143 const double x_right = std::get<0>(
144 psi(0., tol * velocity_left_ + (1. - tol) * velocity_right_));
146 const double norm = (x_right - x_left) * tol;
150 find_velocity = [=,
this](
double x) {
153 return double(velocity_left_);
155 return double(velocity_right_);
159 0.5 * std::tanh(10. * (x - 0.5 * (x_right + x_left)) /
162 velocity_left_ * (0.5 - nu) + velocity_right_ * (nu + 0.5);
164 auto [f, df] = psi(x, v);
166 while (std::abs(f) > norm) {
167 const double v_next = v - f / df;
170 if (std::abs(v_next - v) <
171 tol * 0.5 * (velocity_right_ + velocity_left_)) {
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);
183 const auto [new_f, new_df] = psi(x, v);
195 const auto view = hyperbolic_system_.template view<dim, Number>();
198 const double R_infty = (gamma_ + 1) / (gamma_ - 1);
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());
211 using state_type_1d =
typename Description::
213 const auto state_1d = state_type_1d{
215 Number(rho * (velocity_ + v)),
216 Number(rho * (e + 0.5 * (velocity_ + v) * (velocity_ + v)))}};
218 return view.expand_state(state_1d);
226 Number velocity_left_;
227 Number velocity_right_;
228 Number density_left_;
230 std::function<std::tuple<double, double>(
double,
double)> psi;
231 std::function<double(
double)> find_velocity;
typename Description::template HyperbolicSystemView< dim, Number > View
state_type compute(const dealii::Point< dim > &point, Number t) final
typename View::state_type state_type
BeckerSolution(const HyperbolicSystem &hyperbolic_system, const std::string &subsection)
typename Description::HyperbolicSystem HyperbolicSystem
dealii::Tensor< 1, problem_dimension, Number > state_type
Euler::HyperbolicSystem HyperbolicSystem