ryujin 2.1.1 revision 336b16a72e829721302c626ec7071b92032b8248
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 velocity_expression_ = "3.0";
45 this->add_parameter("velocity expression",
46 velocity_expression_,
47 "A function expression describing the velocity");
48
49 /*
50 * Set up the muparser object with the final flux description from
51 * the parameter file:
52 */
53 const auto set_up_muparser = [this] {
54 using FP = dealii::FunctionParser<dim>;
55 /*
56 * This variant of the constructor initializes the function
57 * parser with support for a time-dependent description involving
58 * a variable »t«:
59 */
60 depth_function_ = std::make_unique<FP>(depth_expression_);
61 velocity_function_ = std::make_unique<FP>(velocity_expression_);
62 };
63
64 set_up_muparser();
65 this->parse_parameters_call_back.connect(set_up_muparser);
66 }
67
68 state_type compute(const dealii::Point<dim> &point, Number t) final
69 {
70 const auto view = hyperbolic_system_.template view<dim, Number>();
71
72 dealii::Tensor<1, 2, Number> primitive;
73
74 depth_function_->set_time(t);
75 primitive[0] = depth_function_->value(point);
76 velocity_function_->set_time(t);
77 primitive[1] = velocity_function_->value(point);
78
79 return view.from_initial_state(primitive);
80 }
81
82 private:
83 const HyperbolicSystem &hyperbolic_system_;
84
85 bool use_primitive_state_functions_;
86
87 std::string depth_expression_;
88 std::string velocity_expression_;
89 std::string bathymetry_expression_;
90
91 std::unique_ptr<dealii::FunctionParser<dim>> depth_function_;
92 std::unique_ptr<dealii::FunctionParser<dim>> velocity_function_;
93 std::unique_ptr<dealii::FunctionParser<dim>> bathymetry_function_;
94 };
95 } // namespace ShallowWaterInitialStates
96} // 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