ryujin 2.1.1 revision a15074459a388761bd8df6bd4ef7e6abe9d8077b
initial_values.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
11#include "offline_data.h"
12
13#include <deal.II/base/parameter_acceptor.h>
14#include <deal.II/base/tensor.h>
15
16#include <functional>
17
18namespace ryujin
19{
36 template <typename Description, int dim, typename Number = double>
37 class InitialValues : public dealii::ParameterAcceptor
38 {
39 public:
44
48 using View =
49 typename Description::template HyperbolicSystemView<dim, Number>;
50
54 static constexpr unsigned int problem_dimension = View::problem_dimension;
55
59 using state_type = typename View::state_type;
60
65
69 static constexpr unsigned int n_precomputed_values =
70 View::n_precomputed_initial_values;
71
75 using precomputed_state_type = typename View::precomputed_state_type;
76
80 InitialValues(const HyperbolicSystem &hyperbolic_system,
81 const OfflineData<dim, Number> &offline_data,
82 const std::string &subsection = "/InitialValues");
83
84
92
93
101 DEAL_II_ALWAYS_INLINE inline state_type
102 initial_state(const dealii::Point<dim> &point, Number t) const
103 {
104 return initial_state_(point, t);
105 }
106
107
112 vector_type interpolate(Number t = 0) const;
113
114
122 DEAL_II_ALWAYS_INLINE inline precomputed_state_type
123 flux_contributions(const dealii::Point<dim> &point) const
124 {
125 return flux_contributions_(point);
126 }
127
128
135
136 private:
141
142 std::string configuration_;
143
144 dealii::Point<dim> initial_position_;
145
146 dealii::Tensor<1, dim> initial_direction_;
147
148 Number perturbation_;
149
151
155
156 dealii::SmartPointer<const HyperbolicSystem> hyperbolic_system_;
157 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
158
160 initial_state_list_type initial_state_list_;
161
162 std::function<state_type(const dealii::Point<dim> &point, Number t)>
163 initial_state_;
164
165 std::function<precomputed_state_type(const dealii::Point<dim> &point)>
166 flux_contributions_;
167
169 };
170
171} /* namespace ryujin */
static constexpr unsigned int n_precomputed_values
MultiComponentVector< Number, problem_dimension > vector_type
InitialValues(const HyperbolicSystem &hyperbolic_system, const OfflineData< dim, Number > &offline_data, const std::string &subsection="/InitialValues")
typename View::precomputed_state_type precomputed_state_type
DEAL_II_ALWAYS_INLINE precomputed_state_type flux_contributions(const dealii::Point< dim > &point) const
typename Description::HyperbolicSystem HyperbolicSystem
typename Description::template HyperbolicSystemView< dim, Number > View
vector_type interpolate(Number t=0) const
DEAL_II_ALWAYS_INLINE state_type initial_state(const dealii::Point< dim > &point, Number t) const
static constexpr unsigned int problem_dimension
MultiComponentVector< Number, n_precomputed_values > interpolate_precomputed_initial_values() const
typename View::state_type state_type
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32