8#include <compile_time_options.h>
12#include <deal.II/distributed/solution_transfer.h>
26 template <
typename Description,
int dim,
typename Number =
double>
57 : offline_data_(&offline_data)
58 , hyperbolic_system_(hyperbolic_system)
59 , solution_transfer_(offline_data.dof_handler())
73 const auto &U = std::get<0>(state_vector);
75 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
76 const auto &affine_constraints = offline_data_->affine_constraints();
79 for (
auto &it : state_)
80 it.reinit(scalar_partitioner);
89 U.extract_component(state_[k], k);
90 affine_constraints.distribute(state_[k]);
91 state_[k].update_ghost_values();
94 std::vector<const ScalarVector *> ptr_state;
95 std::transform(state_.begin(),
97 std::back_inserter(ptr_state),
98 [](
auto &it) { return ⁢ });
99 solution_transfer_.prepare_for_coarsening_and_refinement(ptr_state);
111 auto &U = std::get<0>(state_vector);
113 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
115 U.reinit(offline_data_->vector_partitioner());
118 for (
auto &it : interpolated_state_) {
119 it.reinit(scalar_partitioner);
120 it.zero_out_ghost_values();
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 ⁢ });
128 solution_transfer_.interpolate(ptr_interpolated_state);
135 U.insert_component(interpolated_state_[k], k);
137 U.update_ghost_values();
146 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
149 dealii::parallel::distributed::SolutionTransfer<dim, ScalarVector>
152 std::vector<ScalarVector> state_;
153 std::vector<ScalarVector> interpolated_state_;
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
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Euler::HyperbolicSystem HyperbolicSystem