ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
initial_state_smooth_vortex.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0
3// [LANL Copyright Statement]
4// Copyright (C) 2022 - 2024 by the ryujin authors
5// Copyright (C) 2023 - 2024 by Triad National Security, LLC
6//
7
8#pragma once
9
11
12namespace ryujin
13{
14 namespace ShallowWaterInitialStates
15 {
24 template <typename Description, int dim, typename Number>
25 class SmoothVortex : public InitialState<Description, dim, Number>
26 {
27 public:
29 using View =
30 typename Description::template HyperbolicSystemView<dim, Number>;
31 using state_type = typename View::state_type;
32
33 SmoothVortex(const HyperbolicSystem &hyperbolic_system,
34 const std::string subsection)
35 : InitialState<Description, dim, Number>("smooth vortex", subsection)
36 , hyperbolic_system_(hyperbolic_system)
37 {
38 with_bathymetry = false;
39 this->add_parameter("with bathymetry",
40 with_bathymetry,
41 "If set to true then the initial profile includes "
42 "the bathymetry for the steady vortex. ");
43
44 depth_ = 1.0;
45 this->add_parameter("reference depth", depth_, "Reference water depth");
46
47 mach_number_ = 2.0;
48 this->add_parameter(
49 "mach number", mach_number_, "Mach number of unsteady vortex");
50
51 beta_ = 0.1;
52 this->add_parameter("beta", beta_, "vortex strength beta");
53 }
54
55 state_type compute(const dealii::Point<dim> &point, Number t) final
56 {
57 const auto view = hyperbolic_system_.template view<dim, Number>();
58 const auto gravity = view.gravity();
59
60 dealii::Point<2> point_bar;
61 point_bar[0] = point[0] - mach_number_ * t;
62 point_bar[1] = point[1];
63
64 const Number r_square = Number(point_bar.norm_square());
65
66 /* We assume r0 = 1 */
67 const Number factor =
68 beta_ / Number(2. * M_PI) *
69 exp(Number(0.5) - Number(0.5) * r_square /* /r0^ 2*/);
70
71 Number h =
72 depth_ - Number(1. / (2. * gravity /* * r0^2 */)) * factor * factor;
73
74 if (with_bathymetry)
75 h -= compute_bathymetry(point);
76
77 const Number u = mach_number_ - factor * Number(point_bar[1]);
78 const Number v = factor * Number(point_bar[0]);
79
80 if constexpr (dim == 2)
81 return state_type{{h, h * u, h * v}};
82
83 else {
84 AssertThrow(false, dealii::ExcNotImplemented());
85 __builtin_trap();
86 }
87 }
88
89 auto initial_precomputations(const dealii::Point<dim> &point) ->
91 initial_precomputed_type final
92 {
93 /* Compute bathymetry: */
94 return {compute_bathymetry(point)};
95 }
96
97 private:
98 const HyperbolicSystem &hyperbolic_system_;
99
100 DEAL_II_ALWAYS_INLINE inline Number
101 compute_bathymetry(const dealii::Point<dim> &point) const
102 {
103 const Number r_square = Number(point.norm_square());
104
105 /* We assume r0 = 1 */
106 const Number factor =
107 beta_ / Number(2. * M_PI) *
108 exp(Number(0.5) - Number(0.5) * r_square /* /r0^ 2*/);
109
110 Number bath = 0.;
111 if (with_bathymetry)
112 bath = depth_ / 4. * factor;
113
114 return bath;
115 }
116
117 bool with_bathymetry;
118 Number depth_;
119 Number mach_number_;
120 Number beta_;
121 };
122
123 } // namespace ShallowWaterInitialStates
124} // namespace ryujin
typename Description::HyperbolicSystem HyperbolicSystem
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
typename Description::template HyperbolicSystemView< dim, Number > View
state_type compute(const dealii::Point< dim > &point, Number t) final
SmoothVortex(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32