10#include <compile_time_options.h>
16 namespace ShallowWaterInitialStates
30 template <
typename Description,
int dim,
typename Number>
40 const std::string subsection)
42 , hyperbolic_system_(hyperbolic_system)
46 "free surface radius", a_,
"Radius of the circular free surface");
50 "water height", h_0_,
"Water height at central point");
52 if constexpr (dim == 1) {
56 "paraboloid length", length_,
"Length of 1D paraboloid");
59 this->add_parameter(
"speed", B_,
"The 1D paraboloid speed");
63 this->add_parameter(
"eta", eta_,
"The eta parameter");
69 const auto view = hyperbolic_system_.template view<dim, Number>();
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];
78 Number h, v_x, v_y = 0.;
82 if constexpr (dim == 1) {
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.;
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);
93 const auto term2 = -(B_ * B_ / (4. * g)) * std::exp(-k * t);
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_);
99 auto htilde = h_0_ - compute_bathymetry(point);
100 htilde += term1 + term2 + term3;
102 h = std::max(htilde, Number(0.));
103 v_x = B_ * std::exp(-1. / 2. * k * t) * std::sin(s * t);
106 }
else if constexpr (dim == 2) {
108 const Number &y = point[1];
110 const Number elevation =
111 eta_ * h_0_ / (a_ * a_) *
112 (2. * x * std::cos(omega * t) + 2. * y * std::sin(omega * t));
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);
121 AssertThrow(
false, dealii::ExcNotImplemented());
129 initial_precomputed_type
final
132 return {compute_bathymetry(point)};
138 DEAL_II_ALWAYS_INLINE
inline Number
139 compute_bathymetry(
const dealii::Point<dim> &point)
const
141 if constexpr (dim == 1)
142 return h_0_ / (a_ * a_) *
std::pow(point[0] - 0.5 * length_, 2);
144 return -h_0_ * (Number(1.) - point.norm_square() / (a_ * a_));
typename Description::HyperbolicSystem HyperbolicSystem
typename View::state_type state_type
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