ryujin 2.1.1 revision 0b194b984a74af675d09b5e928529ca8c7b634f2
initial_state_noh.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//
5
6#pragma once
7
9
10namespace ryujin
11{
12 namespace EulerInitialStates
13 {
25 template <typename Description, int dim, typename Number>
26 class Noh : public InitialState<Description, dim, Number>
27 {
28 public:
30 using View =
31 typename Description::template HyperbolicSystemView<dim, Number>;
32 using state_type = typename View::state_type;
33
34 Noh(const HyperbolicSystem &hyperbolic_system,
35 const std::string &subsection)
36 : InitialState<Description, dim, Number>("noh", subsection)
37 , hyperbolic_system_(hyperbolic_system)
38 {
39 gamma_ = 1.4;
40 if constexpr (!View::have_gamma) {
41 this->add_parameter("gamma", gamma_, "The ratio of specific heats");
42 }
43
44 rho0_ = 1.0;
45 this->add_parameter(
46 "reference density", rho0_, "The reference density");
47
48 /*
49 * Exact solution assumes this value is negative, but we are just
50 * switching u0 to -u0 by hand in the formulas.
51 */
52 u0_ = 1.0;
53 this->add_parameter("reference velocity magnitude",
54 u0_,
55 "The reference velocity magnitude");
56
57
58 p0_ = 1.e-12;
59 this->add_parameter(
60 "reference pressure", p0_, "The reference pressure");
61
62 this->parse_parameters_call_back.connect([&]() {
63 if constexpr (View::have_gamma) {
64 const auto view = hyperbolic_system_.template view<dim, Number>();
65 gamma_ = view.gamma();
66 }
67 });
68 }
69
70 /* Compute solution */
71 auto compute(const dealii::Point<dim> &point, Number t)
72 -> state_type final
73 {
74 const auto view = hyperbolic_system_.template view<dim, Number>();
75
76 const auto norm = point.norm();
77 const auto min = 10. * std::numeric_limits<Number>::min();
78
79 /* Initialize primitive variables */
80 Number rho = rho0_;
81 auto vel = -u0_ * point / (norm + min);
82 Number p = p0_;
83
84 /* Define exact solutions */
85 const auto D = u0_ * (gamma_ - 1.) / 2.;
86 const bool in_interior = t == Number(0.) ? false : norm / t < D;
87
88 if (in_interior) {
89 rho = rho0_ * std::pow((gamma_ + 1.) / (gamma_ - 1.), dim);
90 vel = 0. * point;
91 p = 0.5 * rho0_ * u0_ * u0_;
92 p *= std::pow(gamma_ + 1., dim) / std::pow(gamma_ - 1., dim - 1);
93 } else {
94 rho = rho0_ * std::pow(1. + t / (norm + min), dim - 1);
95 }
96
97 /* Set final state */
98 if constexpr (dim == 1)
99 return view.from_initial_state(state_type{{rho, vel[0], p}});
100 else if constexpr (dim == 2)
101 return view.from_initial_state(state_type{{rho, vel[0], vel[1], p}});
102 else
103 return view.from_initial_state(
104 state_type{{rho, vel[0], vel[1], vel[2], p}});
105 }
106
107 private:
108 const HyperbolicSystem &hyperbolic_system_;
109 Number gamma_;
110 Number rho0_;
111 Number u0_;
112 Number p0_;
113 };
114 } // namespace EulerInitialStates
115} // namespace ryujin
Noh(const HyperbolicSystem &hyperbolic_system, const std::string &subsection)
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
auto compute(const dealii::Point< dim > &point, Number t) -> state_type final
typename Description::HyperbolicSystem HyperbolicSystem
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32