ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
initial_state_transient.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0
3// [LANL Copyright Statement]
4// Copyright (C) 2024 by the ryujin authors
5// Copyright (C) 2023 - 2024 by Triad National Security, LLC
6//
7
8#pragma once
9
10#include <compile_time_options.h>
11
13
14namespace ryujin
15{
16 namespace ShallowWaterInitialStates
17 {
23 template <typename Description, int dim, typename Number>
24 class TankExperiments : public InitialState<Description, dim, Number>
25 {
26 public:
28 using View =
29 typename Description::template HyperbolicSystemView<dim, Number>;
30 using state_type = typename View::state_type;
31
32 TankExperiments(const HyperbolicSystem &hyperbolic_system,
33 const std::string subsection)
34 : InitialState<Description, dim, Number>("transient experiments",
35 subsection)
36 , hyperbolic_system_(hyperbolic_system)
37 {
38 dealii::ParameterAcceptor::parse_parameters_call_back.connect(
40
41 state_left_[0] = 1.;
42 state_left_[1] = 0.0;
43 this->add_parameter("flow state left",
44 state_left_,
45 "Initial 1d flow state (h, q) on the left");
46
47 state_right_[0] = 1.;
48 state_right_[1] = 0.0;
49 this->add_parameter("flow state right",
50 state_right_,
51 "Initial 1d flow state (h, q) on the right");
52
53 which_case_ = "G1";
54 this->add_parameter("experimental configuration",
55 which_case_,
56 "Either 'G1', 'G2', 'G3' or 'none' for bathymetry "
57 "configuration");
58 }
59
61 {
62 AssertThrow(
63 which_case_ == "G1" || which_case_ == "G2" || which_case_ == "G3" ||
64 which_case_ == "none",
65 dealii::ExcMessage("Case must be 'G1', 'G2', 'G3' or 'none'"));
66 }
67
68 state_type compute(const dealii::Point<dim> &point, Number /* t */) final
69 {
70 const auto view = hyperbolic_system_.template view<dim, Number>();
71
72 if constexpr (dim == 1) {
73 AssertThrow(false, dealii::ExcNotImplemented());
74 __builtin_trap();
75 }
76
77 const auto temp = point[0] > 1.e-8 ? state_right_ : state_left_;
78 return view.expand_state(temp);
79 }
80
81 auto initial_precomputations(const dealii::Point<dim> &point) ->
83 initial_precomputed_type final
84 {
85 /* Compute bathymetry: */
86 return {compute_bathymetry(point)};
87 }
88
89 private:
90 const HyperbolicSystem &hyperbolic_system_;
91
92 DEAL_II_ALWAYS_INLINE
93 inline Number compute_bathymetry(const dealii::Point<dim> &point) const
94 {
95 const auto &x = point[0], &y = point[1];
96
97 /* Bathymetry base is the same for all configurations */
98 Number bath = 0.;
99 if (x >= 0. && x <= 326. / 100.)
100 bath = -0.00092 * x;
101 else if (x > 326. / 100.)
102 bath = -0.0404 * (x - 326. / 100.) - 0.00092 * 326. / 100.;
103
104 if (which_case_ == "none")
105 return bath;
106
107 /* Initialize obstacle to 0 */
108 Number obstacle = 0.;
109
110
111 // G1 -- rectangular obstacle
112 if (which_case_ == "G1") {
113
114 Number obstacle_length = 16.3 / 100.;
115 Number obstacle_width = 8. / 100.;
116
117 Number xc = 205. / 100. + (16.3 / 2. / 100.); // obstacle center
118
119 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
120 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
121 1.)
122 obstacle = 7. / 100.;
123 } else if (which_case_ == "G2") { // circular bump + rectangle
124
125 // circular bump
126 double xc = 184.5 / 100. + 31. / 2. / 100.;
127 const double radicand =
128 positive_part(1. - std::pow((x - xc) / (31. / 2. / 100.), 2));
129
130 const double semi_circle = 7.3 / 100. * std::sqrt(radicand);
131
132 obstacle = std::max(semi_circle, 0.);
133
134 // rectangular obstacle
135 double obstacle_length = 16.3 / 100.;
136 double obstacle_width = 8. / 100.;
137
138 xc = 235. / 100. + (16.3 / 2. / 100.); // obstacle center
139
140 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
141 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
142 1.)
143 obstacle = 7. / 100.;
144 } else if (which_case_ == "G3") { // narrowing half-circle + rectangle
145
146 // narrowing half-circles canal
147 double xc = 194 / 100. + 31. / 2. / 100.;
148 const double radicand =
149 positive_part(1. - std::pow((x - xc) / (31. / 2. / 100.), 2));
150 const double semi_circle = 7.3 / 100. * std::sqrt(radicand);
151
152 if (y < semi_circle - 24. / 2. / 100. &&
153 std::abs(x - xc) <= 31. / 2. / 100.)
154 obstacle = 21. / 100.;
155
156 if (y > -semi_circle + 24. / 2. / 100. &&
157 std::abs(x - xc) <= 31. / 2. / 100.)
158 obstacle = 21. / 100.;
159
160 // rectangular obstacle
161 double obstacle_length = 16.3 / 100.;
162 double obstacle_width = 8. / 100.;
163
164 xc = 235. / 100. + (16.3 / 2. / 100.); // obstacle center
165
166 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
167 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
168 1.)
169 obstacle = 7. / 100.;
170 }
171
172 return bath + obstacle;
173 }
174
175 dealii::Tensor<1, 2, Number> state_left_;
176 dealii::Tensor<1, 2, Number> state_right_;
177
178 std::string which_case_;
179 std::string flow_type_;
180 };
181
182 } // namespace ShallowWaterInitialStates
183} // namespace ryujin
typename Description::template HyperbolicSystemView< dim, Number > View
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
TankExperiments(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
state_type compute(const dealii::Point< dim > &point, Number) final
typename Description::HyperbolicSystem HyperbolicSystem
T pow(const T x, const T b)
DEAL_II_ALWAYS_INLINE Number positive_part(const Number number)
Definition: simd.h:113
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:34