12#include <deal.II/base/function_parser.h>
16 namespace ShallowWaterInitialStates
24 template <
typename Description,
int dim,
typename Number>
34 const std::string subsection)
36 , hyperbolic_system_(hyperbolic_system)
39 depth_expression_ =
"1.4";
40 this->add_parameter(
"water depth expression",
42 "A function expression describing the water depth");
45 velocity_x_expression_ =
"3.0";
47 "velocity x expression",
48 velocity_x_expression_,
49 "A function expression describing the x-component of the velocity");
51 if constexpr (dim > 1) {
52 velocity_y_expression_ =
"0.0";
53 this->add_parameter(
"velocity y expression",
54 velocity_y_expression_,
55 "A function expression describing the "
56 "y-component of the velocity");
63 const auto set_up_muparser = [
this] {
64 using FP = dealii::FunctionParser<dim>;
70 depth_function_ = std::make_unique<FP>(depth_expression_);
71 velocity_x_function_ = std::make_unique<FP>(velocity_x_expression_);
72 if constexpr (dim > 1)
73 velocity_y_function_ = std::make_unique<FP>(velocity_y_expression_);
77 this->parse_parameters_call_back.connect(set_up_muparser);
82 const auto view = hyperbolic_system_.template view<dim, Number>();
85 depth_function_->set_time(t);
86 full_primitive[0] = depth_function_->value(point);
88 velocity_x_function_->set_time(t);
89 full_primitive[1] = velocity_x_function_->value(point);
91 if constexpr (dim > 1) {
92 velocity_y_function_->set_time(t);
93 full_primitive[2] = velocity_y_function_->value(point);
96 return view.from_primitive_state(full_primitive);
102 std::string depth_expression_;
103 std::string velocity_x_expression_;
104 std::string velocity_y_expression_;
105 std::string bathymetry_expression_;
107 std::unique_ptr<dealii::FunctionParser<dim>> depth_function_;
108 std::unique_ptr<dealii::FunctionParser<dim>> velocity_x_function_;
109 std::unique_ptr<dealii::FunctionParser<dim>> velocity_y_function_;
110 std::unique_ptr<dealii::FunctionParser<dim>> bathymetry_function_;
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
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