10#include <compile_time_options.h>
16 namespace ShallowWaterInitialStates
23 template <
typename Description,
int dim,
typename Number>
33 const std::string subsection)
36 , hyperbolic_system_(hyperbolic_system)
38 dealii::ParameterAcceptor::parse_parameters_call_back.connect(
43 this->add_parameter(
"flow state left",
45 "Initial 1d flow state (h, q) on the left");
48 state_right_[1] = 0.0;
49 this->add_parameter(
"flow state right",
51 "Initial 1d flow state (h, q) on the right");
54 this->add_parameter(
"experimental configuration",
56 "Either 'G1', 'G2', 'G3' or 'none' for bathymetry "
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'"));
70 const auto view = hyperbolic_system_.template view<dim, Number>();
72 if constexpr (dim == 1) {
73 AssertThrow(
false, dealii::ExcNotImplemented());
77 const auto temp = point[0] > 1.e-8 ? state_right_ : state_left_;
78 return view.expand_state(temp);
83 initial_precomputed_type
final
86 return {compute_bathymetry(point)};
93 inline Number compute_bathymetry(
const dealii::Point<dim> &point)
const
95 const auto &x = point[0], &y = point[1];
99 if (x >= 0. && x <= 326. / 100.)
101 else if (x > 326. / 100.)
102 bath = -0.0404 * (x - 326. / 100.) - 0.00092 * 326. / 100.;
104 if (which_case_ ==
"none")
108 Number obstacle = 0.;
112 if (which_case_ ==
"G1") {
114 Number obstacle_length = 16.3 / 100.;
115 Number obstacle_width = 8. / 100.;
117 Number xc = 205. / 100. + (16.3 / 2. / 100.);
119 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
120 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
122 obstacle = 7. / 100.;
123 }
else if (which_case_ ==
"G2") {
126 double xc = 184.5 / 100. + 31. / 2. / 100.;
127 const double radicand =
130 const double semi_circle = 7.3 / 100. * std::sqrt(radicand);
132 obstacle = std::max(semi_circle, 0.);
135 double obstacle_length = 16.3 / 100.;
136 double obstacle_width = 8. / 100.;
138 xc = 235. / 100. + (16.3 / 2. / 100.);
140 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
141 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
143 obstacle = 7. / 100.;
144 }
else if (which_case_ ==
"G3") {
147 double xc = 194 / 100. + 31. / 2. / 100.;
148 const double radicand =
150 const double semi_circle = 7.3 / 100. * std::sqrt(radicand);
152 if (y < semi_circle - 24. / 2. / 100. &&
153 std::abs(x - xc) <= 31. / 2. / 100.)
154 obstacle = 21. / 100.;
156 if (y > -semi_circle + 24. / 2. / 100. &&
157 std::abs(x - xc) <= 31. / 2. / 100.)
158 obstacle = 21. / 100.;
161 double obstacle_length = 16.3 / 100.;
162 double obstacle_width = 8. / 100.;
164 xc = 235. / 100. + (16.3 / 2. / 100.);
166 if (std::abs((x - xc) / obstacle_length + y / obstacle_width) +
167 std::abs((x - xc) / obstacle_length - y / obstacle_width) <=
169 obstacle = 7. / 100.;
172 return bath + obstacle;
175 dealii::Tensor<1, 2, Number> state_left_;
176 dealii::Tensor<1, 2, Number> state_right_;
178 std::string which_case_;
179 std::string flow_type_;
typename Description::template HyperbolicSystemView< dim, Number > View
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
typename View::state_type state_type
TankExperiments(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
state_type compute(const dealii::Point< dim > &point, Number) final
typename Description::HyperbolicSystem HyperbolicSystem
void parse_parameters_callback()
T pow(const T x, const T b)
DEAL_II_ALWAYS_INLINE Number positive_part(const Number number)
Euler::HyperbolicSystem HyperbolicSystem