15#include <deal.II/base/graph_coloring.h>
16#include <deal.II/base/parallel.h>
17#include <deal.II/base/work_stream.h>
18#include <deal.II/dofs/dof_renumbering.h>
19#include <deal.II/dofs/dof_tools.h>
20#include <deal.II/fe/fe_values.h>
21#include <deal.II/grid/grid_tools.h>
22#include <deal.II/lac/dynamic_sparsity_pattern.h>
23#include <deal.II/lac/la_parallel_vector.h>
24#ifdef DEAL_II_WITH_TRILINOS
25#include <deal.II/lac/trilinos_sparse_matrix.h>
28#ifdef FORCE_DEAL_II_SPARSE_MATRIX
29#undef DEAL_II_WITH_TRILINOS
34 using namespace dealii;
37 template <
int dim,
typename Number>
41 const std::string &subsection )
42 : ParameterAcceptor(subsection)
43 , mpi_ensemble_(mpi_ensemble)
44 , discretization_(&discretization)
46 incidence_relaxation_even_ = 0.5;
47 add_parameter(
"incidence matrix relaxation even degree",
48 incidence_relaxation_even_,
49 "Scaling exponent for incidence matrix used for "
50 "discontinuous finite elements with even degree. The default "
51 "value 0.5 scales the jump penalization with (h_i+h_j)^0.5.");
53 incidence_relaxation_odd_ = 0.0;
54 add_parameter(
"incidence matrix relaxation odd degree",
55 incidence_relaxation_odd_,
56 "Scaling exponent for incidence matrix used for "
57 "discontinuous finite elements with even degree. The default "
58 "value of 0.0 sets the jump penalization to a constant 1.");
62 template <
int dim,
typename Number>
71 auto &dof_handler = *dof_handler_;
72 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
74 IndexSet locally_relevant;
75 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
77 affine_constraints_.reinit(locally_relevant);
78 DoFTools::make_hanging_node_constraints(dof_handler, affine_constraints_);
80#ifndef DEAL_II_WITH_TRILINOS
81 AssertThrow(affine_constraints_.n_constraints() == 0,
82 ExcMessage(
"ryujin was built without Trilinos support - no "
83 "hanging node support available"));
91 const auto &periodic_faces =
92 discretization_->triangulation().get_periodic_face_map();
94 for (
const auto &[left, value] : periodic_faces) {
95 const auto &[right, orientation] = value;
97 typename DoFHandler<dim>::cell_iterator dof_cell_left(
98 &left.first->get_triangulation(),
103 typename DoFHandler<dim>::cell_iterator dof_cell_right(
104 &right.first->get_triangulation(),
105 right.first->level(),
106 right.first->index(),
109 if constexpr (std::is_same_v<Number, double>) {
110 DoFTools::make_periodicity_constraints(
111 dof_cell_left->face(left.second),
112 dof_cell_right->face(right.second),
115#
if DEAL_II_VERSION_GTE(9, 6, 0)
123 AssertThrow(
false, dealii::ExcNotImplemented());
128 affine_constraints_.close();
133 const std::vector<IndexSet> &locally_owned_dofs =
134 Utilities::MPI::all_gather(mpi_ensemble_.ensemble_communicator(),
135 dof_handler.locally_owned_dofs());
136 const IndexSet locally_active =
137 dealii::DoFTools::extract_locally_active_dofs(dof_handler);
138 Assert(affine_constraints_.is_consistent_in_parallel(
141 mpi_ensemble_.ensemble_communicator(),
147 sparsity_pattern_.reinit(
148 dof_handler.n_dofs(), dof_handler.n_dofs(), locally_relevant);
150 if (discretization_->have_discontinuous_ansatz()) {
155 dof_handler, sparsity_pattern_, affine_constraints_,
false);
160#ifdef DEAL_II_WITH_TRILINOS
161 DoFTools::make_sparsity_pattern(
162 dof_handler, sparsity_pattern_, affine_constraints_,
false);
171 dof_handler, sparsity_pattern_, affine_constraints_,
false);
181 SparsityTools::distribute_sparsity_pattern(
184 mpi_ensemble_.ensemble_communicator(),
189 template <
int dim,
typename Number>
190 void OfflineData<dim, Number>::setup(
const unsigned int problem_dimension,
191 const unsigned int n_precomputed_values)
194 std::cout <<
"OfflineData<dim, Number>::setup()" << std::endl;
201 const auto &triangulation = discretization_->triangulation();
203 dof_handler_ = std::make_unique<dealii::DoFHandler<dim>>(triangulation);
204 auto &dof_handler = *dof_handler_;
206 dof_handler.distribute_dofs(discretization_->finite_element());
208 n_locally_owned_ = dof_handler.locally_owned_dofs().n_elements();
214 DoFRenumbering::Cuthill_McKee(dof_handler);
225 mpi_ensemble_.ensemble_communicator(),
236 create_constraints_and_sparsity_pattern();
238 dof_handler, sparsity_pattern_, VectorizedArray<Number>::size());
253 mpi_ensemble_.ensemble_communicator(),
255 VectorizedArray<Number>::size());
261 const auto consistent_stride_range [[maybe_unused]] = [&]() {
262 constexpr auto group_size = VectorizedArray<Number>::size();
263 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
264 const auto offset = n_locally_owned_ != 0 ? *locally_owned.begin() : 0;
266 unsigned int group_row_length = 0;
268 for (; i < n_locally_internal_; ++i) {
269 if (i % group_size == 0) {
270 group_row_length = sparsity_pattern_.row_length(offset + i);
272 if (group_row_length != sparsity_pattern_.row_length(offset + i)) {
277 return i / group_size * group_size;
283 const auto mpi_allreduce_logical_or = [&](
const bool local_value) {
284 std::function<bool(
const bool &,
const bool &)> comparator =
285 [](
const bool &left,
const bool &right) ->
bool {
286 return left || right;
288 return Utilities::MPI::all_reduce(
289 local_value, mpi_ensemble_.ensemble_communicator(), comparator);
296 create_constraints_and_sparsity_pattern();
306#if DEAL_II_VERSION_GTE(9, 5, 0)
307 if (mpi_allreduce_logical_or(affine_constraints_.n_constraints() > 0)) {
308 if (mpi_allreduce_logical_or(
309 consistent_stride_range() != n_locally_internal_)) {
319 VectorizedArray<Number>::size());
320 create_constraints_and_sparsity_pattern();
321 n_locally_internal_ = consistent_stride_range();
331 Assert(consistent_stride_range() == n_locally_internal_,
332 dealii::ExcInternalError());
338 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
339 Assert(n_locally_owned_ == locally_owned.n_elements(),
340 dealii::ExcInternalError());
342 IndexSet locally_relevant;
343 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
346 IndexSet additional_dofs(dof_handler.n_dofs());
347 for (
auto &entry : sparsity_pattern_)
348 if (!locally_relevant.is_element(entry.column())) {
349 Assert(locally_owned.is_element(entry.row()), ExcInternalError());
350 additional_dofs.add_index(entry.column());
352 additional_dofs.compress();
353 locally_relevant.add_indices(additional_dofs);
354 locally_relevant.compress();
357 n_locally_relevant_ = locally_relevant.n_elements();
359 scalar_partitioner_ = std::make_shared<dealii::Utilities::MPI::Partitioner>(
360 locally_owned, locally_relevant, mpi_ensemble_.ensemble_communicator());
363 scalar_partitioner_, problem_dimension);
366 scalar_partitioner_, n_precomputed_values);
374 if (mpi_allreduce_logical_or(affine_constraints_.n_constraints() > 0)) {
378 n_export_indices_ = 0;
379 for (
const auto &it : scalar_partitioner_->import_indices())
380 if (it.second <= n_locally_internal_)
381 n_export_indices_ = std::max(n_export_indices_, it.second);
383 constexpr auto simd_length = VectorizedArray<Number>::size();
385 (n_export_indices_ + simd_length - 1) / simd_length * simd_length;
390 unsigned int control = 0;
391 for (
const auto &it : scalar_partitioner_->import_indices())
392 if (it.second <= n_locally_internal_)
393 control = std::max(control, it.second);
395 Assert(control <= n_export_indices_, ExcInternalError());
396 Assert(n_export_indices_ <= n_locally_internal_, ExcInternalError());
405 sparsity_pattern_simd_.reinit(
406 n_locally_internal_, sparsity_pattern_, scalar_partitioner_);
410 template <
int dim,
typename Number>
411 void OfflineData<dim, Number>::create_matrices()
414 std::cout <<
"OfflineData<dim, Number>::create_matrices()" << std::endl;
421 mass_matrix_.reinit(sparsity_pattern_simd_);
422 if (discretization_->have_discontinuous_ansatz())
423 mass_matrix_inverse_.reinit(sparsity_pattern_simd_);
425 lumped_mass_matrix_.reinit(scalar_partitioner_);
426 lumped_mass_matrix_inverse_.reinit(scalar_partitioner_);
428 cij_matrix_.reinit(sparsity_pattern_simd_);
429 if (discretization_->have_discontinuous_ansatz())
430 incidence_matrix_.reinit(sparsity_pattern_simd_);
436 auto &dof_handler = *dof_handler_;
438 measure_of_omega_ = 0.;
440#ifdef DEAL_II_WITH_TRILINOS
443 AffineConstraints<double> affine_constraints_assembly;
445 affine_constraints_assembly.reinit(affine_constraints_.get_local_lines());
446 for (
auto line : affine_constraints_.get_lines()) {
447 affine_constraints_assembly.add_line(line.index);
448 for (
auto entry : line.entries)
449 affine_constraints_assembly.add_entry(
450 line.index, entry.first, entry.second);
451 affine_constraints_assembly.set_inhomogeneity(line.index,
454 affine_constraints_assembly.close();
456 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
457 TrilinosWrappers::SparsityPattern trilinos_sparsity_pattern;
458 trilinos_sparsity_pattern.reinit(locally_owned,
460 mpi_ensemble_.ensemble_communicator());
462 TrilinosWrappers::SparseMatrix mass_matrix_tmp;
463 TrilinosWrappers::SparseMatrix mass_matrix_inverse_tmp;
464 if (discretization_->have_discontinuous_ansatz())
465 mass_matrix_inverse_tmp.reinit(trilinos_sparsity_pattern);
467 std::array<TrilinosWrappers::SparseMatrix, dim> cij_matrix_tmp;
469 mass_matrix_tmp.reinit(trilinos_sparsity_pattern);
470 for (
auto &matrix : cij_matrix_tmp)
471 matrix.reinit(trilinos_sparsity_pattern);
476 AffineConstraints<Number> affine_constraints_assembly;
477 affine_constraints_assembly.copy_from(affine_constraints_);
480 SparsityPattern sparsity_pattern_assembly;
482 DynamicSparsityPattern dsp(n_locally_relevant_, n_locally_relevant_);
483 for (
const auto &entry : sparsity_pattern_) {
484 const auto i = scalar_partitioner_->global_to_local(entry.row());
485 const auto j = scalar_partitioner_->global_to_local(entry.column());
488 sparsity_pattern_assembly.copy_from(dsp);
491 dealii::SparseMatrix<Number> mass_matrix_tmp;
492 dealii::SparseMatrix<Number> mass_matrix_inverse_tmp;
493 if (discretization_->have_discontinuous_ansatz())
494 mass_matrix_inverse_tmp.reinit(sparsity_pattern_assembly);
496 std::array<dealii::SparseMatrix<Number>, dim> cij_matrix_tmp;
498 mass_matrix_tmp.reinit(sparsity_pattern_assembly);
499 for (
auto &matrix : cij_matrix_tmp)
500 matrix.reinit(sparsity_pattern_assembly);
503 const unsigned int dofs_per_cell =
504 discretization_->finite_element().dofs_per_cell;
511 const auto local_assemble_system = [&](
const auto &cell,
516 auto &is_locally_owned = copy.is_locally_owned_;
517 auto &local_dof_indices = copy.local_dof_indices_;
518 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
520 auto &cell_mass_matrix = copy.cell_mass_matrix_;
521 auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
522 auto &cell_cij_matrix = copy.cell_cij_matrix_;
523 auto &interface_cij_matrix = copy.interface_cij_matrix_;
524 auto &cell_measure = copy.cell_measure_;
526 auto &fe_values = scratch.fe_values_;
527 auto &fe_face_values = scratch.fe_face_values_;
528 auto &fe_neighbor_face_values = scratch.fe_neighbor_face_values_;
530#ifdef DEAL_II_WITH_TRILINOS
531 is_locally_owned = cell->is_locally_owned();
539 is_locally_owned = !cell->is_artificial();
541 if (!is_locally_owned)
544 cell_mass_matrix.reinit(dofs_per_cell, dofs_per_cell);
545 for (
auto &matrix : cell_cij_matrix)
546 matrix.reinit(dofs_per_cell, dofs_per_cell);
547 if (discretization_->have_discontinuous_ansatz()) {
548 cell_mass_matrix_inverse.reinit(dofs_per_cell, dofs_per_cell);
549 for (
auto &it : interface_cij_matrix)
550 for (
auto &matrix : it)
551 matrix.reinit(dofs_per_cell, dofs_per_cell);
554 fe_values.reinit(cell);
556 local_dof_indices.resize(dofs_per_cell);
557 cell->get_dof_indices(local_dof_indices);
560 cell_mass_matrix = 0.;
561 for (
auto &matrix : cell_cij_matrix)
563 if (discretization_->have_discontinuous_ansatz()) {
564 cell_mass_matrix_inverse = 0.;
565 for (
auto &it : interface_cij_matrix)
566 for (
auto &matrix : it)
571 for (
unsigned int q : fe_values.quadrature_point_indices()) {
572 const auto JxW = fe_values.JxW(q);
574 if (cell->is_locally_owned())
575 cell_measure += Number(JxW);
577 for (
unsigned int j : fe_values.dof_indices()) {
578 const auto value_JxW = fe_values.shape_value(j, q) * JxW;
579 const auto grad_JxW = fe_values.shape_grad(j, q) * JxW;
581 for (
unsigned int i : fe_values.dof_indices()) {
582 const auto value = fe_values.shape_value(i, q);
584 cell_mass_matrix(i, j) += Number(value * value_JxW);
585 for (
unsigned int d = 0; d < dim; ++d)
586 cell_cij_matrix[d](i, j) += Number((value * grad_JxW)[d]);
596 if (!discretization_->have_discontinuous_ansatz())
599 for (
const auto f_index : cell->face_indices()) {
600 const auto &face = cell->face(f_index);
603 const bool has_neighbor =
604 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
608 neighbor_local_dof_indices[f_index].resize(0);
613 const auto neighbor_cell = cell->neighbor_or_periodic_neighbor(f_index);
614 if (neighbor_cell->is_artificial()) {
617 neighbor_local_dof_indices[f_index].resize(0);
621 fe_face_values.reinit(cell, f_index);
625 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
626 const auto JxW = fe_face_values.JxW(q);
627 const auto &normal = fe_face_values.get_normal_vectors()[q];
629 for (
unsigned int j : fe_face_values.dof_indices()) {
630 const auto value_JxW = fe_face_values.shape_value(j, q) * JxW;
632 for (
unsigned int i : fe_face_values.dof_indices()) {
633 const auto value = fe_face_values.shape_value(i, q);
635 for (
unsigned int d = 0; d < dim; ++d)
636 cell_cij_matrix[d](i, j) -=
637 Number(0.5 * normal[d] * value * value_JxW);
644 const unsigned int f_index_neighbor =
645 cell->has_periodic_neighbor(f_index)
646 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
647 : cell->neighbor_of_neighbor(f_index);
649 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
650 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
652 fe_neighbor_face_values.reinit(neighbor_cell, f_index_neighbor);
654 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
655 const auto JxW = fe_face_values.JxW(q);
656 const auto &normal = fe_face_values.get_normal_vectors()[q];
659 for (
unsigned int j : fe_face_values.dof_indices()) {
660 const auto value_JxW =
661 fe_neighbor_face_values.shape_value(j, q) * JxW;
663 for (
unsigned int i : fe_face_values.dof_indices()) {
664 const auto value = fe_face_values.shape_value(i, q);
666 for (
unsigned int d = 0; d < dim; ++d)
667 interface_cij_matrix[f_index][d](i, j) +=
668 Number(0.5 * normal[d] * value * value_JxW);
678 if (discretization_->have_discontinuous_ansatz()) {
680 cell_mass_matrix_inverse.invert(cell_mass_matrix);
684 const auto copy_local_to_global = [&](
const auto ©) {
685 const auto &is_locally_owned = copy.is_locally_owned_;
686#ifdef DEAL_II_WITH_TRILINOS
687 const auto &local_dof_indices = copy.local_dof_indices_;
688 const auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
695 auto local_dof_indices = copy.local_dof_indices_;
696 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
698 const auto &cell_mass_matrix = copy.cell_mass_matrix_;
699 const auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
700 const auto &cell_cij_matrix = copy.cell_cij_matrix_;
701 const auto &interface_cij_matrix = copy.interface_cij_matrix_;
702 const auto &cell_measure = copy.cell_measure_;
704 if (!is_locally_owned)
707#ifndef DEAL_II_WITH_TRILINOS
709 for (
auto &indices : neighbor_local_dof_indices)
713 affine_constraints_assembly.distribute_local_to_global(
714 cell_mass_matrix, local_dof_indices, mass_matrix_tmp);
716 for (
int k = 0; k < dim; ++k) {
717 affine_constraints_assembly.distribute_local_to_global(
718 cell_cij_matrix[k], local_dof_indices, cij_matrix_tmp[k]);
720 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
721 if (neighbor_local_dof_indices[f_index].size() != 0) {
722 affine_constraints_assembly.distribute_local_to_global(
723 interface_cij_matrix[f_index][k],
725 neighbor_local_dof_indices[f_index],
731 if (discretization_->have_discontinuous_ansatz())
732 affine_constraints_assembly.distribute_local_to_global(
733 cell_mass_matrix_inverse,
735 mass_matrix_inverse_tmp);
737 measure_of_omega_ += cell_measure;
740 WorkStream::run(dof_handler.begin_active(),
742 local_assemble_system,
743 copy_local_to_global,
744 AssemblyScratchData<dim>(*discretization_),
745#ifdef DEAL_II_WITH_TRILINOS
746 AssemblyCopyData<dim, double>());
748 AssemblyCopyData<dim, Number>());
751#ifdef DEAL_II_WITH_TRILINOS
753 for (
auto &it : cij_matrix_tmp)
756 mass_matrix_.read_in(mass_matrix_tmp,
false);
757 if (discretization_->have_discontinuous_ansatz())
758 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
false);
759 cij_matrix_.read_in(cij_matrix_tmp,
false);
761 mass_matrix_.read_in(mass_matrix_tmp,
true);
762 if (discretization_->have_discontinuous_ansatz())
763 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
true);
764 cij_matrix_.read_in(cij_matrix_tmp,
true);
767 mass_matrix_.update_ghost_rows();
768 if (discretization_->have_discontinuous_ansatz())
769 mass_matrix_inverse_.update_ghost_rows();
770 cij_matrix_.update_ghost_rows();
772 measure_of_omega_ = Utilities::MPI::sum(
773 measure_of_omega_, mpi_ensemble_.ensemble_communicator());
780#ifdef DEAL_II_WITH_TRILINOS
783 affine_constraints_assembly.set_zero(one);
785 ScalarVector local_lumped_mass_matrix(scalar_partitioner_);
786 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
789 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
791 lumped_mass_matrix_.local_element(i) =
792 local_lumped_mass_matrix.local_element(i);
793 lumped_mass_matrix_inverse_.local_element(i) =
794 1. / lumped_mass_matrix_.local_element(i);
796 lumped_mass_matrix_.update_ghost_values();
797 lumped_mass_matrix_inverse_.update_ghost_values();
801 dealii::Vector<Number> one(mass_matrix_tmp.m());
803 affine_constraints_assembly.set_zero(one);
805 dealii::Vector<Number> local_lumped_mass_matrix(mass_matrix_tmp.m());
806 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
808 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
810 lumped_mass_matrix_.local_element(i) = local_lumped_mass_matrix(i);
811 lumped_mass_matrix_inverse_.local_element(i) =
812 1. / lumped_mass_matrix_.local_element(i);
814 lumped_mass_matrix_.update_ghost_values();
815 lumped_mass_matrix_inverse_.update_ghost_values();
823 if (discretization_->have_discontinuous_ansatz()) {
824#ifdef DEAL_II_WITH_TRILINOS
825 TrilinosWrappers::SparseMatrix incidence_matrix_tmp;
826 incidence_matrix_tmp.reinit(trilinos_sparsity_pattern);
828 dealii::SparseMatrix<Number> incidence_matrix_tmp;
829 incidence_matrix_tmp.reinit(sparsity_pattern_assembly);
833 const auto local_assemble_system = [&](
const auto &cell,
838 auto &is_locally_owned = copy.is_locally_owned_;
839 auto &local_dof_indices = copy.local_dof_indices_;
840 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
841 auto &interface_incidence_matrix = copy.interface_incidence_matrix_;
842 auto &fe_face_values_nodal = scratch.fe_face_values_nodal_;
843 auto &fe_neighbor_face_values_nodal =
844 scratch.fe_neighbor_face_values_nodal_;
846#ifdef DEAL_II_WITH_TRILINOS
847 is_locally_owned = cell->is_locally_owned();
849 is_locally_owned = !cell->is_artificial();
851 if (!is_locally_owned)
854 for (
auto &matrix : interface_incidence_matrix)
855 matrix.reinit(dofs_per_cell, dofs_per_cell);
857 local_dof_indices.resize(dofs_per_cell);
858 cell->get_dof_indices(local_dof_indices);
861 for (
auto &matrix : interface_incidence_matrix)
864 for (
const auto f_index : cell->face_indices()) {
865 const auto &face = cell->face(f_index);
868 const bool has_neighbor =
869 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
873 neighbor_local_dof_indices[f_index].resize(0);
878 const auto neighbor_cell =
879 cell->neighbor_or_periodic_neighbor(f_index);
880 if (neighbor_cell->is_artificial()) {
883 neighbor_local_dof_indices[f_index].resize(0);
887 const unsigned int f_index_neighbor =
888 cell->has_periodic_neighbor(f_index)
889 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
890 : cell->neighbor_of_neighbor(f_index);
892 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
893 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
895 fe_face_values_nodal.reinit(cell, f_index);
896 fe_neighbor_face_values_nodal.reinit(neighbor_cell, f_index_neighbor);
900 for (
unsigned int q :
901 fe_face_values_nodal.quadrature_point_indices()) {
903 for (
unsigned int j : fe_neighbor_face_values_nodal.dof_indices()) {
904 const auto v_j = fe_neighbor_face_values_nodal.shape_value(j, q);
905 for (
unsigned int i : fe_face_values_nodal.dof_indices()) {
906 const auto v_i = fe_face_values_nodal.shape_value(i, q);
907 constexpr auto eps = std::numeric_limits<Number>::epsilon();
908 if (std::abs(v_i * v_j) > 100. * eps) {
909 const auto &ansatz = discretization_->ansatz();
911 const auto glob_i = local_dof_indices[i];
912 const auto glob_j = neighbor_local_dof_indices[f_index][j];
913 const auto m_i = lumped_mass_matrix_[glob_i];
914 const auto m_j = lumped_mass_matrix_[glob_j];
916 Number(0.5) * (m_i + m_j) / measure_of_omega_;
927 r_ij =
std::pow(hd_ij, incidence_relaxation_even_ / dim);
934 r_ij =
std::pow(hd_ij, incidence_relaxation_odd_ / dim);
937 interface_incidence_matrix[f_index](i, j) += r_ij;
945 const auto copy_local_to_global = [&](
const auto ©) {
946 const auto &is_locally_owned = copy.is_locally_owned_;
947#ifdef DEAL_II_WITH_TRILINOS
948 const auto &local_dof_indices = copy.local_dof_indices_;
949 const auto &neighbor_local_dof_indices =
950 copy.neighbor_local_dof_indices_;
957 auto local_dof_indices = copy.local_dof_indices_;
958 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
960 const auto &interface_incidence_matrix =
961 copy.interface_incidence_matrix_;
963 if (!is_locally_owned)
966#ifndef DEAL_II_WITH_TRILINOS
968 for (
auto &indices : neighbor_local_dof_indices)
972 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
973 if (neighbor_local_dof_indices[f_index].size() != 0) {
974 affine_constraints_assembly.distribute_local_to_global(
975 interface_incidence_matrix[f_index],
977 neighbor_local_dof_indices[f_index],
978 incidence_matrix_tmp);
983 WorkStream::run(dof_handler.begin_active(),
985 local_assemble_system,
986 copy_local_to_global,
987 AssemblyScratchData<dim>(*discretization_),
988#ifdef DEAL_II_WITH_TRILINOS
989 AssemblyCopyData<dim, double>());
991 AssemblyCopyData<dim, Number>());
994#ifdef DEAL_II_WITH_TRILINOS
995 incidence_matrix_.read_in(incidence_matrix_tmp,
false);
997 incidence_matrix_.read_in(incidence_matrix_tmp,
true);
999 incidence_matrix_.update_ghost_rows();
1006 boundary_map_ = construct_boundary_map(
1007 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1009 coupling_boundary_pairs_ = collect_coupling_boundary_pairs(
1010 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1013#ifdef DEBUG_SYMMETRY_CHECK
1018 double total_mass = 0.;
1019 for (
unsigned int i = 0; i < n_locally_owned_; ++i)
1020 total_mass += lumped_mass_matrix_.local_element(i);
1022 Utilities::MPI::sum(total_mass, mpi_ensemble_.ensemble_communicator());
1024 Assert(std::abs(measure_of_omega_ - total_mass) <
1025 1.e-12 * measure_of_omega_,
1027 "Total mass differs from the measure of the domain."));
1033 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1035 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1036 if (row_length == 1)
1040 mass_matrix_.get_entry(i, 0) - lumped_mass_matrix_.local_element(i);
1043 constexpr auto simd_length = VectorizedArray<Number>::size();
1044 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1045 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1046 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1048 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1050 const auto m_ij = mass_matrix_.get_entry(i, col_idx);
1051 if (discretization_->have_discontinuous_ansatz()) {
1054 Assert(std::abs(m_ij) > -1.e-12, dealii::ExcInternalError());
1056 Assert(std::abs(m_ij) > 1.e-12, dealii::ExcInternalError());
1060 const auto m_ji = mass_matrix_.get_transposed_entry(i, col_idx);
1061 if (std::abs(m_ij - m_ji) >= 1.e-12) {
1063 std::stringstream ss;
1064 ss <<
"m_ij matrix is not symmetric: " << m_ij <<
" <-> " << m_ji;
1065 Assert(
false, dealii::ExcMessage(ss.str()));
1069 Assert(std::abs(sum) < 1.e-12, dealii::ExcInternalError());
1076 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1078 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1079 if (row_length == 1)
1082 auto sum = cij_matrix_.get_tensor(i, 0);
1085 constexpr auto simd_length = VectorizedArray<Number>::size();
1086 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1087 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1088 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1090 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1092 const auto c_ij = cij_matrix_.get_tensor(i, col_idx);
1093 Assert(c_ij.norm() > 1.e-12, dealii::ExcInternalError());
1096 const auto c_ji = cij_matrix_.get_transposed_tensor(i, col_idx);
1097 if ((c_ij + c_ji).norm() >= 1.e-12) {
1101 CouplingDescription coupling{i, col_idx, j};
1102 const auto it = std::find(coupling_boundary_pairs_.begin(),
1103 coupling_boundary_pairs_.end(),
1105 if (it == coupling_boundary_pairs_.end()) {
1106 std::stringstream ss;
1107 ss <<
"c_ij matrix is not anti-symmetric: " << c_ij <<
" <-> "
1109 Assert(
false, dealii::ExcMessage(ss.str()));
1114 Assert(sum.norm() < 1.e-12, dealii::ExcInternalError());
1120 template <
int dim,
typename Number>
1121 void OfflineData<dim, Number>::create_multigrid_data()
1124 std::cout <<
"OfflineData<dim, Number>::create_multigrid_data()"
1128 auto &dof_handler = *dof_handler_;
1130 dof_handler.distribute_mg_dofs();
1132 const auto n_levels = dof_handler.get_triangulation().n_global_levels();
1134 AffineConstraints<float> level_constraints;
1137 level_boundary_map_.resize(n_levels);
1138 level_lumped_mass_matrix_.resize(n_levels);
1140 for (
unsigned int level = 0; level < n_levels; ++level) {
1143 IndexSet relevant_dofs;
1144 dealii::DoFTools::extract_locally_relevant_level_dofs(
1145 dof_handler, level, relevant_dofs);
1146 const auto partitioner = std::make_shared<Utilities::MPI::Partitioner>(
1147 dof_handler.locally_owned_mg_dofs(level),
1149 mpi_ensemble_.ensemble_communicator());
1150 level_lumped_mass_matrix_[level].reinit(partitioner);
1151 std::vector<types::global_dof_index> dof_indices(
1152 dof_handler.get_fe().dofs_per_cell);
1153 dealii::Vector<Number> mass_values(dof_handler.get_fe().dofs_per_cell);
1154 FEValues<dim> fe_values(discretization_->mapping(),
1155 discretization_->finite_element(),
1156 discretization_->quadrature(),
1157 update_values | update_JxW_values);
1158 for (
const auto &cell : dof_handler.cell_iterators_on_level(level))
1161 if (cell->is_locally_owned_on_level()) {
1162 fe_values.reinit(cell);
1163 for (
unsigned int i = 0; i < mass_values.size(); ++i) {
1165 for (
unsigned int q = 0; q < fe_values.n_quadrature_points; ++q)
1166 sum += fe_values.shape_value(i, q) * fe_values.JxW(q);
1167 mass_values(i) = sum;
1169 cell->get_mg_dof_indices(dof_indices);
1170 level_constraints.distribute_local_to_global(
1171 mass_values, dof_indices, level_lumped_mass_matrix_[level]);
1177 level_boundary_map_[level] = construct_boundary_map(
1178 dof_handler.begin_mg(level), dof_handler.end_mg(level), *partitioner);
1183 template <
int dim,
typename Number>
1184 template <
typename ITERATOR1,
typename ITERATOR2>
1186 const ITERATOR1 &begin,
1187 const ITERATOR2 &end,
1188 const Utilities::MPI::Partitioner &partitioner)
const -> BoundaryMap
1191 std::cout <<
"OfflineData<dim, Number>::construct_boundary_map()"
1199 using BoundaryData = std::tuple<dealii::Tensor<1, dim, Number> ,
1202 dealii::types::boundary_id ,
1203 dealii::Point<dim>> ;
1204 std::multimap<unsigned int, BoundaryData> preliminary_map;
1206 std::vector<dealii::types::global_dof_index> local_dof_indices;
1208 const dealii::QGauss<dim - 1> face_quadrature(3);
1209 dealii::FEFaceValues<dim> fe_face_values(discretization_->mapping(),
1210 discretization_->finite_element(),
1212 dealii::update_normal_vectors |
1213 dealii::update_values |
1214 dealii::update_JxW_values);
1216 const unsigned int dofs_per_cell =
1217 discretization_->finite_element().dofs_per_cell;
1219 const auto support_points =
1220 discretization_->finite_element().get_unit_support_points();
1222 for (
auto cell = begin; cell != end; ++cell) {
1230 if ((cell->is_active() && cell->is_artificial()) ||
1231 (!cell->is_active() && cell->is_artificial_on_level()))
1234 local_dof_indices.resize(dofs_per_cell);
1235 cell->get_active_or_mg_dof_indices(local_dof_indices);
1237 for (
auto f : cell->face_indices()) {
1238 const auto face = cell->face(f);
1239 const auto id = face->boundary_id();
1241 if (!face->at_boundary())
1252 fe_face_values.reinit(cell, f);
1254 for (
unsigned int j : fe_face_values.dof_indices()) {
1255 if (!discretization_->finite_element().has_support_on_face(j, f))
1258 Number boundary_mass = 0.;
1259 dealii::Tensor<1, dim, Number> normal;
1261 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
1262 const auto JxW = fe_face_values.JxW(q);
1263 const auto phi_i = fe_face_values.shape_value(j, q);
1265 boundary_mass += phi_i * JxW;
1266 normal += phi_i * fe_face_values.normal_vector(q) * JxW;
1269 const auto global_index = local_dof_indices[j];
1270 const auto index = partitioner.global_to_local(global_index);
1273 if (index >= n_locally_owned_)
1277 const unsigned int row_length =
1278 sparsity_pattern_simd_.row_length(index);
1279 if (row_length == 1)
1282 Point<dim> position =
1283 discretization_->mapping().transform_unit_to_real_cell(
1284 cell, support_points[j]);
1290 preliminary_map.insert(
1291 {index, {normal, boundary_mass, boundary_mass, id, position}});
1307 std::multimap<unsigned int, BoundaryData> filtered_map;
1308 std::set<dealii::types::global_dof_index> boundary_dofs;
1309 for (
auto entry : preliminary_map) {
1310 bool inserted =
false;
1311 const auto range = filtered_map.equal_range(entry.first);
1312 for (
auto it = range.first; it != range.second; ++it) {
1317 new_point] = entry.second;
1318 auto &[normal, normal_mass, boundary_mass, id, point] = it->second;
1323 Assert(point.distance(new_point) < 1.0e-14, dealii::ExcInternalError());
1325 if (normal * new_normal / normal.norm() / new_normal.norm() > 0.50) {
1330 normal += new_normal;
1331 boundary_mass += new_boundary_mass;
1335 }
else if constexpr (dim == 2) {
1351 filtered_map.insert(entry);
1358 BoundaryMap boundary_map;
1360 std::begin(filtered_map),
1361 std::end(filtered_map),
1362 std::back_inserter(boundary_map),
1364 auto index = it.first;
1365 const auto &[normal, normal_mass, boundary_mass, id, point] =
1368 const auto new_normal_mass =
1369 normal.norm() + std::numeric_limits<Number>::epsilon();
1370 const auto new_normal = normal / new_normal_mass;
1372 return {index, new_normal, new_normal_mass, boundary_mass, id, point};
1375 return boundary_map;
1379 template <
int dim,
typename Number>
1380 template <
typename ITERATOR1,
typename ITERATOR2>
1382 const ITERATOR1 &begin,
1383 const ITERATOR2 &end,
1384 const Utilities::MPI::Partitioner &partitioner)
const
1385 -> CouplingBoundaryPairs
1388 std::cout <<
"OfflineData<dim, Number>::collect_coupling_boundary_pairs()"
1398 std::set<unsigned int> locally_relevant_boundary_indices;
1400 std::vector<dealii::types::global_dof_index> local_dof_indices;
1402 const auto &finite_element = discretization_->finite_element();
1403 const unsigned int dofs_per_cell = finite_element.dofs_per_cell;
1404 local_dof_indices.resize(dofs_per_cell);
1406 for (
auto cell = begin; cell != end; ++cell) {
1409 if (cell->is_artificial())
1412 cell->get_active_or_mg_dof_indices(local_dof_indices);
1414 for (
auto f : cell->face_indices()) {
1415 const auto face = cell->face(f);
1416 const auto id = face->boundary_id();
1418 if (!face->at_boundary())
1425 for (
unsigned int j = 0; j < dofs_per_cell; ++j) {
1427 if (!finite_element.has_support_on_face(j, f))
1430 const auto global_index = local_dof_indices[j];
1431 const auto index = partitioner.global_to_local(global_index);
1434 if (index >= n_locally_relevant_)
1437 locally_relevant_boundary_indices.insert(index);
1446 CouplingBoundaryPairs result;
1448 for (
const auto i : locally_relevant_boundary_indices) {
1451 if (i >= n_locally_owned_)
1454 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1457 if (row_length == 1)
1460 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1461 constexpr auto simd_length = VectorizedArray<Number>::size();
1463 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1464 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1467 if (locally_relevant_boundary_indices.count(j) != 0) {
1468 result.push_back({i, col_idx, j});
std::tuple< unsigned int, dealii::Tensor< 1, dim, Number >, Number, Number, dealii::types::boundary_id, dealii::Point< dim > > BoundaryDescription
OfflineData(const MPIEnsemble &mpi_ensemble, const Discretization< dim > &discretization, const std::string &subsection="/OfflineData")
void make_extended_sparsity_pattern(const dealii::DoFHandler< dim > &dof_handler, SPARSITY &dsp, const dealii::AffineConstraints< Number > &affine_constraints, bool keep_constrained)
void transform_to_local_range(const dealii::Utilities::MPI::Partitioner &partitioner, dealii::AffineConstraints< Number > &affine_constraints)
unsigned int internal_range(dealii::DoFHandler< dim > &dof_handler, const dealii::DynamicSparsityPattern &sparsity, const std::size_t group_size)
unsigned int export_indices_first(dealii::DoFHandler< dim > &dof_handler, const MPI_Comm &mpi_communicator, const unsigned int n_locally_internal, const std::size_t group_size)
unsigned int inconsistent_strides_last(dealii::DoFHandler< dim > &dof_handler, const dealii::DynamicSparsityPattern &sparsity, const unsigned int n_locally_internal, const std::size_t group_size)
void make_extended_sparsity_pattern_dg(const dealii::DoFHandler< dim > &dof_handler, SPARSITY &dsp, const dealii::AffineConstraints< Number > &affine_constraints, bool keep_constrained)
std::shared_ptr< const dealii::Utilities::MPI::Partitioner > create_vector_partitioner(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &scalar_partitioner, const unsigned int n_components)
T pow(const T x, const T b)
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
DEAL_II_ALWAYS_INLINE FT add(const FT &flux_left_ij, const FT &flux_right_ij)