ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
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
10#include <compile_time_options.h>
11
13
14#include <deal.II/base/function_parser.h>
15
16namespace ryujin
17{
18 namespace ShallowWaterInitialStates
19 {
26 template <typename Description, int dim, typename Number>
27 class Function : public InitialState<Description, dim, Number>
28 {
29 public:
31 using View =
32 typename Description::template HyperbolicSystemView<dim, Number>;
33 using state_type = typename View::state_type;
34
35 Function(const HyperbolicSystem &hyperbolic_system,
36 const std::string subsection)
37 : InitialState<Description, dim, Number>("function", subsection)
38 , hyperbolic_system_(hyperbolic_system)
39 {
40
41 depth_expression_ = "1.4";
42 this->add_parameter(
43 "water elevation expression",
44 depth_expression_,
45 "A function expression describing the water elevation. When "
46 "bathymetry is 0, this reduces to the water depth.");
47
48 bathymetry_expression_ = "0.";
49 this->add_parameter("bathymetry expression",
50 bathymetry_expression_,
51 "A function expression describing the bathymetry");
52
53
54 velocity_x_expression_ = "3.0";
55 this->add_parameter(
56 "velocity x expression",
57 velocity_x_expression_,
58 "A function expression describing the x-component of the velocity");
59
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");
66 }
67
68 /*
69 * Set up the muparser object with the final flux description from
70 * the parameter file:
71 */
72 const auto set_up_muparser = [this] {
73 using FP = dealii::FunctionParser<dim>;
74 /*
75 * This variant of the constructor initializes the function
76 * parser with support for a time-dependent description involving
77 * a variable »t«:
78 */
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_);
84 };
85
86 set_up_muparser();
87 this->parse_parameters_call_back.connect(set_up_muparser);
88 }
89
90 state_type compute(const dealii::Point<dim> &point, Number t) final
91 {
92 const auto view = hyperbolic_system_.template view<dim, Number>();
93 state_type full_primitive;
94
95 /* Compute bathymetry */
96 const Number z = compute_bathymetry(point);
97
98 depth_function_->set_time(t);
99 full_primitive[0] =
100 std::max(Number(depth_function_->value(point)) - z, Number(0.));
101
102 velocity_x_function_->set_time(t);
103 full_primitive[1] = velocity_x_function_->value(point);
104
105 if constexpr (dim > 1) {
106 velocity_y_function_->set_time(t);
107 full_primitive[2] = velocity_y_function_->value(point);
108 }
109
110 return view.from_primitive_state(full_primitive);
111 }
112
113 auto initial_precomputations(const dealii::Point<dim> &point) ->
115 initial_precomputed_type final
116 {
117 /* Compute bathymetry: */
118 return {compute_bathymetry(point)};
119 }
120
121 private:
122 const HyperbolicSystem &hyperbolic_system_;
123
124 std::string depth_expression_;
125 std::string velocity_x_expression_;
126 std::string velocity_y_expression_;
127 std::string bathymetry_expression_;
128
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_;
133
134 DEAL_II_ALWAYS_INLINE inline Number
135 compute_bathymetry(const dealii::Point<dim> &point) const
136 {
137 bathymetry_function_->set_time(0.);
138 return bathymetry_function_->value(point);
139 }
140 };
141 } // namespace ShallowWaterInitialStates
142} // namespace ryujin
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
Definition: description.h:34