ryujin 2.1.1 revision 955e869188d49b3c97ca7b1cf4fd9ceb0e6f46ef
initial_values.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include "initial_values.h"
9
10#include <deal.II/numerics/vector_tools.h>
11#include <deal.II/numerics/vector_tools.templates.h>
12
13#include <random>
14
15namespace ryujin
16{
17 using namespace dealii;
18
19 template <typename Description, int dim, typename Number>
21 const MPIEnsemble &mpi_ensemble,
22 const OfflineData<dim, Number> &offline_data,
23 const HyperbolicSystem &hyperbolic_system,
24 const ParabolicSystem &parabolic_system,
25 const std::string &subsection)
26 : ParameterAcceptor(subsection)
27 , mpi_ensemble_(mpi_ensemble)
28 , offline_data_(&offline_data)
29 , hyperbolic_system_(&hyperbolic_system)
30 , parabolic_system_(&parabolic_system)
31 {
32 ParameterAcceptor::parse_parameters_call_back.connect(std::bind(
34 this));
35
36 configuration_ = "uniform";
37 add_parameter("configuration",
38 configuration_,
39 "The initial state configuration. Valid names are given by "
40 "any of the subsections defined below.");
41
42 initial_direction_[0] = 1.;
43 add_parameter(
44 "direction",
45 initial_direction_,
46 "Initial direction of initial configuration (Galilei transform)");
47
48 initial_position_[0] = 1.;
49 add_parameter(
50 "position",
51 initial_position_,
52 "Initial position of initial configuration (Galilei transform)");
53
54 perturbation_ = 0.;
55 add_parameter("perturbation",
56 perturbation_,
57 "Add a random perturbation of the specified magnitude to the "
58 "initial state.");
59
60 /*
61 * And finally populate the initial state list with all initial state
62 * configurations defined in the InitialStateLibrary namespace:
63 */
65 initial_state_list_,
66 *hyperbolic_system_,
67 *parabolic_system_,
68 subsection);
69 }
70
71 namespace
72 {
76 template <int dim>
77 inline DEAL_II_ALWAYS_INLINE dealii::Point<dim>
78 affine_transform(const dealii::Tensor<1, dim> initial_direction,
79 const dealii::Point<dim> initial_position,
80 const dealii::Point<dim> x)
81 {
82 auto direction = x - initial_position;
83
84 /* Roll third component of initial_direction onto xy-plane: */
85 if constexpr (dim == 3) {
86 auto n_x = initial_direction[0];
87 auto n_z = initial_direction[2];
88 const auto norm = std::sqrt(n_x * n_x + n_z * n_z);
89 n_x /= norm;
90 n_z /= norm;
91 auto new_direction = direction;
92 if (norm > 1.0e-14) {
93 new_direction[0] = n_x * direction[0] + n_z * direction[2];
94 new_direction[2] = -n_z * direction[0] + n_x * direction[2];
95 }
96 direction = new_direction;
97 }
98
99 /* Roll second component of initial_direction onto x-axis: */
100 if constexpr (dim >= 2) {
101 auto n_x = initial_direction[0];
102 auto n_y = initial_direction[1];
103 const auto norm = std::sqrt(n_x * n_x + n_y * n_y);
104 n_x /= norm;
105 n_y /= norm;
106 auto new_direction = direction;
107 if (norm > 1.0e-14) {
108 new_direction[0] = n_x * direction[0] + n_y * direction[1];
109 new_direction[1] = -n_y * direction[0] + n_x * direction[1];
110 }
111 direction = new_direction;
112 }
113
114 return Point<dim>() + direction;
115 }
116
117
121 template <int dim, typename Number>
122 inline DEAL_II_ALWAYS_INLINE dealii::Tensor<1, dim, Number>
123 affine_transform_vector(const dealii::Tensor<1, dim> initial_direction,
124 dealii::Tensor<1, dim, Number> direction)
125 {
126 if constexpr (dim >= 2) {
127 auto n_x = initial_direction[0];
128 auto n_y = initial_direction[1];
129 const auto norm = std::sqrt(n_x * n_x + n_y * n_y);
130 n_x /= norm;
131 n_y /= norm;
132 auto new_direction = direction;
133 if (norm > 1.0e-14) {
134 new_direction[0] = n_x * direction[0] - n_y * direction[1];
135 new_direction[1] = n_y * direction[0] + n_x * direction[1];
136 }
137 direction = new_direction;
138 }
139
140 if constexpr (dim == 3) {
141 auto n_x = initial_direction[0];
142 auto n_z = initial_direction[2];
143 const auto norm = std::sqrt(n_x * n_x + n_z * n_z);
144 n_x /= norm;
145 n_z /= norm;
146 auto new_direction = direction;
147 if (norm > 1.0e-14) {
148 new_direction[0] = n_x * direction[0] - n_z * direction[2];
149 new_direction[2] = n_z * direction[0] + n_x * direction[2];
150 }
151 direction = new_direction;
152 }
153
154 return direction;
155 }
156 } /* namespace */
157
158
159 template <typename Description, int dim, typename Number>
161 {
162 /* First, let's normalize the direction: */
163
164 AssertThrow(initial_direction_.norm() != 0.,
165 ExcMessage("Initial direction is set to the zero vector."));
166 initial_direction_ /= initial_direction_.norm();
167
168 /* Populate std::function object: */
169
170 {
171 bool initialized = false;
172 for (auto &it : initial_state_list_)
173 if (it->name() == configuration_) {
174 initial_state_ = [this, &it](const dealii::Point<dim> &point,
175 Number t) {
176 const auto transformed_point =
177 affine_transform(initial_direction_, initial_position_, point);
178 auto state = it->compute(transformed_point, t);
179 const auto view = hyperbolic_system_->template view<dim, Number>();
180 state =
181 view.apply_galilei_transform(state, [&](const auto &momentum) {
182 return affine_transform_vector(initial_direction_, momentum);
183 });
184 return state;
185 };
186
187 initial_precomputed_ = [this, &it](const dealii::Point<dim> &point) {
188 const auto transformed_point =
189 affine_transform(initial_direction_, initial_position_, point);
190 return it->initial_precomputations(transformed_point);
191 };
192
193 initialized = true;
194 break;
195 }
196
197 AssertThrow(
198 initialized,
199 ExcMessage(
200 "Could not find an initial state description with name \"" +
201 configuration_ + "\""));
202 }
203
204 /* Add a random perturbation to the original function object: */
205
206 if (perturbation_ != 0.) {
207 initial_state_ = [old_state = this->initial_state_,
208 perturbation = this->perturbation_](
209 const dealii::Point<dim> &point, Number t) {
210 auto state = old_state(point, t);
211
212 if (t > 0.)
213 return state;
214
215 static std::default_random_engine generator =
216 std::default_random_engine(std::random_device()());
217 static std::uniform_real_distribution<Number> distribution(-1., 1.);
218 static auto draw = std::bind(distribution, generator);
219 for (unsigned int i = 0; i < problem_dimension; ++i)
220 state[i] *= (Number(1.) + perturbation * draw());
221
222 return state;
223 };
224 }
225 }
226
227
228 template <typename Description, int dim, typename Number>
230 Number t) const -> HyperbolicVector
231 {
232#ifdef DEBUG_OUTPUT
233 std::cout << "InitialValues<dim, Number>::"
234 << "interpolate_hyperbolic_vector(t = " << t << ")" << std::endl;
235#endif
236
238 U.reinit(offline_data_->hyperbolic_vector_partitioner());
239
241
242 const auto callable = [&](const auto &p) { return initial_state(p, t); };
243
244 ScalarVector temp;
245 const auto scalar_partitioner = offline_data_->scalar_partitioner();
246 temp.reinit(scalar_partitioner);
247
248 for (unsigned int d = 0; d < problem_dimension; ++d) {
249 VectorTools::interpolate(offline_data_->dof_handler(),
250 to_function<dim, Number>(callable, d),
251 temp);
252 U.insert_component(temp, d);
253 }
254
255 U.update_ghost_values();
256
257 return U;
258 }
259
260
261 template <typename Description, int dim, typename Number>
264 {
265#ifdef DEBUG_OUTPUT
266 std::cout << "InitialValues<dim, Number>::"
267 << "interpolate_initial_precomputed_vector()" << std::endl;
268#endif
269
270 const auto scalar_partitioner = offline_data_->scalar_partitioner();
271
272 InitialPrecomputedVector precomputed;
273 precomputed.reinit_with_scalar_partitioner(scalar_partitioner);
274
275 if constexpr (n_initial_precomputed_values == 0)
276 return precomputed;
277
279
280 const auto callable = [&](const auto &p) { return initial_precomputed(p); };
281
282 ScalarVector temp;
283 temp.reinit(scalar_partitioner);
284
285 for (unsigned int d = 0; d < n_initial_precomputed_values; ++d) {
286 VectorTools::interpolate(offline_data_->dof_handler(),
287 to_function<dim, Number>(callable, d),
288 temp);
289 precomputed.insert_component(temp, d);
290 }
291
292 precomputed.update_ghost_values();
293 return precomputed;
294 }
295
296} /* namespace ryujin */
static void populate_initial_state_list(initial_state_list_type &initial_state_list, const HyperbolicSystem &h, const ParabolicSystem &p, const std::string &s)
InitialValues(const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem &parabolic_system, const std::string &subsection="/InitialValues")
InitialPrecomputedVector interpolate_initial_precomputed_vector() const
typename View::InitialPrecomputedVector InitialPrecomputedVector
typename Description::HyperbolicSystem HyperbolicSystem
typename Description::ParabolicSystem ParabolicSystem
HyperbolicVector interpolate_hyperbolic_vector(Number t=0) const
typename View::HyperbolicVector HyperbolicVector
Vectors::ScalarVector< Number > ScalarVector
Definition: offline_data.h:58
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
Definition: state_vector.h:31