10#include <deal.II/numerics/vector_tools.h>
11#include <deal.II/numerics/vector_tools.templates.h>
17 using namespace dealii;
19 template <
typename Description,
int dim,
typename Number>
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_(¶bolic_system)
32 ParameterAcceptor::parse_parameters_call_back.connect(std::bind(
36 configuration_ =
"uniform";
37 add_parameter(
"configuration",
39 "The initial state configuration. Valid names are given by "
40 "any of the subsections defined below.");
42 initial_direction_[0] = 1.;
46 "Initial direction of initial configuration (Galilei transform)");
48 initial_position_[0] = 1.;
52 "Initial position of initial configuration (Galilei transform)");
55 add_parameter(
"perturbation",
57 "Add a random perturbation of the specified magnitude to the "
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)
82 auto direction = x - initial_position;
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);
91 auto new_direction = direction;
93 new_direction[0] = n_x * direction[0] + n_z * direction[2];
94 new_direction[2] = -n_z * direction[0] + n_x * direction[2];
96 direction = new_direction;
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);
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];
111 direction = new_direction;
114 return Point<dim>() + direction;
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)
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);
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];
137 direction = new_direction;
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);
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];
151 direction = new_direction;
159 template <
typename Description,
int dim,
typename Number>
164 AssertThrow(initial_direction_.norm() != 0.,
165 ExcMessage(
"Initial direction is set to the zero vector."));
166 initial_direction_ /= initial_direction_.norm();
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,
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>();
181 view.apply_galilei_transform(state, [&](
const auto &momentum) {
182 return affine_transform_vector(initial_direction_, momentum);
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);
200 "Could not find an initial state description with name \"" +
201 configuration_ +
"\""));
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);
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());
228 template <
typename Description,
int dim,
typename Number>
233 std::cout <<
"InitialValues<dim, Number>::"
234 <<
"interpolate_hyperbolic_vector(t = " << t <<
")" << std::endl;
238 U.reinit(offline_data_->hyperbolic_vector_partitioner());
242 const auto callable = [&](
const auto &p) {
return initial_state(p, t); };
245 const auto scalar_partitioner = offline_data_->scalar_partitioner();
246 temp.reinit(scalar_partitioner);
248 for (
unsigned int d = 0; d < problem_dimension; ++d) {
249 VectorTools::interpolate(offline_data_->dof_handler(),
250 to_function<dim, Number>(callable, d),
252 U.insert_component(temp, d);
255 U.update_ghost_values();
261 template <
typename Description,
int dim,
typename Number>
266 std::cout <<
"InitialValues<dim, Number>::"
267 <<
"interpolate_initial_precomputed_vector()" << std::endl;
270 const auto scalar_partitioner = offline_data_->scalar_partitioner();
273 precomputed.reinit_with_scalar_partitioner(scalar_partitioner);
275 if constexpr (n_initial_precomputed_values == 0)
280 const auto callable = [&](
const auto &p) {
return initial_precomputed(p); };
283 temp.reinit(scalar_partitioner);
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),
289 precomputed.insert_component(temp, d);
292 precomputed.update_ghost_values();
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 ¶bolic_system, const std::string &subsection="/InitialValues")
InitialPrecomputedVector interpolate_initial_precomputed_vector() const
typename View::InitialPrecomputedVector InitialPrecomputedVector
void parse_parameters_callback()
typename Description::HyperbolicSystem HyperbolicSystem
typename Description::ParabolicSystem ParabolicSystem
HyperbolicVector interpolate_hyperbolic_vector(Number t=0) const
typename View::HyperbolicVector HyperbolicVector
Vectors::ScalarVector< Number > ScalarVector
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector