ryujin 2.1.1 revision 5df067c50838bd7aa4d14c27f2b7e7974dc646f3
solution_transfer.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
10#include "offline_data.h"
11
12#include <deal.II/distributed/solution_transfer.h>
13
14namespace ryujin
15{
16
26 template <typename Description, int dim, typename Number = double>
27 class SolutionTransfer final
28 {
29 public:
34
36
37 using View =
38 typename Description::template HyperbolicSystemView<dim, Number>;
39
40 static constexpr auto problem_dimension = View::problem_dimension;
41
43
45
47
51
56 const HyperbolicSystem &hyperbolic_system)
57 : offline_data_(&offline_data)
58 , hyperbolic_system_(hyperbolic_system)
59 , solution_transfer_(offline_data.dof_handler())
60 {
61 }
62
71 void prepare_for_interpolation(const StateVector &state_vector)
72 {
73 const auto &U = std::get<0>(state_vector);
74
75 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
76 const auto &affine_constraints = offline_data_->affine_constraints();
77
78 state_.resize(problem_dimension);
79 for (auto &it : state_)
80 it.reinit(scalar_partitioner);
81
82 /*
83 * We need to copy over to an auxiliary state vector formed by a
84 * ScalarVector for each component because dealii::SolutionTransfer
85 * cannot work on our StateVector or MultiComponentVector
86 */
87
88 for (unsigned int k = 0; k < problem_dimension; ++k) {
89 U.extract_component(state_[k], k);
90 affine_constraints.distribute(state_[k]);
91 state_[k].update_ghost_values();
92 }
93
94 std::vector<const ScalarVector *> ptr_state;
95 std::transform(state_.begin(),
96 state_.end(),
97 std::back_inserter(ptr_state),
98 [](auto &it) { return &it; });
99 solution_transfer_.prepare_for_coarsening_and_refinement(ptr_state);
100 }
101
109 void interpolate(StateVector &state_vector)
110 {
111 auto &U = std::get<0>(state_vector);
112
113 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
114
115 U.reinit(offline_data_->vector_partitioner());
116
117 interpolated_state_.resize(problem_dimension);
118 for (auto &it : interpolated_state_) {
119 it.reinit(scalar_partitioner);
120 it.zero_out_ghost_values();
121 }
122
123 std::vector<ScalarVector *> ptr_interpolated_state;
124 std::transform(interpolated_state_.begin(),
125 interpolated_state_.end(),
126 std::back_inserter(ptr_interpolated_state),
127 [](auto &it) { return &it; });
128 solution_transfer_.interpolate(ptr_interpolated_state);
129
130 /*
131 * Read back from interpolated_state_ vector:
132 */
133
134 for (unsigned int k = 0; k < problem_dimension; ++k) {
135 U.insert_component(interpolated_state_[k], k);
136 }
137 U.update_ghost_values();
138 }
139
140 private:
142
146 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
147 const HyperbolicSystem &hyperbolic_system_;
148
149 dealii::parallel::distributed::SolutionTransfer<dim, ScalarVector>
150 solution_transfer_;
151
152 std::vector<ScalarVector> state_;
153 std::vector<ScalarVector> interpolated_state_;
155 };
156
157} /* namespace ryujin */
SolutionTransfer(const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system)
Vectors::ScalarVector< Number > ScalarVector
void interpolate(StateVector &state_vector)
typename Description::template HyperbolicSystemView< dim, Number > View
static constexpr auto problem_dimension
void prepare_for_interpolation(const StateVector &state_vector)
typename Description::HyperbolicSystem HyperbolicSystem
typename View::StateVector StateVector
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
Definition: state_vector.h:25
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Definition: state_vector.h:45
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32