ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
initial_state_function.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0
3// [LANL Copyright Statement]
4// Copyright (C) 2023 - 2024 by the ryujin authors
5// Copyright (C) 2023 - 2024 by Triad National Security, LLC
6//
7
8#pragma once
9
11
12#include <deal.II/base/function_parser.h>
13
14namespace ryujin
15{
16 namespace ShallowWaterInitialStates
17 {
24 template <typename Description, int dim, typename Number>
25 class Function : public InitialState<Description, dim, Number>
26 {
27 public:
29 using View =
30 typename Description::template HyperbolicSystemView<dim, Number>;
31 using state_type = typename View::state_type;
32
33 Function(const HyperbolicSystem &hyperbolic_system,
34 const std::string subsection)
35 : InitialState<Description, dim, Number>("function", subsection)
36 , hyperbolic_system_(hyperbolic_system)
37 {
38
39 depth_expression_ = "1.4";
40 this->add_parameter("water depth expression",
41 depth_expression_,
42 "A function expression describing the water depth");
43
44
45 velocity_x_expression_ = "3.0";
46 this->add_parameter(
47 "velocity x expression",
48 velocity_x_expression_,
49 "A function expression describing the x-component of the velocity");
50
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");
57 }
58
59 /*
60 * Set up the muparser object with the final flux description from
61 * the parameter file:
62 */
63 const auto set_up_muparser = [this] {
64 using FP = dealii::FunctionParser<dim>;
65 /*
66 * This variant of the constructor initializes the function
67 * parser with support for a time-dependent description involving
68 * a variable »t«:
69 */
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_);
74 };
75
76 set_up_muparser();
77 this->parse_parameters_call_back.connect(set_up_muparser);
78 }
79
80 state_type compute(const dealii::Point<dim> &point, Number t) final
81 {
82 const auto view = hyperbolic_system_.template view<dim, Number>();
83 state_type full_primitive;
84
85 depth_function_->set_time(t);
86 full_primitive[0] = depth_function_->value(point);
87
88 velocity_x_function_->set_time(t);
89 full_primitive[1] = velocity_x_function_->value(point);
90
91 if constexpr (dim > 1) {
92 velocity_y_function_->set_time(t);
93 full_primitive[2] = velocity_y_function_->value(point);
94 }
95
96 return view.from_primitive_state(full_primitive);
97 }
98
99 private:
100 const HyperbolicSystem &hyperbolic_system_;
101
102 std::string depth_expression_;
103 std::string velocity_x_expression_;
104 std::string velocity_y_expression_;
105 std::string bathymetry_expression_;
106
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_;
111 };
112 } // namespace ShallowWaterInitialStates
113} // namespace ryujin
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
Definition: description.h:32