8#include <compile_time_options.h>
12#if DEAL_II_VERSION_GTE(9, 6, 0)
16#include <deal.II/base/config.h>
17#include <deal.II/distributed/tria.h>
18#include <deal.II/dofs/dof_accessor.h>
19#include <deal.II/dofs/dof_tools.h>
20#if DEAL_II_VERSION_GTE(9, 6, 0)
21#include <deal.II/grid/cell_status.h>
23#include <deal.II/grid/tria_accessor.h>
24#include <deal.II/grid/tria_iterator.h>
25#include <deal.II/lac/block_vector.h>
26#include <deal.II/lac/la_parallel_block_vector.h>
27#include <deal.II/lac/la_parallel_vector.h>
28#include <deal.II/lac/petsc_block_vector.h>
29#include <deal.II/lac/petsc_vector.h>
30#include <deal.II/lac/trilinos_parallel_block_vector.h>
31#include <deal.II/lac/trilinos_vector.h>
32#include <deal.II/lac/vector.h>
33#include <deal.II/matrix_free/fe_point_evaluation.h>
38 template <
typename Description,
int dim,
typename Number>
44 : mpi_ensemble_(mpi_ensemble)
45 , offline_data_(&offline_data)
46 , hyperbolic_system_(&hyperbolic_system)
47 , parabolic_system_(¶bolic_system)
48 , handle_(dealii::numbers::invalid_unsigned_int)
59 template <
typename state_type>
61 pack_state_values(
const std::vector<state_type> &state_values)
63 std::vector<char> buffer(
sizeof(
state_type) * state_values.size());
64 std::memcpy(buffer.data(), state_values.data(), buffer.size());
72 template <
typename state_type>
73 std::vector<state_type> unpack_state_values(
74 const boost::iterator_range<std::vector<char>::const_iterator>
77 const std::size_t n_bytes = data_range.size();
78 Assert(n_bytes %
sizeof(state_type) == 0, dealii::ExcInternalError());
79 std::vector<state_type> state_values(n_bytes /
sizeof(state_type));
80 std::memcpy(state_values.data(),
82 state_values.size() *
sizeof(state_type));
88 template <
typename Description,
int dim,
typename Number>
89 inline DEAL_II_ALWAYS_INLINE
auto
90 SolutionTransfer<Description, dim, Number>::get_tensor(
91 const HyperbolicVector &U,
const dealii::types::global_dof_index global_i)
94 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
95 const auto &affine_constraints = offline_data_->affine_constraints();
96 const auto local_i = scalar_partitioner->global_to_local(global_i);
97 if (affine_constraints.is_constrained(global_i)) {
99 const auto &line = *affine_constraints.get_constraint_entries(global_i);
100 for (
const auto &[global_k, c_k] : line) {
101 const auto local_k = scalar_partitioner->global_to_local(global_k);
102 result += c_k * U.get_tensor(local_k);
106 return U.get_tensor(local_i);
111 template <
typename Description,
int dim,
typename Number>
112 inline DEAL_II_ALWAYS_INLINE
void
113 SolutionTransfer<Description, dim, Number>::add_tensor(
115 const state_type &new_U_i,
116 const dealii::types::global_dof_index global_i)
118 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
119 const auto local_i = scalar_partitioner->global_to_local(global_i);
120 U.add_tensor(new_U_i, local_i);
124 template <
typename Description,
int dim,
typename Number>
126 const StateVector &old_state_vector [[maybe_unused]])
130 <<
"SolutionTransfer<Description, dim, Number>::prepare_projection()"
134#if !DEAL_II_VERSION_GTE(9, 6, 0)
138 "The SolutionTransfer class needs deal.II version 9.6.0 or newer"));
141 AssertThrow(have_distributed_triangulation<dim>,
143 "The SolutionTransfer class is not implemented for a "
144 "distributed::shared::Triangulation which we use in 1D"));
146 const auto &discretization = offline_data_->discretization();
147 auto &triangulation = *discretization.triangulation_;
149 Assert(handle_ == dealii::numbers::invalid_unsigned_int,
151 "You can only add one solution per SolutionTransfer object."));
157 handle_ = triangulation.register_data_attach(
158 [
this, &old_state_vector](
const auto cell,
159 const dealii::CellStatus status) {
160 const auto &dof_handler = offline_data_->dof_handler();
161 const auto dof_cell =
typename dealii::DoFHandler<dim>::cell_iterator(
162 &cell->get_triangulation(),
167 const auto &U = std::get<0>(old_state_vector);
173 const auto n_dofs_per_cell = dof_handler.get_fe().n_dofs_per_cell();
174 std::vector<state_type> state_values(n_dofs_per_cell);
177 case dealii::CellStatus::cell_will_persist:
179 case dealii::CellStatus::cell_will_be_refined: {
185 Assert(dof_cell->is_active(), dealii::ExcInternalError());
186 std::vector<dealii::types::global_dof_index> dof_indices(
188 dof_cell->get_dof_indices(dof_indices);
190 std::transform(std::begin(dof_indices),
191 std::end(dof_indices),
192 std::begin(state_values),
194 const auto U_i = get_tensor(U, i);
199 case dealii::CellStatus::children_will_be_coarsened: {
206 Assert(dof_cell->has_children(), dealii::ExcInternalError());
208 const auto &discretization = offline_data_->discretization();
209 const auto &finite_element = discretization.finite_element();
210 const auto &mapping = discretization.mapping();
211 const auto &quadrature = discretization.quadrature();
213 dealii::FEValues<dim> fe_values(
217 dealii::update_values | dealii::update_JxW_values |
218 dealii::update_quadrature_points);
220 const auto polynomial_space =
221 dealii::internal::FEPointEvaluation::get_polynomial_space(
224 std::vector<dealii::Point<dim, Number>> unit_points(
230 std::vector<dealii::Point<dim>> unit_points_temp(
231 std::is_same_v<Number, float> ? quadrature.size() : 0);
235 std::vector<state_type> state_values_quad(quadrature.size());
236 std::vector<state_type> local_rhs(n_dofs_per_cell);
238 std::vector<dealii::types::global_dof_index> dof_indices(
241 for (
unsigned int child = 0; child < dof_cell->n_children();
243 const auto child_cell = dof_cell->child(child);
244 Assert(child_cell->is_active(), dealii::ExcInternalError());
246 fe_values.reinit(child_cell);
248 if constexpr (std::is_same_v<Number, float>) {
249 mapping.transform_points_real_to_unit_cell(
251 fe_values.get_quadrature_points(),
253 std::transform(std::begin(unit_points_temp),
254 std::end(unit_points_temp),
255 std::begin(unit_points),
256 [](
const auto &x) {
return x; });
258 mapping.transform_points_real_to_unit_cell(
259 dof_cell, fe_values.get_quadrature_points(), unit_points);
262 child_cell->get_dof_indices(dof_indices);
264 for (
auto &it : state_values_quad)
267 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
268 const auto U_i = get_tensor(U, dof_indices[i]);
269 for (
unsigned int q = 0; q < quadrature.size(); ++q) {
270 state_values_quad[q] += U_i * fe_values.shape_value(i, q);
274 for (
unsigned int q = 0; q < quadrature.size(); ++q)
275 state_values_quad[q] *= fe_values.JxW(q);
277 for (
unsigned int q = 0; q < quadrature.size(); ++q) {
278 const unsigned int n_shapes = polynomial_space.size();
279 AssertIndexRange(n_shapes, 10);
280 dealii::ndarray<Number, 10, 2, dim> shapes;
282 std::array<Number, dim> point;
283 for (
unsigned int d = 0; d < dim; ++d)
284 point[d] = unit_points[q][d];
285 for (
unsigned int i = 0; i < n_shapes; ++i)
286 polynomial_space[i].values_of_array(point, 1, &shapes[i][0]);
288 Assert(finite_element.degree == 1, dealii::ExcNotImplemented());
294 state_type>(shapes.data(),
296 state_values_quad[q],
305 fe_values.reinit(dof_cell);
307 dealii::FullMatrix<double> mij(n_dofs_per_cell, n_dofs_per_cell);
308 dealii::Vector<double> mi(n_dofs_per_cell);
309 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
310 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j) {
312 for (
unsigned int q = 0; q < quadrature.size(); ++q)
313 sum += fe_values.shape_value(i, q) *
314 fe_values.shape_value(j, q) * fe_values.JxW(q);
322 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
323 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j) {
324 state_values[i] += mij(i, j) * local_rhs[j];
329 case dealii::CellStatus::cell_invalid:
330 Assert(
false, dealii::ExcInternalError());
335 return pack_state_values(state_values);
342 template <
typename Description,
int dim,
typename Number>
347 std::cout <<
"SolutionTransfer<Description, dim, Number>::project()"
351#if !DEAL_II_VERSION_GTE(9, 6, 0)
355 "The SolutionTransfer class needs deal.II version 9.6.0 or newer"));
359 AssertThrow(have_distributed_triangulation<dim>,
361 "The SolutionTransfer class is not implemented for a "
362 "distributed::shared::Triangulation which we use in 1D"));
364 const auto &scalar_partitioner = offline_data_->scalar_partitioner();
365 const auto &affine_constraints = offline_data_->affine_constraints();
367 const auto &discretization = offline_data_->discretization();
368 auto &triangulation = *discretization.triangulation_;
371 handle_ != dealii::numbers::invalid_unsigned_int,
373 "Cannot project() a state vector without valid handle. "
374 "prepare_projection() or set_handle() have to be called first."));
381 projected_mass.reinit(offline_data_->scalar_partitioner());
383 projected_state.reinit(offline_data_->hyperbolic_vector_partitioner());
385 triangulation.notify_ready_to_unpack(
387 [
this, &projected_mass, &projected_state](
389 const dealii::CellStatus status,
390 const auto &data_range) {
391 const auto &dof_handler = offline_data_->dof_handler();
392 const auto dof_cell =
typename dealii::DoFHandler<dim>::cell_iterator(
393 &cell->get_triangulation(),
402 const auto n_dofs_per_cell = dof_handler.get_fe().n_dofs_per_cell();
403 std::vector<dealii::types::global_dof_index> dof_indices(
406 const auto state_values = unpack_state_values<state_type>(data_range);
409 case dealii::CellStatus::cell_will_persist:
411 case dealii::CellStatus::children_will_be_coarsened: {
417 Assert(dof_cell->is_active(), dealii::ExcInternalError());
418 dof_cell->get_dof_indices(dof_indices);
420 const auto &discretization = offline_data_->discretization();
421 const auto &finite_element = discretization.finite_element();
422 const auto &mapping = discretization.mapping();
423 const auto &quadrature = discretization.quadrature();
425 dealii::FEValues<dim> fe_values(mapping,
428 dealii::update_values |
429 dealii::update_JxW_values);
431 fe_values.reinit(dof_cell);
433 dealii::Vector<double> mi(n_dofs_per_cell);
434 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
436 for (
unsigned int q = 0; q < quadrature.size(); ++q)
437 sum += fe_values.shape_value(i, q) * fe_values.JxW(q);
441 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
442 const auto global_i = dof_indices[i];
443 add_tensor(projected_state, mi(i) * state_values[i], global_i);
444 projected_mass(global_i) += mi(i);
449 case dealii::CellStatus::cell_will_be_refined: {
455 Assert(dof_cell->has_children(), dealii::ExcInternalError());
457 const auto &discretization = offline_data_->discretization();
458 const auto &finite_element = discretization.finite_element();
459 const auto &mapping = discretization.mapping();
460 const auto &quadrature = discretization.quadrature();
462 dealii::FEValues<dim> fe_values(
466 dealii::update_values | dealii::update_JxW_values |
467 dealii::update_quadrature_points);
469 const auto polynomial_space =
470 dealii::internal::FEPointEvaluation::get_polynomial_space(
472 std::vector<dealii::Point<dim, Number>> unit_points(
478 std::vector<dealii::Point<dim>> unit_points_temp(
479 std::is_same_v<Number, float> ? quadrature.size() : 0);
481 dealii::FullMatrix<double> mij(n_dofs_per_cell, n_dofs_per_cell);
482 dealii::Vector<double> mi(n_dofs_per_cell);
483 std::vector<state_type> local_rhs(n_dofs_per_cell);
485 for (
unsigned int child = 0; child < dof_cell->n_children();
487 const auto child_cell = dof_cell->child(child);
489 Assert(child_cell->is_active(), dealii::ExcInternalError());
490 child_cell->get_dof_indices(dof_indices);
494 fe_values.reinit(child_cell);
496 if constexpr (std::is_same_v<Number, float>) {
497 mapping.transform_points_real_to_unit_cell(
499 fe_values.get_quadrature_points(),
501 std::transform(std::begin(unit_points_temp),
502 std::end(unit_points_temp),
503 std::begin(unit_points),
504 [](
const auto &x) {
return x; });
506 mapping.transform_points_real_to_unit_cell(
507 dof_cell, fe_values.get_quadrature_points(), unit_points);
510 for (
auto &it : local_rhs)
513 for (
unsigned int q = 0; q < quadrature.size(); ++q) {
514 Assert(finite_element.degree == 1, dealii::ExcNotImplemented());
516 dealii::internal::evaluate_tensor_product_value(
518 make_const_array_view(state_values),
521 coefficient *= fe_values.JxW(q);
523 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
524 local_rhs[i] += coefficient * fe_values.shape_value(i, q);
531 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
532 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j) {
534 for (
unsigned int q = 0; q < quadrature.size(); ++q)
535 sum += fe_values.shape_value(i, q) *
536 fe_values.shape_value(j, q) * fe_values.JxW(q);
544 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i) {
546 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j) {
547 U_i += mij(i, j) * local_rhs[j];
550 const auto global_i = dof_indices[i];
551 add_tensor(projected_state, mi(i) * U_i, global_i);
552 projected_mass(global_i) += mi(i);
557 case dealii::CellStatus::cell_invalid:
558 Assert(
false, dealii::ExcInternalError());
569 auto &new_U = std::get<0>(new_state_vector);
570 const auto n_locally_owned = offline_data_->n_locally_owned();
574 const auto update_new_state_vector = [&]() {
578 for (
unsigned int local_i = 0; local_i < n_locally_owned; ++local_i) {
579 const auto global_i = scalar_partitioner->local_to_global(local_i);
580 if (affine_constraints.is_constrained(global_i))
583 const auto U_i = projected_state.get_tensor(local_i);
584 const auto m_i = projected_mass.local_element(local_i);
585 new_U.write_tensor(U_i / m_i, local_i);
587 new_U.update_ghost_values();
590 update_new_state_vector();
601 for (
const auto &line : affine_constraints.get_lines()) {
602 const auto global_i = line.index;
603 const auto local_i = scalar_partitioner->global_to_local(global_i);
606 if (local_i >= n_locally_owned)
610 const auto m_i_star = projected_mass.local_element(local_i);
611 const auto U_i_star = projected_state.get_tensor(local_i) / m_i_star;
615 for (
const auto &[global_k, c_k] : line.entries) {
616 const auto local_k = scalar_partitioner->global_to_local(global_k);
617 U_i_interp += c_k * new_U.get_tensor(local_k);
621 const auto defect = U_i_star - U_i_interp;
622 for (
const auto &[global_k, c_k] : line.entries) {
623 const auto local_k = scalar_partitioner->global_to_local(global_k);
624 const auto U_j = new_U.get_tensor(local_k);
626 projected_state.add_tensor(c_k * m_i_star * (U_j + defect), local_k);
627 projected_mass.local_element(local_k) += c_k * m_i_star;
631 update_new_state_vector();
637 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
638 for (
unsigned int local_i = 0; local_i < n_locally_owned; ++local_i) {
639 const auto global_i = scalar_partitioner->local_to_global(local_i);
640 if (affine_constraints.is_constrained(global_i))
643 const auto m_i = projected_mass.local_element(local_i);
644 const auto m_i_reference = lumped_mass_matrix.local_element(local_i);
645 Assert(std::abs(m_i - m_i_reference) < 1.e-10,
647 "SolutionTransfer::projection(): something went wrong. Final "
648 "masses do not agree with those computed in OfflineData."));
typename View::HyperbolicVector HyperbolicVector
Vectors::ScalarVector< Number > ScalarVector
typename View::state_type state_type
typename Description::ParabolicSystem ParabolicSystem
SolutionTransfer(const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system)
void prepare_projection(const StateVector &old_state_vector)
typename Description::HyperbolicSystem HyperbolicSystem
void project(StateVector &new_state_vector)
typename View::StateVector StateVector
void integrate_tensor_product_value(const dealii::ndarray< Number, 2, dim > *shapes, const unsigned int n_shapes, const Number2 &value, Number2 *values, const dealii::Point< dim, Number > &p, const bool do_add)
DEAL_II_ALWAYS_INLINE FT add(const FT &flux_left_ij, const FT &flux_right_ij)