10#include <deal.II/base/function_parser.h>
14 namespace EulerInitialStates
21 template <
typename Description,
int dim,
typename Number>
31 const std::string subsection)
33 , hyperbolic_system_(hyperbolic_system)
36 density_expression_ =
"1.4";
37 this->add_parameter(
"density expression",
39 "A function expression describing the density");
41 velocity_x_expression_ =
"3.0";
43 "velocity x expression",
44 velocity_x_expression_,
45 "A function expression describing the x-component of the velocity");
47 if constexpr (dim > 1) {
48 velocity_y_expression_ =
"0.0";
49 this->add_parameter(
"velocity y expression",
50 velocity_y_expression_,
51 "A function expression describing the "
52 "y-component of the velocity");
55 if constexpr (dim > 2) {
56 velocity_z_expression_ =
"0.0";
57 this->add_parameter(
"velocity z expression",
58 velocity_z_expression_,
59 "A function expression describing the "
60 "z-component of the velocity");
63 pressure_expression_ =
"1.0";
64 this->add_parameter(
"pressure expression",
66 "A function expression describing the pressure");
72 const auto set_up_muparser = [
this] {
73 using FP = dealii::FunctionParser<dim>;
79 density_function_ = std::make_unique<FP>(density_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 if constexpr (dim > 2)
84 velocity_z_function_ = std::make_unique<FP>(velocity_z_expression_);
85 pressure_function_ = std::make_unique<FP>(pressure_expression_);
89 this->parse_parameters_call_back.connect(set_up_muparser);
94 const auto view = hyperbolic_system_.template view<dim, Number>();
97 density_function_->set_time(t);
98 full_primitive_state[0] = density_function_->value(point);
100 velocity_x_function_->set_time(t);
101 full_primitive_state[1] = velocity_x_function_->value(point);
103 if constexpr (dim > 1) {
104 velocity_y_function_->set_time(t);
105 full_primitive_state[2] = velocity_y_function_->value(point);
107 if constexpr (dim > 2) {
108 velocity_z_function_->set_time(t);
109 full_primitive_state[3] = velocity_z_function_->value(point);
112 pressure_function_->set_time(t);
113 full_primitive_state[1 + dim] = pressure_function_->value(point);
115 return view.from_primitive_state(full_primitive_state);
121 std::string density_expression_;
122 std::string velocity_x_expression_;
123 std::string velocity_y_expression_;
124 std::string velocity_z_expression_;
125 std::string pressure_expression_;
127 std::unique_ptr<dealii::FunctionParser<dim>> density_function_;
128 std::unique_ptr<dealii::FunctionParser<dim>> velocity_x_function_;
129 std::unique_ptr<dealii::FunctionParser<dim>> velocity_y_function_;
130 std::unique_ptr<dealii::FunctionParser<dim>> velocity_z_function_;
131 std::unique_ptr<dealii::FunctionParser<dim>> pressure_function_;
typename Description::HyperbolicSystem HyperbolicSystem
state_type compute(const dealii::Point< dim > &point, Number t) final
Function(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
typename Description::template HyperbolicSystemView< dim, Number > View
typename View::state_type state_type
Euler::HyperbolicSystem HyperbolicSystem