15 namespace ShallowWaterInitialStates
23 template <
typename Description,
int dim,
typename Number>
33 const std::string subsection)
36 , hyperbolic_system_(hyperbolic_system)
38 well_balancing_validation =
false;
40 "well balancing validation",
41 well_balancing_validation,
42 "If set to true then the initial profile is returned for all "
44 "(t>0); otherwise a constant inflow is computed for t>0 suitable "
45 "for prescribing Dirichlet conditions at the inflow boundary.");
49 this->add_parameter(
"left water depth",
51 "Depth of water to the left of pseudo-dam");
53 this->add_parameter(
"right water depth",
55 "Depth of water to the right of pseudo-dam");
58 this->add_parameter(
"cone magnitude",
60 "To modify magnitude of cone heights");
65 const auto view = hyperbolic_system_.template view<dim, Number>();
67 const Number x = point[0];
71 if (t <= 1.e-10 || well_balancing_validation) {
72 Number h = x < 16. ? left_depth : right_depth;
73 h = std::max(h - compute_bathymetry(point), Number(0.));
79 const auto &h = left_depth;
80 const auto a = view.speed_of_sound(
state_type{{h, Number(0.)}});
86 initial_precomputed_type
final
89 return {compute_bathymetry(point)};
95 DEAL_II_ALWAYS_INLINE
inline Number
96 compute_bathymetry(
const dealii::Point<dim> &point)
const
98 if constexpr (dim == 1) {
100 const Number &x = point[0];
102 Number z3 = 3. - 3. / 10. * std::sqrt(
std::pow(x - 47.5, 2));
103 return cone_magnitude * std::max({z3, Number(0.)});
106 const Number &x = point[0];
107 const Number &y = point[1];
121 return cone_magnitude * std::max({z1, z2, z3, Number(0.)});
124 bool well_balancing_validation;
127 Number cone_magnitude;
typename Description::HyperbolicSystem HyperbolicSystem
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
state_type compute(const dealii::Point< dim > &point, Number t) final
ThreeBumpsDamBreak(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem