ryujin 2.1.1 revision 9dcb748690310d6a540ebb8b066d1a0834fc7604
initial_state_sloping_ramp_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 {
20 template <int dim, typename Number, typename state_type>
21 class SlopingRampDamBreak : public InitialState<dim, Number, state_type, 1>
22 {
23 public:
24 SlopingRampDamBreak(const HyperbolicSystem &hyperbolic_system,
25 const std::string s)
26 : InitialState<dim, Number, state_type, 1>("sloping ramp", s)
27 , hyperbolic_system(hyperbolic_system)
28 {
29 inflow = false;
30 this->add_parameter(
31 "inflow",
32 inflow,
33 "If set to true then a constant inflow is computed for t>0 "
34 "suitable for prescribing Dirichlet conditions at the inflow "
35 "boundary.");
36
37
38 left_depth = 1.875;
39 this->add_parameter("left water depth",
40 left_depth,
41 "Depth of water to the left of pseudo-dam");
42 right_depth = 0.;
43 this->add_parameter("right water depth",
44 right_depth,
45 "Depth of water to the right of pseudo-dam");
46
47 ramp_slope = 1.;
48 this->add_parameter(
49 "ramp slope", ramp_slope, "To modify slope of ramp");
50
51 flat_length = 5.;
52 this->add_parameter("length of flat part",
53 flat_length,
54 "To modify length of flat portion");
55 }
56
57 state_type compute(const dealii::Point<dim> &point, Number t) final
58 {
59 const Number x = point[0];
60
61 /* Initial state: */
62
63 if (t <= 1.e-10 || !inflow) {
64 Number h = x < 0 ? left_depth : right_depth;
65 h = std::max(h - compute_bathymetry(point), Number(0.));
66 return hyperbolic_system.template expand_state<dim>(
68 }
69
70 /* For t > 0 prescribe constant inflow Dirichlet data on the left: */
71
72 const auto &h = left_depth;
73 const auto a = hyperbolic_system.speed_of_sound(
74 HyperbolicSystem ::state_type<1, Number>{{h, Number(0.)}});
75 return hyperbolic_system.template expand_state<dim>(
77 }
78
79 auto initial_precomputations(const dealii::Point<dim> &point) ->
81 final
82 {
83 /* Compute bathymetry: */
84 return {compute_bathymetry(point)};
85 }
86
87 private:
88 const HyperbolicSystem &hyperbolic_system;
89
90 DEAL_II_ALWAYS_INLINE inline Number
91 compute_bathymetry(const dealii::Point<dim> &point) const
92 {
93
94 const Number &x = point[0];
95 return std::max(ramp_slope * (x - flat_length), Number(0.));
96 }
97
98 bool inflow;
99 Number left_depth;
100 Number right_depth;
101 Number flat_length;
102 Number ramp_slope;
103 };
104
105 } // namespace ShallowWater
106} // 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
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< dim, Number, state_type, 1 >::precomputed_type final
state_type compute(const dealii::Point< dim > &point, Number t) final
SlopingRampDamBreak(const HyperbolicSystem &hyperbolic_system, const std::string s)