ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
initial_state_ritter_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
11
12namespace ryujin
13{
14 namespace ShallowWaterInitialStates
15 {
22 template <typename Description, int dim, typename Number>
23 class RitterDamBreak : public InitialState<Description, dim, Number>
24 {
25 public:
27 using View =
28 typename Description::template HyperbolicSystemView<dim, Number>;
29 using state_type = typename View::state_type;
30
31 RitterDamBreak(const HyperbolicSystem &hyperbolic_system,
32 const std::string subsection)
33 : InitialState<Description, dim, Number>("ritter dam break",
34 subsection)
35 , hyperbolic_system_(hyperbolic_system)
36 {
37 dealii::ParameterAcceptor::parse_parameters_call_back.connect(
39
40 t_initial_ = 0.1;
41 this->add_parameter("time initial",
42 t_initial_,
43 "Time at which initial state is prescribed");
44
45 left_depth = 0.005;
46 this->add_parameter("left water depth",
47 left_depth,
48 "Depth of water to the left of pseudo-dam (x<0)");
49 }
50
52 {
53 AssertThrow(t_initial_ > 0.,
54 dealii::ExcMessage("Expansion must be computed at an "
55 "initial time greater than 0."));
56 }
57
58 state_type compute(const dealii::Point<dim> &point, Number t) final
59 {
60 const auto view = hyperbolic_system_.template view<dim, Number>();
61 const auto g = view.gravity();
62
63 const auto x = point[0];
64
65 const Number aL = std::sqrt(g * left_depth);
66 const Number xA = -(t + t_initial_) * aL;
67 const Number xB = Number(2.) * (t + t_initial_) * aL;
68
69 const Number tmp = aL - x / (2. * (t + t_initial_));
70
71 const Number h_expansion = 4. / (9. * g) * tmp * tmp;
72 const Number v_expansion = 2. / 3. * (x / (t + t_initial_) + aL);
73
74 if (x <= xA)
75 return state_type{{left_depth, Number(0.)}};
76 else if (x <= xB)
77 return state_type{{h_expansion, h_expansion * v_expansion}};
78 else
79 return state_type{{Number(0.), Number(0.)}};
80 }
81
82 /* Default bathymetry of 0 */
83
84 private:
85 const HyperbolicSystem &hyperbolic_system_;
86
87 Number t_initial_;
88 Number left_depth;
89 };
90
91 } // namespace ShallowWaterInitialStates
92} // namespace ryujin
typename Description::template HyperbolicSystemView< dim, Number > View
state_type compute(const dealii::Point< dim > &point, Number t) final
RitterDamBreak(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32