12#include <deal.II/base/function_parser.h>
16 namespace ShallowWaterInitialStates
25 template <
typename Description,
int dim,
typename Number>
36 const std::string subsection)
38 , hyperbolic_system_(hyperbolic_system)
39 , geotiff_reader_(subsection +
"/geotiff")
41 height_expression_ =
"1.4";
43 "water height expression",
45 "A function expression describing the initial total water height");
47 velocity_expression_ =
"0.0";
49 "velocity expression",
51 "A function expression describing the initial velocity");
53 const auto set_up = [
this] {
54 using FP = dealii::FunctionParser<dim>;
60 height_function_ = std::make_unique<FP>(height_expression_);
61 velocity_function_ = std::make_unique<FP>(velocity_expression_);
65 this->parse_parameters_call_back.connect(set_up);
70 const auto z = geotiff_reader_.compute_height(point);
72 dealii::Tensor<1, 2, Number> primitive;
74 height_function_->set_time(t);
75 primitive[0] = std::max(0., height_function_->value(point) - z);
77 velocity_function_->set_time(t);
78 primitive[1] = velocity_function_->value(point);
80 const auto view = hyperbolic_system_.template view<dim, Number>();
81 return view.from_initial_state(primitive);
86 initial_precomputed_type
final
89 return {
static_cast<Number
>(geotiff_reader_.compute_height(point))};
98 std::string height_expression_;
99 std::string velocity_expression_;
103 std::unique_ptr<dealii::FunctionParser<dim>> height_function_;
104 std::unique_ptr<dealii::FunctionParser<dim>> velocity_function_;
typename Description::HyperbolicSystem HyperbolicSystem
auto initial_precomputations(const dealii::Point< dim > &point) -> typename InitialState< Description, dim, Number >::initial_precomputed_type final
typename View::state_type state_type
state_type compute(const dealii::Point< dim > &point, Number t) final
GeoTIFF(const HyperbolicSystem &hyperbolic_system, const std::string subsection)
typename Description::template HyperbolicSystemView< dim, Number > View
Euler::HyperbolicSystem HyperbolicSystem