14 namespace ShallowWaterInitialStates
24 template <
typename Description,
int dim,
typename Number>
34 const std::string subsection)
36 , hyperbolic_system_(hyperbolic_system)
38 with_bathymetry =
false;
39 this->add_parameter(
"with bathymetry",
41 "If set to true then the initial profile includes "
42 "the bathymetry for the steady vortex. ");
45 this->add_parameter(
"reference depth", depth_,
"Reference water depth");
49 "mach number", mach_number_,
"Mach number of unsteady vortex");
52 this->add_parameter(
"beta", beta_,
"vortex strength beta");
57 const auto view = hyperbolic_system_.template view<dim, Number>();
58 const auto gravity = view.gravity();
60 dealii::Point<2> point_bar;
61 point_bar[0] = point[0] - mach_number_ * t;
62 point_bar[1] = point[1];
64 const Number r_square = Number(point_bar.norm_square());
68 beta_ / Number(2. * M_PI) *
69 exp(Number(0.5) - Number(0.5) * r_square );
72 depth_ - Number(1. / (2. * gravity )) * factor * factor;
75 h -= compute_bathymetry(point);
77 const Number u = mach_number_ - factor * Number(point_bar[1]);
78 const Number v = factor * Number(point_bar[0]);
80 if constexpr (dim == 2)
84 AssertThrow(
false, dealii::ExcNotImplemented());
91 initial_precomputed_type
final
94 return {compute_bathymetry(point)};
100 DEAL_II_ALWAYS_INLINE
inline Number
101 compute_bathymetry(
const dealii::Point<dim> &point)
const
103 const Number r_square = Number(point.norm_square());
106 const Number factor =
107 beta_ / Number(2. * M_PI) *
108 exp(Number(0.5) - Number(0.5) * r_square );
112 bath = depth_ / 4. * factor;
117 bool with_bathymetry;
typename Description::HyperbolicSystem HyperbolicSystem
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
typename Description::template HyperbolicSystemView< dim, Number > View
state_type compute(const dealii::Point< dim > &point, Number t) final
typename View::state_type state_type
SmoothVortex(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
Euler::HyperbolicSystem HyperbolicSystem