10#include <compile_time_options.h>
14#include <deal.II/base/function_parser.h>
18 namespace ShallowWaterInitialStates
26 template <
typename Description,
int dim,
typename Number>
36 const std::string subsection)
38 , hyperbolic_system_(hyperbolic_system)
41 depth_expression_ =
"1.4";
43 "water elevation expression",
45 "A function expression describing the water elevation. When "
46 "bathymetry is 0, this reduces to the water depth.");
48 bathymetry_expression_ =
"0.";
49 this->add_parameter(
"bathymetry expression",
50 bathymetry_expression_,
51 "A function expression describing the bathymetry");
54 velocity_x_expression_ =
"3.0";
56 "velocity x expression",
57 velocity_x_expression_,
58 "A function expression describing the x-component of the velocity");
60 if constexpr (dim > 1) {
61 velocity_y_expression_ =
"0.0";
62 this->add_parameter(
"velocity y expression",
63 velocity_y_expression_,
64 "A function expression describing the "
65 "y-component of the velocity");
72 const auto set_up_muparser = [
this] {
73 using FP = dealii::FunctionParser<dim>;
79 depth_function_ = std::make_unique<FP>(depth_expression_);
80 velocity_x_function_ = std::make_unique<FP>(velocity_x_expression_);
81 if constexpr (dim > 1)
82 velocity_y_function_ = std::make_unique<FP>(velocity_y_expression_);
83 bathymetry_function_ = std::make_unique<FP>(bathymetry_expression_);
87 this->parse_parameters_call_back.connect(set_up_muparser);
92 const auto view = hyperbolic_system_.template view<dim, Number>();
96 const Number z = compute_bathymetry(point);
98 depth_function_->set_time(t);
100 std::max(Number(depth_function_->value(point)) - z, Number(0.));
102 velocity_x_function_->set_time(t);
103 full_primitive[1] = velocity_x_function_->value(point);
105 if constexpr (dim > 1) {
106 velocity_y_function_->set_time(t);
107 full_primitive[2] = velocity_y_function_->value(point);
110 return view.from_primitive_state(full_primitive);
115 initial_precomputed_type
final
118 return {compute_bathymetry(point)};
124 std::string depth_expression_;
125 std::string velocity_x_expression_;
126 std::string velocity_y_expression_;
127 std::string bathymetry_expression_;
129 std::unique_ptr<dealii::FunctionParser<dim>> depth_function_;
130 std::unique_ptr<dealii::FunctionParser<dim>> velocity_x_function_;
131 std::unique_ptr<dealii::FunctionParser<dim>> velocity_y_function_;
132 std::unique_ptr<dealii::FunctionParser<dim>> bathymetry_function_;
134 DEAL_II_ALWAYS_INLINE
inline Number
135 compute_bathymetry(
const dealii::Point<dim> &point)
const
137 bathymetry_function_->set_time(0.);
138 return bathymetry_function_->value(point);
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
state_type compute(const dealii::Point< dim > &point, Number t) final
typename Description::HyperbolicSystem HyperbolicSystem
Function(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
Euler::HyperbolicSystem HyperbolicSystem