ryujin 2.1.1 revision 9dcb748690310d6a540ebb8b066d1a0834fc7604
initial_state_three_bumps_dam_break.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: MIT
3// Copyright (C) 2020 - 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include "hyperbolic_system.h"
9#include <initial_state.h>
10
11namespace ryujin
12{
13 namespace ShallowWater
14 {
21 template <int dim, typename Number, typename state_type>
22 class ThreeBumpsDamBreak : public InitialState<dim, Number, state_type, 1>
23 {
24 public:
25 ThreeBumpsDamBreak(const HyperbolicSystem &hyperbolic_system,
26 const std::string s)
27 : InitialState<dim, Number, state_type, 1>("three bumps dam break", s)
28 , hyperbolic_system(hyperbolic_system)
29 {
30 well_balancing_validation = false;
31 this->add_parameter(
32 "well balancing validation",
33 well_balancing_validation,
34 "If set to true then the initial profile is returned for all "
35 "times "
36 "(t>0); otherwise a constant inflow is computed for t>0 suitable "
37 "for prescribing Dirichlet conditions at the inflow boundary.");
38
39
40 left_depth = 1.875;
41 this->add_parameter("left water depth",
42 left_depth,
43 "Depth of water to the left of pseudo-dam");
44 right_depth = 0.;
45 this->add_parameter("right water depth",
46 right_depth,
47 "Depth of water to the right of pseudo-dam");
48
49 cone_magnitude = 1.;
50 this->add_parameter("cone magnitude",
51 cone_magnitude,
52 "To modify magnitude of cone heights");
53 }
54
55 state_type compute(const dealii::Point<dim> &point, Number t) final
56 {
57 const Number x = point[0];
58
59 /* Initial state: */
60
61 if (t <= 1.e-10 || well_balancing_validation) {
62 Number h = x < 0 ? left_depth : right_depth;
63 h = std::max(h - compute_bathymetry(point), Number(0.));
64 return hyperbolic_system.template expand_state<dim>(
66 }
67
68 /* For t > 0 prescribe constant inflow Dirichlet data on the left: */
69
70 const auto &h = left_depth;
71 const auto a = hyperbolic_system.speed_of_sound(
72 HyperbolicSystem ::state_type<1, Number>{{h, Number(0.)}});
73 return hyperbolic_system.template expand_state<dim>(
75 }
76
77 auto initial_precomputations(const dealii::Point<dim> &point) ->
79 final
80 {
81 /* Compute bathymetry: */
82 return {compute_bathymetry(point)};
83 }
84
85 private:
86 const HyperbolicSystem &hyperbolic_system;
87
88 DEAL_II_ALWAYS_INLINE inline Number
89 compute_bathymetry(const dealii::Point<dim> &point) const
90 {
91 if constexpr (dim == 1) {
92 /* When dim = 1, we only have one cone */
93 const Number &x = point[0];
94
95 Number z3 = 3. - 3. / 10. * std::sqrt(std::pow(x - 47.5, 2));
96 return cone_magnitude * std::max({z3, Number(0.)});
97 }
98
99 const Number &x = point[0];
100 const Number &y = point[1];
101
102 Number z1 =
103 1. -
104 1. / 8. * std::sqrt(std::pow(x - 30., 2) + std::pow(y - 6., 2));
105
106 Number z2 =
107 1. -
108 1. / 8. * std::sqrt(std::pow(x - 30., 2) + std::pow(y - 24., 2));
109
110 Number z3 =
111 3. -
112 3. / 10. * std::sqrt(std::pow(x - 47.5, 2) + std::pow(y - 15., 2));
113
114 return cone_magnitude * std::max({z1, z2, z3, Number(0.)});
115 }
116
117 bool well_balancing_validation;
118 Number left_depth;
119 Number right_depth;
120 Number cone_magnitude;
121 };
122
123 } // namespace ShallowWater
124} // namespace ryujin
typename HyperbolicSystemView::state_type state_type
Number speed_of_sound(const dealii::Tensor< 1, problem_dim, Number > &U) const
dealii::Tensor< 1, problem_dimension< dim >, Number > state_type
state_type compute(const dealii::Point< dim > &point, Number t) final
ThreeBumpsDamBreak(const HyperbolicSystem &hyperbolic_system, const std::string s)
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< dim, Number, state_type, 1 >::precomputed_type final
T pow(const T x, const T b)