ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
initial_state_paraboloid.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 <compile_time_options.h>
11
13
14namespace ryujin
15{
16 namespace ShallowWaterInitialStates
17 {
30 template <typename Description, int dim, typename Number>
31 class Paraboloid : public InitialState<Description, dim, Number>
32 {
33 public:
35 using View =
36 typename Description::template HyperbolicSystemView<dim, Number>;
37 using state_type = typename View::state_type;
38
39 Paraboloid(const HyperbolicSystem &hyperbolic_system,
40 const std::string subsection)
41 : InitialState<Description, dim, Number>("paraboloid", subsection)
42 , hyperbolic_system_(hyperbolic_system)
43 {
44 a_ = 1.;
45 this->add_parameter(
46 "free surface radius", a_, "Radius of the circular free surface");
47
48 h_0_ = 0.1;
49 this->add_parameter(
50 "water height", h_0_, "Water height at central point");
51
52 if constexpr (dim == 1) {
53
54 length_ = 10000.;
55 this->add_parameter(
56 "paraboloid length", length_, "Length of 1D paraboloid");
57
58 B_ = 2.;
59 this->add_parameter("speed", B_, "The 1D paraboloid speed");
60
61 } else {
62 eta_ = 0.5;
63 this->add_parameter("eta", eta_, "The eta parameter");
64 }
65 }
66
67 state_type compute(const dealii::Point<dim> &point, Number t) final
68 {
69 const auto view = hyperbolic_system_.template view<dim, Number>();
70
71 /* Common quantities */
72 const auto z = compute_bathymetry(point);
73 const auto g = view.gravity();
74 const Number omega = std::sqrt(2. * g * h_0_) / a_;
75 const Number &x = point[0];
76
77 /* Initialize primitive variables */
78 Number h, v_x, v_y = 0.;
79
80 /* Define slightly different profiles for each dimension */
81
82 if constexpr (dim == 1) {
83
84 const Number k = view.manning_friction_coefficient();
85 const Number p = std::sqrt(8. * g * h_0_) / a_;
86 const Number s = std::sqrt(p * p - k * k) / 2.;
87
88 auto term1 =
89 (a_ * a_ * B_ * B_) / (8. * g * g * h_0_) * std::exp(-k * t);
90 term1 *= (1. / 4. * k * k - s * s) * std::cos(2. * s * t) -
91 s * k * sin(2. * s * t);
92
93 const auto term2 = -(B_ * B_ / (4. * g)) * std::exp(-k * t);
94
95 auto term3 = -(B_ / g) * std::exp(-1. / 2. * k * t);
96 term3 *= (s * std::cos(s * t) + 1. / 2. * k * std::sin(s * t)) *
97 (point[0] - 1. / 2. * length_);
98
99 auto htilde = h_0_ - compute_bathymetry(point);
100 htilde += term1 + term2 + term3;
101
102 h = std::max(htilde, Number(0.));
103 v_x = B_ * std::exp(-1. / 2. * k * t) * std::sin(s * t);
104
105 return state_type{{h, h * v_x}};
106 } else if constexpr (dim == 2) {
107
108 const Number &y = point[1];
109
110 const Number elevation =
111 eta_ * h_0_ / (a_ * a_) *
112 (2. * x * std::cos(omega * t) + 2. * y * std::sin(omega * t));
113
114 h = std::max(elevation - z, Number(0.));
115 v_x = -eta_ * omega * std::sin(omega * t);
116 v_y = eta_ * omega * std::cos(omega * t);
117
118 return state_type{{h, h * v_x, h * v_y}};
119
120 } else {
121 AssertThrow(false, dealii::ExcNotImplemented());
122 __builtin_trap();
123 }
124 }
125
126
127 auto initial_precomputations(const dealii::Point<dim> &point) ->
129 initial_precomputed_type final
130 {
131 /* Compute bathymetry: */
132 return {compute_bathymetry(point)};
133 }
134
135 private:
136 const HyperbolicSystem &hyperbolic_system_;
137
138 DEAL_II_ALWAYS_INLINE inline Number
139 compute_bathymetry(const dealii::Point<dim> &point) const
140 {
141 if constexpr (dim == 1)
142 return h_0_ / (a_ * a_) * std::pow(point[0] - 0.5 * length_, 2);
143 else
144 return -h_0_ * (Number(1.) - point.norm_square() / (a_ * a_));
145 }
146
147 Number a_;
148 Number h_0_;
149 Number eta_;
150 Number length_;
151 Number B_;
152 };
153
154 } // namespace ShallowWaterInitialStates
155} // namespace ryujin
typename Description::HyperbolicSystem HyperbolicSystem
Paraboloid(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
state_type compute(const dealii::Point< dim > &point, Number t) final
typename Description::template HyperbolicSystemView< dim, Number > View
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
T pow(const T x, const T b)
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:34