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#if DEAL_II_VERSION_GTE(9, 6, 0)
75 const auto locally_relevant =
76 DoFTools::extract_locally_relevant_dofs(dof_handler);
78 IndexSet locally_relevant;
79 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
82#if DEAL_II_VERSION_GTE(9, 6, 0)
83 affine_constraints_.reinit(locally_owned, locally_relevant);
85 affine_constraints_.reinit(locally_relevant);
87 DoFTools::make_hanging_node_constraints(dof_handler, affine_constraints_);
89#ifndef DEAL_II_WITH_TRILINOS
90 AssertThrow(affine_constraints_.n_constraints() == 0,
91 ExcMessage(
"ryujin was built without Trilinos support - no "
92 "hanging node support available"));
100 const auto &periodic_faces =
101 discretization_->triangulation().get_periodic_face_map();
103 for (
const auto &[left, value] : periodic_faces) {
104 const auto &[right, orientation] = value;
106 typename DoFHandler<dim>::cell_iterator dof_cell_left(
107 &left.first->get_triangulation(),
112 typename DoFHandler<dim>::cell_iterator dof_cell_right(
113 &right.first->get_triangulation(),
114 right.first->level(),
115 right.first->index(),
118 if constexpr (std::is_same_v<Number, double>) {
119 DoFTools::make_periodicity_constraints(
120 dof_cell_left->face(left.second),
121 dof_cell_right->face(right.second),
124#
if DEAL_II_VERSION_GTE(9, 6, 0)
132 AssertThrow(
false, dealii::ExcNotImplemented());
137 affine_constraints_.close();
142 const std::vector<IndexSet> &locally_owned_dofs =
143 Utilities::MPI::all_gather(mpi_ensemble_.ensemble_communicator(),
144 dof_handler.locally_owned_dofs());
145 const IndexSet locally_active =
146 dealii::DoFTools::extract_locally_active_dofs(dof_handler);
147 Assert(affine_constraints_.is_consistent_in_parallel(
150 mpi_ensemble_.ensemble_communicator(),
156 sparsity_pattern_.reinit(
157 dof_handler.n_dofs(), dof_handler.n_dofs(), locally_relevant);
159 if (discretization_->have_discontinuous_ansatz()) {
164 dof_handler, sparsity_pattern_, affine_constraints_,
false);
169#ifdef DEAL_II_WITH_TRILINOS
170 DoFTools::make_sparsity_pattern(
171 dof_handler, sparsity_pattern_, affine_constraints_,
false);
180 dof_handler, sparsity_pattern_, affine_constraints_,
false);
190 SparsityTools::distribute_sparsity_pattern(
193 mpi_ensemble_.ensemble_communicator(),
198 template <
int dim,
typename Number>
199 void OfflineData<dim, Number>::setup(
const unsigned int problem_dimension,
200 const unsigned int n_precomputed_values)
203 std::cout <<
"OfflineData<dim, Number>::setup()" << std::endl;
210 const auto &triangulation = discretization_->triangulation();
212 dof_handler_ = std::make_unique<dealii::DoFHandler<dim>>(triangulation);
213 auto &dof_handler = *dof_handler_;
215 dof_handler.distribute_dofs(discretization_->finite_element());
217 n_locally_owned_ = dof_handler.locally_owned_dofs().n_elements();
223 DoFRenumbering::Cuthill_McKee(dof_handler);
234 mpi_ensemble_.ensemble_communicator(),
245 create_constraints_and_sparsity_pattern();
247 dof_handler, sparsity_pattern_, VectorizedArray<Number>::size());
262 mpi_ensemble_.ensemble_communicator(),
264 VectorizedArray<Number>::size());
270 const auto consistent_stride_range [[maybe_unused]] = [&]() {
271 constexpr auto group_size = VectorizedArray<Number>::size();
272 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
273 const auto offset = n_locally_owned_ != 0 ? *locally_owned.begin() : 0;
275 unsigned int group_row_length = 0;
277 for (; i < n_locally_internal_; ++i) {
278 if (i % group_size == 0) {
279 group_row_length = sparsity_pattern_.row_length(offset + i);
281 if (group_row_length != sparsity_pattern_.row_length(offset + i)) {
286 return i / group_size * group_size;
292 const auto mpi_allreduce_logical_or = [&](
const bool local_value) {
293 std::function<bool(
const bool &,
const bool &)> comparator =
294 [](
const bool &left,
const bool &right) ->
bool {
295 return left || right;
297 return Utilities::MPI::all_reduce(
298 local_value, mpi_ensemble_.ensemble_communicator(), comparator);
305 create_constraints_and_sparsity_pattern();
315#if DEAL_II_VERSION_GTE(9, 5, 0)
316 if (mpi_allreduce_logical_or(affine_constraints_.n_constraints() > 0)) {
317 if (mpi_allreduce_logical_or(
318 consistent_stride_range() != n_locally_internal_)) {
328 VectorizedArray<Number>::size());
329 create_constraints_and_sparsity_pattern();
330 n_locally_internal_ = consistent_stride_range();
340 Assert(consistent_stride_range() == n_locally_internal_,
341 dealii::ExcInternalError());
347 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
348 Assert(n_locally_owned_ == locally_owned.n_elements(),
349 dealii::ExcInternalError());
351#if DEAL_II_VERSION_GTE(9, 6, 0)
352 auto locally_relevant =
353 DoFTools::extract_locally_relevant_dofs(dof_handler);
355 IndexSet locally_relevant;
356 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
360 IndexSet additional_dofs(dof_handler.n_dofs());
361 for (
auto &entry : sparsity_pattern_)
362 if (!locally_relevant.is_element(entry.column())) {
363 Assert(locally_owned.is_element(entry.row()), ExcInternalError());
364 additional_dofs.add_index(entry.column());
366 additional_dofs.compress();
367 locally_relevant.add_indices(additional_dofs);
368 locally_relevant.compress();
371 n_locally_relevant_ = locally_relevant.n_elements();
373 scalar_partitioner_ = std::make_shared<dealii::Utilities::MPI::Partitioner>(
374 locally_owned, locally_relevant, mpi_ensemble_.ensemble_communicator());
377 scalar_partitioner_, problem_dimension);
380 scalar_partitioner_, n_precomputed_values);
388 if (mpi_allreduce_logical_or(affine_constraints_.n_constraints() > 0)) {
392 n_export_indices_ = 0;
393 for (
const auto &it : scalar_partitioner_->import_indices())
394 if (it.second <= n_locally_internal_)
395 n_export_indices_ = std::max(n_export_indices_, it.second);
397 constexpr auto simd_length = VectorizedArray<Number>::size();
399 (n_export_indices_ + simd_length - 1) / simd_length * simd_length;
404 unsigned int control = 0;
405 for (
const auto &it : scalar_partitioner_->import_indices())
406 if (it.second <= n_locally_internal_)
407 control = std::max(control, it.second);
409 Assert(control <= n_export_indices_, ExcInternalError());
410 Assert(n_export_indices_ <= n_locally_internal_, ExcInternalError());
419 sparsity_pattern_simd_.reinit(
420 n_locally_internal_, sparsity_pattern_, scalar_partitioner_);
424 template <
int dim,
typename Number>
425 void OfflineData<dim, Number>::create_matrices()
428 std::cout <<
"OfflineData<dim, Number>::create_matrices()" << std::endl;
435 mass_matrix_.reinit(sparsity_pattern_simd_);
436 if (discretization_->have_discontinuous_ansatz())
437 mass_matrix_inverse_.reinit(sparsity_pattern_simd_);
439 lumped_mass_matrix_.reinit(scalar_partitioner_);
440 lumped_mass_matrix_inverse_.reinit(scalar_partitioner_);
442 cij_matrix_.reinit(sparsity_pattern_simd_);
443 if (discretization_->have_discontinuous_ansatz())
444 incidence_matrix_.reinit(sparsity_pattern_simd_);
450 auto &dof_handler = *dof_handler_;
452 measure_of_omega_ = 0.;
454#ifdef DEAL_II_WITH_TRILINOS
457 AffineConstraints<double> affine_constraints_assembly;
459#if DEAL_II_VERSION_GTE(9, 6, 0)
460 affine_constraints_assembly.reinit(affine_constraints_.get_local_lines(),
461 affine_constraints_.get_local_lines());
463 affine_constraints_assembly.reinit(affine_constraints_.get_local_lines());
465 for (
auto line : affine_constraints_.get_lines()) {
466 affine_constraints_assembly.add_line(line.index);
467 for (
auto entry : line.entries)
468 affine_constraints_assembly.add_entry(
469 line.index, entry.first, entry.second);
470 affine_constraints_assembly.set_inhomogeneity(line.index,
473 affine_constraints_assembly.close();
475 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
476 TrilinosWrappers::SparsityPattern trilinos_sparsity_pattern;
477 trilinos_sparsity_pattern.reinit(locally_owned,
479 mpi_ensemble_.ensemble_communicator());
481 TrilinosWrappers::SparseMatrix mass_matrix_tmp;
482 TrilinosWrappers::SparseMatrix mass_matrix_inverse_tmp;
483 if (discretization_->have_discontinuous_ansatz())
484 mass_matrix_inverse_tmp.reinit(trilinos_sparsity_pattern);
486 std::array<TrilinosWrappers::SparseMatrix, dim> cij_matrix_tmp;
488 mass_matrix_tmp.reinit(trilinos_sparsity_pattern);
489 for (
auto &matrix : cij_matrix_tmp)
490 matrix.reinit(trilinos_sparsity_pattern);
495 AffineConstraints<Number> affine_constraints_assembly;
496 affine_constraints_assembly.copy_from(affine_constraints_);
499 SparsityPattern sparsity_pattern_assembly;
501 DynamicSparsityPattern dsp(n_locally_relevant_, n_locally_relevant_);
502 for (
const auto &entry : sparsity_pattern_) {
503 const auto i = scalar_partitioner_->global_to_local(entry.row());
504 const auto j = scalar_partitioner_->global_to_local(entry.column());
507 sparsity_pattern_assembly.copy_from(dsp);
510 dealii::SparseMatrix<Number> mass_matrix_tmp;
511 dealii::SparseMatrix<Number> mass_matrix_inverse_tmp;
512 if (discretization_->have_discontinuous_ansatz())
513 mass_matrix_inverse_tmp.reinit(sparsity_pattern_assembly);
515 std::array<dealii::SparseMatrix<Number>, dim> cij_matrix_tmp;
517 mass_matrix_tmp.reinit(sparsity_pattern_assembly);
518 for (
auto &matrix : cij_matrix_tmp)
519 matrix.reinit(sparsity_pattern_assembly);
527 const auto local_assemble_system = [&](
const auto &cell,
532 auto &is_locally_owned = copy.is_locally_owned_;
533 auto &local_dof_indices = copy.local_dof_indices_;
534 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
536 auto &cell_mass_matrix = copy.cell_mass_matrix_;
537 auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
538 auto &cell_cij_matrix = copy.cell_cij_matrix_;
539 auto &interface_cij_matrix = copy.interface_cij_matrix_;
540 auto &cell_measure = copy.cell_measure_;
542 auto &hp_fe_values = scratch.hp_fe_values_;
543 auto &hp_fe_face_values = scratch.hp_fe_face_values_;
544 auto &hp_fe_neighbor_face_values = scratch.hp_fe_neighbor_face_values_;
546#ifdef DEAL_II_WITH_TRILINOS
547 is_locally_owned = cell->is_locally_owned();
555 is_locally_owned = !cell->is_artificial();
557 if (!is_locally_owned)
560 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
562 cell_mass_matrix.reinit(dofs_per_cell, dofs_per_cell);
563 for (
auto &matrix : cell_cij_matrix)
564 matrix.reinit(dofs_per_cell, dofs_per_cell);
565 if (discretization_->have_discontinuous_ansatz()) {
566 cell_mass_matrix_inverse.reinit(dofs_per_cell, dofs_per_cell);
567 for (
auto &it : interface_cij_matrix)
568 for (
auto &matrix : it)
569 matrix.reinit(dofs_per_cell, dofs_per_cell);
572 hp_fe_values.reinit(cell);
573 const auto &fe_values = hp_fe_values.get_present_fe_values();
575 local_dof_indices.resize(dofs_per_cell);
576 cell->get_dof_indices(local_dof_indices);
580 cell_mass_matrix = 0.;
581 for (
auto &matrix : cell_cij_matrix)
583 if (discretization_->have_discontinuous_ansatz()) {
584 cell_mass_matrix_inverse = 0.;
585 for (
auto &it : interface_cij_matrix)
586 for (
auto &matrix : it)
591 for (
unsigned int q : fe_values.quadrature_point_indices()) {
592 const auto JxW = fe_values.JxW(q);
594 if (cell->is_locally_owned())
595 cell_measure += Number(JxW);
597 for (
unsigned int j : fe_values.dof_indices()) {
598 const auto value_JxW = fe_values.shape_value(j, q) * JxW;
599 const auto grad_JxW = fe_values.shape_grad(j, q) * JxW;
601 for (
unsigned int i : fe_values.dof_indices()) {
602 const auto value = fe_values.shape_value(i, q);
604 cell_mass_matrix(i, j) += Number(value * value_JxW);
605 for (
unsigned int d = 0; d < dim; ++d)
606 cell_cij_matrix[d](i, j) += Number((value * grad_JxW)[d]);
616 if (!discretization_->have_discontinuous_ansatz())
619 for (
const auto f_index : cell->face_indices()) {
620 const auto &face = cell->face(f_index);
623 const bool has_neighbor =
624 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
628 neighbor_local_dof_indices[f_index].resize(0);
633 const auto neighbor_cell = cell->neighbor_or_periodic_neighbor(f_index);
634 if (neighbor_cell->is_artificial()) {
637 neighbor_local_dof_indices[f_index].resize(0);
641 hp_fe_face_values.reinit(cell, f_index);
642 const auto &fe_face_values = hp_fe_face_values.get_present_fe_values();
646 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
647 const auto JxW = fe_face_values.JxW(q);
648 const auto &normal = fe_face_values.get_normal_vectors()[q];
650 for (
unsigned int j : fe_face_values.dof_indices()) {
651 const auto value_JxW = fe_face_values.shape_value(j, q) * JxW;
653 for (
unsigned int i : fe_face_values.dof_indices()) {
654 const auto value = fe_face_values.shape_value(i, q);
656 for (
unsigned int d = 0; d < dim; ++d)
657 cell_cij_matrix[d](i, j) -=
658 Number(0.5 * normal[d] * value * value_JxW);
665 const unsigned int f_index_neighbor =
666 cell->has_periodic_neighbor(f_index)
667 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
668 : cell->neighbor_of_neighbor(f_index);
670 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
671 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
673 hp_fe_neighbor_face_values.reinit(neighbor_cell, f_index_neighbor);
674 const auto &fe_neighbor_face_values =
675 hp_fe_neighbor_face_values.get_present_fe_values();
677 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
678 const auto JxW = fe_face_values.JxW(q);
679 const auto &normal = fe_face_values.get_normal_vectors()[q];
682 for (
unsigned int j : fe_face_values.dof_indices()) {
683 const auto value_JxW =
684 fe_neighbor_face_values.shape_value(j, q) * JxW;
686 for (
unsigned int i : fe_face_values.dof_indices()) {
687 const auto value = fe_face_values.shape_value(i, q);
689 for (
unsigned int d = 0; d < dim; ++d)
690 interface_cij_matrix[f_index][d](i, j) +=
691 Number(0.5 * normal[d] * value * value_JxW);
701 if (discretization_->have_discontinuous_ansatz()) {
703 cell_mass_matrix_inverse.invert(cell_mass_matrix);
707 const auto copy_local_to_global = [&](
const auto ©) {
708 const auto &is_locally_owned = copy.is_locally_owned_;
709#ifdef DEAL_II_WITH_TRILINOS
710 const auto &local_dof_indices = copy.local_dof_indices_;
711 const auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
718 auto local_dof_indices = copy.local_dof_indices_;
719 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
721 const auto &cell_mass_matrix = copy.cell_mass_matrix_;
722 const auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
723 const auto &cell_cij_matrix = copy.cell_cij_matrix_;
724 const auto &interface_cij_matrix = copy.interface_cij_matrix_;
725 const auto &cell_measure = copy.cell_measure_;
727 if (!is_locally_owned)
730#ifndef DEAL_II_WITH_TRILINOS
732 for (
auto &indices : neighbor_local_dof_indices)
736 affine_constraints_assembly.distribute_local_to_global(
737 cell_mass_matrix, local_dof_indices, mass_matrix_tmp);
739 for (
int k = 0; k < dim; ++k) {
740 affine_constraints_assembly.distribute_local_to_global(
741 cell_cij_matrix[k], local_dof_indices, cij_matrix_tmp[k]);
743 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
744 if (neighbor_local_dof_indices[f_index].size() != 0) {
745 affine_constraints_assembly.distribute_local_to_global(
746 interface_cij_matrix[f_index][k],
748 neighbor_local_dof_indices[f_index],
754 if (discretization_->have_discontinuous_ansatz())
755 affine_constraints_assembly.distribute_local_to_global(
756 cell_mass_matrix_inverse,
758 mass_matrix_inverse_tmp);
760 measure_of_omega_ += cell_measure;
763 WorkStream::run(dof_handler.begin_active(),
765 local_assemble_system,
766 copy_local_to_global,
767 AssemblyScratchData<dim>(*discretization_),
768#ifdef DEAL_II_WITH_TRILINOS
769 AssemblyCopyData<dim, double>());
771 AssemblyCopyData<dim, Number>());
774#ifdef DEAL_II_WITH_TRILINOS
776 for (
auto &it : cij_matrix_tmp)
779 mass_matrix_.read_in(mass_matrix_tmp,
false);
780 if (discretization_->have_discontinuous_ansatz())
781 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
false);
782 cij_matrix_.read_in(cij_matrix_tmp,
false);
784 mass_matrix_.read_in(mass_matrix_tmp,
true);
785 if (discretization_->have_discontinuous_ansatz())
786 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
true);
787 cij_matrix_.read_in(cij_matrix_tmp,
true);
790 mass_matrix_.update_ghost_rows();
791 if (discretization_->have_discontinuous_ansatz())
792 mass_matrix_inverse_.update_ghost_rows();
793 cij_matrix_.update_ghost_rows();
795 measure_of_omega_ = Utilities::MPI::sum(
796 measure_of_omega_, mpi_ensemble_.ensemble_communicator());
803#ifdef DEAL_II_WITH_TRILINOS
806 affine_constraints_assembly.set_zero(one);
808 ScalarVector local_lumped_mass_matrix(scalar_partitioner_);
809 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
812 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
814 lumped_mass_matrix_.local_element(i) =
815 local_lumped_mass_matrix.local_element(i);
816 lumped_mass_matrix_inverse_.local_element(i) =
817 1. / lumped_mass_matrix_.local_element(i);
819 lumped_mass_matrix_.update_ghost_values();
820 lumped_mass_matrix_inverse_.update_ghost_values();
824 dealii::Vector<Number> one(mass_matrix_tmp.m());
826 affine_constraints_assembly.set_zero(one);
828 dealii::Vector<Number> local_lumped_mass_matrix(mass_matrix_tmp.m());
829 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
831 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
833 lumped_mass_matrix_.local_element(i) = local_lumped_mass_matrix(i);
834 lumped_mass_matrix_inverse_.local_element(i) =
835 1. / lumped_mass_matrix_.local_element(i);
837 lumped_mass_matrix_.update_ghost_values();
838 lumped_mass_matrix_inverse_.update_ghost_values();
846 if (discretization_->have_discontinuous_ansatz()) {
847#ifdef DEAL_II_WITH_TRILINOS
848 TrilinosWrappers::SparseMatrix incidence_matrix_tmp;
849 incidence_matrix_tmp.reinit(trilinos_sparsity_pattern);
851 dealii::SparseMatrix<Number> incidence_matrix_tmp;
852 incidence_matrix_tmp.reinit(sparsity_pattern_assembly);
856 const auto local_assemble_system = [&](
const auto &cell,
861 auto &is_locally_owned = copy.is_locally_owned_;
862 auto &local_dof_indices = copy.local_dof_indices_;
863 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
864 auto &interface_incidence_matrix = copy.interface_incidence_matrix_;
865 auto &hp_fe_face_values_nodal = scratch.hp_fe_face_values_nodal_;
866 auto &hp_fe_neighbor_face_values_nodal =
867 scratch.hp_fe_neighbor_face_values_nodal_;
869#ifdef DEAL_II_WITH_TRILINOS
870 is_locally_owned = cell->is_locally_owned();
872 is_locally_owned = !cell->is_artificial();
874 if (!is_locally_owned)
877 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
879 for (
auto &matrix : interface_incidence_matrix)
880 matrix.reinit(dofs_per_cell, dofs_per_cell);
882 local_dof_indices.resize(dofs_per_cell);
883 cell->get_dof_indices(local_dof_indices);
886 for (
auto &matrix : interface_incidence_matrix)
889 for (
const auto f_index : cell->face_indices()) {
890 const auto &face = cell->face(f_index);
893 const bool has_neighbor =
894 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
898 neighbor_local_dof_indices[f_index].resize(0);
903 const auto neighbor_cell =
904 cell->neighbor_or_periodic_neighbor(f_index);
905 if (neighbor_cell->is_artificial()) {
908 neighbor_local_dof_indices[f_index].resize(0);
912 const unsigned int f_index_neighbor =
913 cell->has_periodic_neighbor(f_index)
914 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
915 : cell->neighbor_of_neighbor(f_index);
917 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
918 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
920 hp_fe_face_values_nodal.reinit(cell, f_index);
921 const auto &fe_face_values_nodal =
922 hp_fe_face_values_nodal.get_present_fe_values();
923 hp_fe_neighbor_face_values_nodal.reinit(neighbor_cell,
925 const auto &fe_neighbor_face_values_nodal =
926 hp_fe_neighbor_face_values_nodal.get_present_fe_values();
930 for (
unsigned int q :
931 fe_face_values_nodal.quadrature_point_indices()) {
933 for (
unsigned int j : fe_neighbor_face_values_nodal.dof_indices()) {
934 const auto v_j = fe_neighbor_face_values_nodal.shape_value(j, q);
935 for (
unsigned int i : fe_face_values_nodal.dof_indices()) {
936 const auto v_i = fe_face_values_nodal.shape_value(i, q);
937 constexpr auto eps = std::numeric_limits<Number>::epsilon();
938 if (std::abs(v_i * v_j) > 100. * eps) {
939 const auto &ansatz = discretization_->ansatz();
941 const auto glob_i = local_dof_indices[i];
942 const auto glob_j = neighbor_local_dof_indices[f_index][j];
943 const auto m_i = lumped_mass_matrix_[glob_i];
944 const auto m_j = lumped_mass_matrix_[glob_j];
946 Number(0.5) * (m_i + m_j) / measure_of_omega_;
957 r_ij =
std::pow(hd_ij, incidence_relaxation_even_ / dim);
964 r_ij =
std::pow(hd_ij, incidence_relaxation_odd_ / dim);
967 interface_incidence_matrix[f_index](i, j) += r_ij;
975 const auto copy_local_to_global = [&](
const auto ©) {
976 const auto &is_locally_owned = copy.is_locally_owned_;
977#ifdef DEAL_II_WITH_TRILINOS
978 const auto &local_dof_indices = copy.local_dof_indices_;
979 const auto &neighbor_local_dof_indices =
980 copy.neighbor_local_dof_indices_;
987 auto local_dof_indices = copy.local_dof_indices_;
988 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
990 const auto &interface_incidence_matrix =
991 copy.interface_incidence_matrix_;
993 if (!is_locally_owned)
996#ifndef DEAL_II_WITH_TRILINOS
998 for (
auto &indices : neighbor_local_dof_indices)
1002 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
1003 if (neighbor_local_dof_indices[f_index].size() != 0) {
1004 affine_constraints_assembly.distribute_local_to_global(
1005 interface_incidence_matrix[f_index],
1007 neighbor_local_dof_indices[f_index],
1008 incidence_matrix_tmp);
1013 WorkStream::run(dof_handler.begin_active(),
1015 local_assemble_system,
1016 copy_local_to_global,
1017 AssemblyScratchData<dim>(*discretization_),
1018#ifdef DEAL_II_WITH_TRILINOS
1019 AssemblyCopyData<dim, double>());
1021 AssemblyCopyData<dim, Number>());
1024#ifdef DEAL_II_WITH_TRILINOS
1025 incidence_matrix_.read_in(incidence_matrix_tmp,
false);
1027 incidence_matrix_.read_in(incidence_matrix_tmp,
true);
1029 incidence_matrix_.update_ghost_rows();
1036 boundary_map_ = construct_boundary_map(
1037 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1039 coupling_boundary_pairs_ = collect_coupling_boundary_pairs(
1040 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1043#ifdef DEBUG_SYMMETRY_CHECK
1048 double total_mass = 0.;
1049 for (
unsigned int i = 0; i < n_locally_owned_; ++i)
1050 total_mass += lumped_mass_matrix_.local_element(i);
1052 Utilities::MPI::sum(total_mass, mpi_ensemble_.ensemble_communicator());
1054 Assert(std::abs(measure_of_omega_ - total_mass) <
1055 1.e-12 * measure_of_omega_,
1057 "Total mass differs from the measure of the domain."));
1063 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1065 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1066 if (row_length == 1)
1070 mass_matrix_.get_entry(i, 0) - lumped_mass_matrix_.local_element(i);
1073 constexpr auto simd_length = VectorizedArray<Number>::size();
1074 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1075 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1076 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1078 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1080 const auto m_ij = mass_matrix_.get_entry(i, col_idx);
1081 if (discretization_->have_discontinuous_ansatz()) {
1084 Assert(std::abs(m_ij) > -1.e-12, dealii::ExcInternalError());
1086 Assert(std::abs(m_ij) > 1.e-12, dealii::ExcInternalError());
1090 const auto m_ji = mass_matrix_.get_transposed_entry(i, col_idx);
1091 if (std::abs(m_ij - m_ji) >= 1.e-12) {
1093 std::stringstream ss;
1094 ss <<
"m_ij matrix is not symmetric: " << m_ij <<
" <-> " << m_ji;
1095 Assert(
false, dealii::ExcMessage(ss.str()));
1099 Assert(std::abs(sum) < 1.e-12, dealii::ExcInternalError());
1106 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1108 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1109 if (row_length == 1)
1112 auto sum = cij_matrix_.get_tensor(i, 0);
1115 constexpr auto simd_length = VectorizedArray<Number>::size();
1116 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1117 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1118 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1120 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1122 const auto c_ij = cij_matrix_.get_tensor(i, col_idx);
1123 Assert(c_ij.norm() > 1.e-12, dealii::ExcInternalError());
1126 const auto c_ji = cij_matrix_.get_transposed_tensor(i, col_idx);
1127 if ((c_ij + c_ji).norm() >= 1.e-12) {
1131 CouplingDescription coupling{i, col_idx, j};
1132 const auto it = std::find(coupling_boundary_pairs_.begin(),
1133 coupling_boundary_pairs_.end(),
1135 if (it == coupling_boundary_pairs_.end()) {
1136 std::stringstream ss;
1137 ss <<
"c_ij matrix is not anti-symmetric: " << c_ij <<
" <-> "
1139 Assert(
false, dealii::ExcMessage(ss.str()));
1144 Assert(sum.norm() < 1.e-12, dealii::ExcInternalError());
1150 template <
int dim,
typename Number>
1151 void OfflineData<dim, Number>::create_multigrid_data()
1154 std::cout <<
"OfflineData<dim, Number>::create_multigrid_data()"
1158 auto &dof_handler = *dof_handler_;
1160 Assert(!dof_handler.has_hp_capabilities(), dealii::ExcInternalError());
1162 dof_handler.distribute_mg_dofs();
1164 const auto n_levels = dof_handler.get_triangulation().n_global_levels();
1166 AffineConstraints<float> level_constraints;
1169 level_boundary_map_.resize(n_levels);
1170 level_lumped_mass_matrix_.resize(n_levels);
1172 for (
unsigned int level = 0; level < n_levels; ++level) {
1175#if DEAL_II_VERSION_GTE(9, 6, 0)
1176 const auto relevant_dofs =
1177 dealii::DoFTools::extract_locally_relevant_level_dofs(dof_handler,
1180 IndexSet relevant_dofs;
1181 dealii::DoFTools::extract_locally_relevant_level_dofs(
1182 dof_handler, level, relevant_dofs);
1185 const auto partitioner = std::make_shared<Utilities::MPI::Partitioner>(
1186 dof_handler.locally_owned_mg_dofs(level),
1188 mpi_ensemble_.ensemble_communicator());
1189 level_lumped_mass_matrix_[level].reinit(partitioner);
1190 std::vector<types::global_dof_index> dof_indices(
1191 dof_handler.get_fe().dofs_per_cell);
1192 dealii::Vector<Number> mass_values(dof_handler.get_fe().dofs_per_cell);
1193 dealii::hp::FEValues<dim> hp_fe_values(discretization_->mapping(),
1194 discretization_->finite_element(),
1195 discretization_->quadrature(),
1196 update_values | update_JxW_values);
1197 for (
const auto &cell : dof_handler.cell_iterators_on_level(level))
1200 if (cell->is_locally_owned_on_level()) {
1201 hp_fe_values.reinit(cell);
1202 const auto &fe_values = hp_fe_values.get_present_fe_values();
1203 for (
unsigned int i = 0; i < mass_values.size(); ++i) {
1205 for (
unsigned int q = 0; q < fe_values.n_quadrature_points; ++q)
1206 sum += fe_values.shape_value(i, q) * fe_values.JxW(q);
1207 mass_values(i) = sum;
1209 cell->get_mg_dof_indices(dof_indices);
1210 level_constraints.distribute_local_to_global(
1211 mass_values, dof_indices, level_lumped_mass_matrix_[level]);
1217 level_boundary_map_[level] = construct_boundary_map(
1218 dof_handler.begin_mg(level), dof_handler.end_mg(level), *partitioner);
1223 template <
int dim,
typename Number>
1224 template <
typename ITERATOR1,
typename ITERATOR2>
1226 const ITERATOR1 &begin,
1227 const ITERATOR2 &end,
1228 const Utilities::MPI::Partitioner &partitioner)
const -> BoundaryMap
1231 std::cout <<
"OfflineData<dim, Number>::construct_boundary_map()"
1239 using BoundaryData = std::tuple<dealii::Tensor<1, dim, Number> ,
1242 dealii::types::boundary_id ,
1243 dealii::Point<dim>> ;
1244 std::multimap<unsigned int, BoundaryData> preliminary_map;
1246 std::vector<dealii::types::global_dof_index> local_dof_indices;
1248 dealii::hp::FEFaceValues<dim> hp_fe_face_values(
1249 discretization_->mapping(),
1250 discretization_->finite_element(),
1251 discretization_->face_quadrature(),
1252 dealii::update_normal_vectors | dealii::update_values |
1253 dealii::update_JxW_values);
1255 for (
auto cell = begin; cell != end; ++cell) {
1263 if ((cell->is_active() && cell->is_artificial()) ||
1264 (!cell->is_active() && cell->is_artificial_on_level()))
1267 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1268 const auto &support_points = cell->get_fe().get_unit_support_points();
1270 local_dof_indices.resize(dofs_per_cell);
1271 cell->get_active_or_mg_dof_indices(local_dof_indices);
1273 for (
auto f : cell->face_indices()) {
1274 const auto face = cell->face(f);
1275 const auto id = face->boundary_id();
1279 if (!face->at_boundary())
1290 hp_fe_face_values.reinit(cell, f);
1291 const auto &fe_face_values = hp_fe_face_values.get_present_fe_values();
1292 const auto &mapping =
1293 hp_fe_face_values.get_mapping_collection()[cell->active_fe_index()];
1295 for (
unsigned int j : fe_face_values.dof_indices()) {
1296 if (!cell->get_fe().has_support_on_face(j, f))
1299 Number boundary_mass = 0.;
1300 dealii::Tensor<1, dim, Number> normal;
1302 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
1303 const auto JxW = fe_face_values.JxW(q);
1304 const auto phi_i = fe_face_values.shape_value(j, q);
1306 boundary_mass += phi_i * JxW;
1307 normal += phi_i * fe_face_values.normal_vector(q) * JxW;
1310 const auto global_index = local_dof_indices[j];
1311 const auto index = partitioner.global_to_local(global_index);
1314 if (index >= n_locally_owned_)
1318 const unsigned int row_length =
1319 sparsity_pattern_simd_.row_length(index);
1320 if (row_length == 1)
1323 Point<dim> position =
1324 mapping.transform_unit_to_real_cell(cell, support_points[j]);
1330 preliminary_map.insert(
1331 {index, {normal, boundary_mass, boundary_mass, id, position}});
1347 std::multimap<unsigned int, BoundaryData> filtered_map;
1348 std::set<dealii::types::global_dof_index> boundary_dofs;
1349 for (
auto entry : preliminary_map) {
1350 bool inserted =
false;
1351 const auto range = filtered_map.equal_range(entry.first);
1352 for (
auto it = range.first; it != range.second; ++it) {
1357 new_point] = entry.second;
1358 auto &[normal, normal_mass, boundary_mass, id, point] = it->second;
1363 Assert(point.distance(new_point) < 1.0e-14, dealii::ExcInternalError());
1365 if (normal * new_normal / normal.norm() / new_normal.norm() > 0.50) {
1370 normal += new_normal;
1371 boundary_mass += new_boundary_mass;
1375 }
else if constexpr (dim == 2) {
1391 filtered_map.insert(entry);
1398 BoundaryMap boundary_map;
1400 std::begin(filtered_map),
1401 std::end(filtered_map),
1402 std::back_inserter(boundary_map),
1404 auto index = it.first;
1405 const auto &[normal, normal_mass, boundary_mass, id, point] =
1408 const auto new_normal_mass =
1409 normal.norm() + std::numeric_limits<Number>::epsilon();
1410 const auto new_normal = normal / new_normal_mass;
1412 return {index, new_normal, new_normal_mass, boundary_mass, id, point};
1415 return boundary_map;
1419 template <
int dim,
typename Number>
1420 template <
typename ITERATOR1,
typename ITERATOR2>
1422 const ITERATOR1 &begin,
1423 const ITERATOR2 &end,
1424 const Utilities::MPI::Partitioner &partitioner)
const
1425 -> CouplingBoundaryPairs
1428 std::cout <<
"OfflineData<dim, Number>::collect_coupling_boundary_pairs()"
1438 std::set<unsigned int> locally_relevant_boundary_indices;
1440 std::vector<dealii::types::global_dof_index> local_dof_indices;
1442 for (
auto cell = begin; cell != end; ++cell) {
1445 if (cell->is_artificial())
1448 const auto &finite_element = cell->get_fe();
1449 const unsigned int dofs_per_cell = finite_element.dofs_per_cell;
1450 local_dof_indices.resize(dofs_per_cell);
1451 cell->get_active_or_mg_dof_indices(local_dof_indices);
1453 for (
auto f : cell->face_indices()) {
1454 const auto face = cell->face(f);
1455 const auto id = face->boundary_id();
1457 if (!face->at_boundary())
1464 for (
unsigned int j = 0; j < dofs_per_cell; ++j) {
1466 if (!cell->get_fe().has_support_on_face(j, f))
1469 const auto global_index = local_dof_indices[j];
1470 const auto index = partitioner.global_to_local(global_index);
1473 if (index >= n_locally_relevant_)
1476 locally_relevant_boundary_indices.insert(index);
1485 CouplingBoundaryPairs result;
1487 for (
const auto i : locally_relevant_boundary_indices) {
1490 if (i >= n_locally_owned_)
1493 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1496 if (row_length == 1)
1499 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1500 constexpr auto simd_length = VectorizedArray<Number>::size();
1502 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1503 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1506 if (locally_relevant_boundary_indices.count(j) != 0) {
1507 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)