ryujin 2.1.1 revision 336b16a72e829721302c626ec7071b92032b8248
initial_state_three_bumps_dam_break.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
10#include "hyperbolic_system.h"
12
13namespace ryujin
14{
15 namespace ShallowWaterInitialStates
16 {
23 template <typename Description, int dim, typename Number>
24 class ThreeBumpsDamBreak : 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 ThreeBumpsDamBreak(const HyperbolicSystem &hyperbolic_system,
33 const std::string subsection)
34 : InitialState<Description, dim, Number>("three bumps dam break",
35 subsection)
36 , hyperbolic_system_(hyperbolic_system)
37 {
38 well_balancing_validation = false;
39 this->add_parameter(
40 "well balancing validation",
41 well_balancing_validation,
42 "If set to true then the initial profile is returned for all "
43 "times "
44 "(t>0); otherwise a constant inflow is computed for t>0 suitable "
45 "for prescribing Dirichlet conditions at the inflow boundary.");
46
47
48 left_depth = 1.875;
49 this->add_parameter("left water depth",
50 left_depth,
51 "Depth of water to the left of pseudo-dam");
52 right_depth = 0.;
53 this->add_parameter("right water depth",
54 right_depth,
55 "Depth of water to the right of pseudo-dam");
56
57 cone_magnitude = 1.;
58 this->add_parameter("cone magnitude",
59 cone_magnitude,
60 "To modify magnitude of cone heights");
61 }
62
63 state_type compute(const dealii::Point<dim> &point, Number t) final
64 {
65 const auto view = hyperbolic_system_.template view<dim, Number>();
66
67 const Number x = point[0];
68
69 /* Initial state: */
70
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.));
74 return state_type{{h, 0.}};
75 }
76
77 /* For t > 0 prescribe constant inflow Dirichlet data on the left: */
78
79 const auto &h = left_depth;
80 const auto a = view.speed_of_sound(state_type{{h, Number(0.)}});
81 return state_type{{h, h * a}};
82 }
83
84 auto initial_precomputations(const dealii::Point<dim> &point) ->
86 initial_precomputed_type final
87 {
88 /* Compute bathymetry: */
89 return {compute_bathymetry(point)};
90 }
91
92 private:
93 const HyperbolicSystem &hyperbolic_system_;
94
95 DEAL_II_ALWAYS_INLINE inline Number
96 compute_bathymetry(const dealii::Point<dim> &point) const
97 {
98 if constexpr (dim == 1) {
99 /* When dim = 1, we only have one cone */
100 const Number &x = point[0];
101
102 Number z3 = 3. - 3. / 10. * std::sqrt(std::pow(x - 47.5, 2));
103 return cone_magnitude * std::max({z3, Number(0.)});
104 }
105
106 const Number &x = point[0];
107 const Number &y = point[1];
108
109 Number z1 =
110 1. -
111 1. / 8. * std::sqrt(std::pow(x - 30., 2) + std::pow(y - 6., 2));
112
113 Number z2 =
114 1. -
115 1. / 8. * std::sqrt(std::pow(x - 30., 2) + std::pow(y - 24., 2));
116
117 Number z3 =
118 3. -
119 3. / 10. * std::sqrt(std::pow(x - 47.5, 2) + std::pow(y - 15., 2));
120
121 return cone_magnitude * std::max({z1, z2, z3, Number(0.)});
122 }
123
124 bool well_balancing_validation;
125 Number left_depth;
126 Number right_depth;
127 Number cone_magnitude;
128 };
129
130 } // namespace ShallowWaterInitialStates
131} // namespace ryujin
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 Description::template HyperbolicSystemView< dim, Number > View
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32