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_);
412 mass_matrix_.reinit(sparsity_pattern_simd_);
413 if (discretization_->have_discontinuous_ansatz())
414 mass_matrix_inverse_.reinit(sparsity_pattern_simd_);
416 lumped_mass_matrix_.reinit(scalar_partitioner_);
417 lumped_mass_matrix_inverse_.reinit(scalar_partitioner_);
419 cij_matrix_.reinit(sparsity_pattern_simd_);
420 if (discretization_->have_discontinuous_ansatz())
421 incidence_matrix_.reinit(sparsity_pattern_simd_);
425 template <
int dim,
typename Number>
426 void OfflineData<dim, Number>::assemble()
429 std::cout <<
"OfflineData<dim, Number>::assemble()" << std::endl;
432 auto &dof_handler = *dof_handler_;
434 measure_of_omega_ = 0.;
436#ifdef DEAL_II_WITH_TRILINOS
439 AffineConstraints<double> affine_constraints_assembly;
441 affine_constraints_assembly.reinit(affine_constraints_.get_local_lines());
442 for (
auto line : affine_constraints_.get_lines()) {
443 affine_constraints_assembly.add_line(line.index);
444 for (
auto entry : line.entries)
445 affine_constraints_assembly.add_entry(
446 line.index, entry.first, entry.second);
447 affine_constraints_assembly.set_inhomogeneity(line.index,
450 affine_constraints_assembly.close();
452 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
453 TrilinosWrappers::SparsityPattern trilinos_sparsity_pattern;
454 trilinos_sparsity_pattern.reinit(locally_owned,
456 mpi_ensemble_.ensemble_communicator());
458 TrilinosWrappers::SparseMatrix mass_matrix_tmp;
459 TrilinosWrappers::SparseMatrix mass_matrix_inverse_tmp;
460 if (discretization_->have_discontinuous_ansatz())
461 mass_matrix_inverse_tmp.reinit(trilinos_sparsity_pattern);
463 std::array<TrilinosWrappers::SparseMatrix, dim> cij_matrix_tmp;
465 mass_matrix_tmp.reinit(trilinos_sparsity_pattern);
466 for (
auto &matrix : cij_matrix_tmp)
467 matrix.reinit(trilinos_sparsity_pattern);
472 AffineConstraints<Number> affine_constraints_assembly;
473 affine_constraints_assembly.copy_from(affine_constraints_);
476 SparsityPattern sparsity_pattern_assembly;
478 DynamicSparsityPattern dsp(n_locally_relevant_, n_locally_relevant_);
479 for (
const auto &entry : sparsity_pattern_) {
480 const auto i = scalar_partitioner_->global_to_local(entry.row());
481 const auto j = scalar_partitioner_->global_to_local(entry.column());
484 sparsity_pattern_assembly.copy_from(dsp);
487 dealii::SparseMatrix<Number> mass_matrix_tmp;
488 dealii::SparseMatrix<Number> mass_matrix_inverse_tmp;
489 if (discretization_->have_discontinuous_ansatz())
490 mass_matrix_inverse_tmp.reinit(sparsity_pattern_assembly);
492 std::array<dealii::SparseMatrix<Number>, dim> cij_matrix_tmp;
494 mass_matrix_tmp.reinit(sparsity_pattern_assembly);
495 for (
auto &matrix : cij_matrix_tmp)
496 matrix.reinit(sparsity_pattern_assembly);
499 const unsigned int dofs_per_cell =
500 discretization_->finite_element().dofs_per_cell;
507 const auto local_assemble_system = [&](
const auto &cell,
512 auto &is_locally_owned = copy.is_locally_owned_;
513 auto &local_dof_indices = copy.local_dof_indices_;
514 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
516 auto &cell_mass_matrix = copy.cell_mass_matrix_;
517 auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
518 auto &cell_cij_matrix = copy.cell_cij_matrix_;
519 auto &interface_cij_matrix = copy.interface_cij_matrix_;
520 auto &cell_measure = copy.cell_measure_;
522 auto &fe_values = scratch.fe_values_;
523 auto &fe_face_values = scratch.fe_face_values_;
524 auto &fe_neighbor_face_values = scratch.fe_neighbor_face_values_;
526#ifdef DEAL_II_WITH_TRILINOS
527 is_locally_owned = cell->is_locally_owned();
535 is_locally_owned = !cell->is_artificial();
537 if (!is_locally_owned)
540 cell_mass_matrix.reinit(dofs_per_cell, dofs_per_cell);
541 for (
auto &matrix : cell_cij_matrix)
542 matrix.reinit(dofs_per_cell, dofs_per_cell);
543 if (discretization_->have_discontinuous_ansatz()) {
544 cell_mass_matrix_inverse.reinit(dofs_per_cell, dofs_per_cell);
545 for (
auto &it : interface_cij_matrix)
546 for (
auto &matrix : it)
547 matrix.reinit(dofs_per_cell, dofs_per_cell);
550 fe_values.reinit(cell);
552 local_dof_indices.resize(dofs_per_cell);
553 cell->get_dof_indices(local_dof_indices);
556 cell_mass_matrix = 0.;
557 for (
auto &matrix : cell_cij_matrix)
559 if (discretization_->have_discontinuous_ansatz()) {
560 cell_mass_matrix_inverse = 0.;
561 for (
auto &it : interface_cij_matrix)
562 for (
auto &matrix : it)
567 for (
unsigned int q : fe_values.quadrature_point_indices()) {
568 const auto JxW = fe_values.JxW(q);
570 if (cell->is_locally_owned())
571 cell_measure += Number(JxW);
573 for (
unsigned int j : fe_values.dof_indices()) {
574 const auto value_JxW = fe_values.shape_value(j, q) * JxW;
575 const auto grad_JxW = fe_values.shape_grad(j, q) * JxW;
577 for (
unsigned int i : fe_values.dof_indices()) {
578 const auto value = fe_values.shape_value(i, q);
580 cell_mass_matrix(i, j) += Number(value * value_JxW);
581 for (
unsigned int d = 0; d < dim; ++d)
582 cell_cij_matrix[d](i, j) += Number((value * grad_JxW)[d]);
592 if (!discretization_->have_discontinuous_ansatz())
595 for (
const auto f_index : cell->face_indices()) {
596 const auto &face = cell->face(f_index);
599 const bool has_neighbor =
600 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
604 neighbor_local_dof_indices[f_index].resize(0);
609 const auto neighbor_cell = cell->neighbor_or_periodic_neighbor(f_index);
610 if (neighbor_cell->is_artificial()) {
613 neighbor_local_dof_indices[f_index].resize(0);
617 fe_face_values.reinit(cell, f_index);
621 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
622 const auto JxW = fe_face_values.JxW(q);
623 const auto &normal = fe_face_values.get_normal_vectors()[q];
625 for (
unsigned int j : fe_face_values.dof_indices()) {
626 const auto value_JxW = fe_face_values.shape_value(j, q) * JxW;
628 for (
unsigned int i : fe_face_values.dof_indices()) {
629 const auto value = fe_face_values.shape_value(i, q);
631 for (
unsigned int d = 0; d < dim; ++d)
632 cell_cij_matrix[d](i, j) -=
633 Number(0.5 * normal[d] * value * value_JxW);
640 const unsigned int f_index_neighbor =
641 cell->has_periodic_neighbor(f_index)
642 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
643 : cell->neighbor_of_neighbor(f_index);
645 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
646 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
648 fe_neighbor_face_values.reinit(neighbor_cell, f_index_neighbor);
650 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
651 const auto JxW = fe_face_values.JxW(q);
652 const auto &normal = fe_face_values.get_normal_vectors()[q];
655 for (
unsigned int j : fe_face_values.dof_indices()) {
656 const auto value_JxW =
657 fe_neighbor_face_values.shape_value(j, q) * JxW;
659 for (
unsigned int i : fe_face_values.dof_indices()) {
660 const auto value = fe_face_values.shape_value(i, q);
662 for (
unsigned int d = 0; d < dim; ++d)
663 interface_cij_matrix[f_index][d](i, j) +=
664 Number(0.5 * normal[d] * value * value_JxW);
674 if (discretization_->have_discontinuous_ansatz()) {
676 cell_mass_matrix_inverse.invert(cell_mass_matrix);
680 const auto copy_local_to_global = [&](
const auto ©) {
681 const auto &is_locally_owned = copy.is_locally_owned_;
682#ifdef DEAL_II_WITH_TRILINOS
683 const auto &local_dof_indices = copy.local_dof_indices_;
684 const auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
691 auto local_dof_indices = copy.local_dof_indices_;
692 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
694 const auto &cell_mass_matrix = copy.cell_mass_matrix_;
695 const auto &cell_mass_matrix_inverse = copy.cell_mass_matrix_inverse_;
696 const auto &cell_cij_matrix = copy.cell_cij_matrix_;
697 const auto &interface_cij_matrix = copy.interface_cij_matrix_;
698 const auto &cell_measure = copy.cell_measure_;
700 if (!is_locally_owned)
703#ifndef DEAL_II_WITH_TRILINOS
705 for (
auto &indices : neighbor_local_dof_indices)
709 affine_constraints_assembly.distribute_local_to_global(
710 cell_mass_matrix, local_dof_indices, mass_matrix_tmp);
712 for (
int k = 0; k < dim; ++k) {
713 affine_constraints_assembly.distribute_local_to_global(
714 cell_cij_matrix[k], local_dof_indices, cij_matrix_tmp[k]);
716 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
717 if (neighbor_local_dof_indices[f_index].size() != 0) {
718 affine_constraints_assembly.distribute_local_to_global(
719 interface_cij_matrix[f_index][k],
721 neighbor_local_dof_indices[f_index],
727 if (discretization_->have_discontinuous_ansatz())
728 affine_constraints_assembly.distribute_local_to_global(
729 cell_mass_matrix_inverse,
731 mass_matrix_inverse_tmp);
733 measure_of_omega_ += cell_measure;
736 WorkStream::run(dof_handler.begin_active(),
738 local_assemble_system,
739 copy_local_to_global,
740 AssemblyScratchData<dim>(*discretization_),
741#ifdef DEAL_II_WITH_TRILINOS
742 AssemblyCopyData<dim, double>());
744 AssemblyCopyData<dim, Number>());
747#ifdef DEAL_II_WITH_TRILINOS
749 for (
auto &it : cij_matrix_tmp)
752 mass_matrix_.read_in(mass_matrix_tmp,
false);
753 if (discretization_->have_discontinuous_ansatz())
754 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
false);
755 cij_matrix_.read_in(cij_matrix_tmp,
false);
757 mass_matrix_.read_in(mass_matrix_tmp,
true);
758 if (discretization_->have_discontinuous_ansatz())
759 mass_matrix_inverse_.read_in(mass_matrix_inverse_tmp,
true);
760 cij_matrix_.read_in(cij_matrix_tmp,
true);
763 mass_matrix_.update_ghost_rows();
764 if (discretization_->have_discontinuous_ansatz())
765 mass_matrix_inverse_.update_ghost_rows();
766 cij_matrix_.update_ghost_rows();
768 measure_of_omega_ = Utilities::MPI::sum(
769 measure_of_omega_, mpi_ensemble_.ensemble_communicator());
776#ifdef DEAL_II_WITH_TRILINOS
779 affine_constraints_assembly.set_zero(one);
781 ScalarVector local_lumped_mass_matrix(scalar_partitioner_);
782 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
785 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
787 lumped_mass_matrix_.local_element(i) =
788 local_lumped_mass_matrix.local_element(i);
789 lumped_mass_matrix_inverse_.local_element(i) =
790 1. / lumped_mass_matrix_.local_element(i);
792 lumped_mass_matrix_.update_ghost_values();
793 lumped_mass_matrix_inverse_.update_ghost_values();
797 dealii::Vector<Number> one(mass_matrix_tmp.m());
799 affine_constraints_assembly.set_zero(one);
801 dealii::Vector<Number> local_lumped_mass_matrix(mass_matrix_tmp.m());
802 mass_matrix_tmp.vmult(local_lumped_mass_matrix, one);
804 for (
unsigned int i = 0; i < scalar_partitioner_->locally_owned_size();
806 lumped_mass_matrix_.local_element(i) = local_lumped_mass_matrix(i);
807 lumped_mass_matrix_inverse_.local_element(i) =
808 1. / lumped_mass_matrix_.local_element(i);
810 lumped_mass_matrix_.update_ghost_values();
811 lumped_mass_matrix_inverse_.update_ghost_values();
819 if (discretization_->have_discontinuous_ansatz()) {
820#ifdef DEAL_II_WITH_TRILINOS
821 TrilinosWrappers::SparseMatrix incidence_matrix_tmp;
822 incidence_matrix_tmp.reinit(trilinos_sparsity_pattern);
824 dealii::SparseMatrix<Number> incidence_matrix_tmp;
825 incidence_matrix_tmp.reinit(sparsity_pattern_assembly);
829 const auto local_assemble_system = [&](
const auto &cell,
834 auto &is_locally_owned = copy.is_locally_owned_;
835 auto &local_dof_indices = copy.local_dof_indices_;
836 auto &neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
837 auto &interface_incidence_matrix = copy.interface_incidence_matrix_;
838 auto &fe_face_values_nodal = scratch.fe_face_values_nodal_;
839 auto &fe_neighbor_face_values_nodal =
840 scratch.fe_neighbor_face_values_nodal_;
842#ifdef DEAL_II_WITH_TRILINOS
843 is_locally_owned = cell->is_locally_owned();
845 is_locally_owned = !cell->is_artificial();
847 if (!is_locally_owned)
850 for (
auto &matrix : interface_incidence_matrix)
851 matrix.reinit(dofs_per_cell, dofs_per_cell);
853 local_dof_indices.resize(dofs_per_cell);
854 cell->get_dof_indices(local_dof_indices);
857 for (
auto &matrix : interface_incidence_matrix)
860 for (
const auto f_index : cell->face_indices()) {
861 const auto &face = cell->face(f_index);
864 const bool has_neighbor =
865 !face->at_boundary() || cell->has_periodic_neighbor(f_index);
869 neighbor_local_dof_indices[f_index].resize(0);
874 const auto neighbor_cell =
875 cell->neighbor_or_periodic_neighbor(f_index);
876 if (neighbor_cell->is_artificial()) {
879 neighbor_local_dof_indices[f_index].resize(0);
883 const unsigned int f_index_neighbor =
884 cell->has_periodic_neighbor(f_index)
885 ? cell->periodic_neighbor_of_periodic_neighbor(f_index)
886 : cell->neighbor_of_neighbor(f_index);
888 neighbor_local_dof_indices[f_index].resize(dofs_per_cell);
889 neighbor_cell->get_dof_indices(neighbor_local_dof_indices[f_index]);
891 fe_face_values_nodal.reinit(cell, f_index);
892 fe_neighbor_face_values_nodal.reinit(neighbor_cell, f_index_neighbor);
896 for (
unsigned int q :
897 fe_face_values_nodal.quadrature_point_indices()) {
899 for (
unsigned int j : fe_neighbor_face_values_nodal.dof_indices()) {
900 const auto v_j = fe_neighbor_face_values_nodal.shape_value(j, q);
901 for (
unsigned int i : fe_face_values_nodal.dof_indices()) {
902 const auto v_i = fe_face_values_nodal.shape_value(i, q);
903 constexpr auto eps = std::numeric_limits<Number>::epsilon();
904 if (std::abs(v_i * v_j) > 100. * eps) {
905 const auto &ansatz = discretization_->ansatz();
907 const auto glob_i = local_dof_indices[i];
908 const auto glob_j = neighbor_local_dof_indices[f_index][j];
909 const auto m_i = lumped_mass_matrix_[glob_i];
910 const auto m_j = lumped_mass_matrix_[glob_j];
912 Number(0.5) * (m_i + m_j) / measure_of_omega_;
923 r_ij =
std::pow(hd_ij, incidence_relaxation_even_ / dim);
930 r_ij =
std::pow(hd_ij, incidence_relaxation_odd_ / dim);
933 interface_incidence_matrix[f_index](i, j) += r_ij;
941 const auto copy_local_to_global = [&](
const auto ©) {
942 const auto &is_locally_owned = copy.is_locally_owned_;
943#ifdef DEAL_II_WITH_TRILINOS
944 const auto &local_dof_indices = copy.local_dof_indices_;
945 const auto &neighbor_local_dof_indices =
946 copy.neighbor_local_dof_indices_;
953 auto local_dof_indices = copy.local_dof_indices_;
954 auto neighbor_local_dof_indices = copy.neighbor_local_dof_indices_;
956 const auto &interface_incidence_matrix =
957 copy.interface_incidence_matrix_;
959 if (!is_locally_owned)
962#ifndef DEAL_II_WITH_TRILINOS
964 for (
auto &indices : neighbor_local_dof_indices)
968 for (
unsigned int f_index = 0; f_index < copy.n_faces; ++f_index) {
969 if (neighbor_local_dof_indices[f_index].size() != 0) {
970 affine_constraints_assembly.distribute_local_to_global(
971 interface_incidence_matrix[f_index],
973 neighbor_local_dof_indices[f_index],
974 incidence_matrix_tmp);
979 WorkStream::run(dof_handler.begin_active(),
981 local_assemble_system,
982 copy_local_to_global,
983 AssemblyScratchData<dim>(*discretization_),
984#ifdef DEAL_II_WITH_TRILINOS
985 AssemblyCopyData<dim, double>());
987 AssemblyCopyData<dim, Number>());
990#ifdef DEAL_II_WITH_TRILINOS
991 incidence_matrix_.read_in(incidence_matrix_tmp,
false);
993 incidence_matrix_.read_in(incidence_matrix_tmp,
true);
995 incidence_matrix_.update_ghost_rows();
1002 boundary_map_ = construct_boundary_map(
1003 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1005 coupling_boundary_pairs_ = collect_coupling_boundary_pairs(
1006 dof_handler.begin_active(), dof_handler.end(), *scalar_partitioner_);
1014 double total_mass = 0.;
1015 for (
unsigned int i = 0; i < n_locally_owned_; ++i)
1016 total_mass += lumped_mass_matrix_.local_element(i);
1018 Utilities::MPI::sum(total_mass, mpi_ensemble_.ensemble_communicator());
1020 Assert(std::abs(measure_of_omega_ - total_mass) <
1021 1.e-12 * measure_of_omega_,
1023 "Total mass differs from the measure of the domain."));
1029 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1031 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1032 if (row_length == 1)
1036 mass_matrix_.get_entry(i, 0) - lumped_mass_matrix_.local_element(i);
1039 constexpr auto simd_length = VectorizedArray<Number>::size();
1040 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1041 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1042 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1044 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1046 const auto m_ij = mass_matrix_.get_entry(i, col_idx);
1047 if (discretization_->have_discontinuous_ansatz()) {
1050 Assert(std::abs(m_ij) > -1.e-12, dealii::ExcInternalError());
1052 Assert(std::abs(m_ij) > 1.e-12, dealii::ExcInternalError());
1056 const auto m_ji = mass_matrix_.get_transposed_entry(i, col_idx);
1057 if (std::abs(m_ij - m_ji) >= 1.e-12) {
1059 std::stringstream ss;
1060 ss <<
"m_ij matrix is not symmetric: " << m_ij <<
" <-> " << m_ji;
1061 Assert(
false, dealii::ExcMessage(ss.str()));
1065 Assert(std::abs(sum) < 1.e-12, dealii::ExcInternalError());
1072 for (
unsigned int i = 0; i < n_locally_owned_; ++i) {
1074 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1075 if (row_length == 1)
1078 auto sum = cij_matrix_.get_tensor(i, 0);
1081 constexpr auto simd_length = VectorizedArray<Number>::size();
1082 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1083 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1084 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1086 Assert(j < n_locally_relevant_, dealii::ExcInternalError());
1088 const auto c_ij = cij_matrix_.get_tensor(i, col_idx);
1089 Assert(c_ij.norm() > 1.e-12, dealii::ExcInternalError());
1092 const auto c_ji = cij_matrix_.get_transposed_tensor(i, col_idx);
1093 if ((c_ij + c_ji).norm() >= 1.e-12) {
1097 CouplingDescription coupling{i, col_idx, j};
1098 const auto it = std::find(coupling_boundary_pairs_.begin(),
1099 coupling_boundary_pairs_.end(),
1101 if (it == coupling_boundary_pairs_.end()) {
1102 std::stringstream ss;
1103 ss <<
"c_ij matrix is not anti-symmetric: " << c_ij <<
" <-> "
1105 Assert(
false, dealii::ExcMessage(ss.str()));
1110 Assert(sum.norm() < 1.e-12, dealii::ExcInternalError());
1116 template <
int dim,
typename Number>
1117 void OfflineData<dim, Number>::create_multigrid_data()
1120 std::cout <<
"OfflineData<dim, Number>::create_multigrid_data()"
1124 auto &dof_handler = *dof_handler_;
1126 dof_handler.distribute_mg_dofs();
1128 const auto n_levels = dof_handler.get_triangulation().n_global_levels();
1130 AffineConstraints<float> level_constraints;
1133 level_boundary_map_.resize(n_levels);
1134 level_lumped_mass_matrix_.resize(n_levels);
1136 for (
unsigned int level = 0; level < n_levels; ++level) {
1139 IndexSet relevant_dofs;
1140 dealii::DoFTools::extract_locally_relevant_level_dofs(
1141 dof_handler, level, relevant_dofs);
1142 const auto partitioner = std::make_shared<Utilities::MPI::Partitioner>(
1143 dof_handler.locally_owned_mg_dofs(level),
1145 mpi_ensemble_.ensemble_communicator());
1146 level_lumped_mass_matrix_[level].reinit(partitioner);
1147 std::vector<types::global_dof_index> dof_indices(
1148 dof_handler.get_fe().dofs_per_cell);
1149 dealii::Vector<Number> mass_values(dof_handler.get_fe().dofs_per_cell);
1150 FEValues<dim> fe_values(discretization_->mapping(),
1151 discretization_->finite_element(),
1152 discretization_->quadrature(),
1153 update_values | update_JxW_values);
1154 for (
const auto &cell : dof_handler.cell_iterators_on_level(level))
1157 if (cell->is_locally_owned_on_level()) {
1158 fe_values.reinit(cell);
1159 for (
unsigned int i = 0; i < mass_values.size(); ++i) {
1161 for (
unsigned int q = 0; q < fe_values.n_quadrature_points; ++q)
1162 sum += fe_values.shape_value(i, q) * fe_values.JxW(q);
1163 mass_values(i) = sum;
1165 cell->get_mg_dof_indices(dof_indices);
1166 level_constraints.distribute_local_to_global(
1167 mass_values, dof_indices, level_lumped_mass_matrix_[level]);
1173 level_boundary_map_[level] = construct_boundary_map(
1174 dof_handler.begin_mg(level), dof_handler.end_mg(level), *partitioner);
1179 template <
int dim,
typename Number>
1180 template <
typename ITERATOR1,
typename ITERATOR2>
1182 const ITERATOR1 &begin,
1183 const ITERATOR2 &end,
1184 const Utilities::MPI::Partitioner &partitioner)
const -> BoundaryMap
1187 std::cout <<
"OfflineData<dim, Number>::construct_boundary_map()"
1195 using BoundaryData = std::tuple<dealii::Tensor<1, dim, Number> ,
1198 dealii::types::boundary_id ,
1199 dealii::Point<dim>> ;
1200 std::multimap<unsigned int, BoundaryData> preliminary_map;
1202 std::vector<dealii::types::global_dof_index> local_dof_indices;
1204 const dealii::QGauss<dim - 1> face_quadrature(3);
1205 dealii::FEFaceValues<dim> fe_face_values(discretization_->mapping(),
1206 discretization_->finite_element(),
1208 dealii::update_normal_vectors |
1209 dealii::update_values |
1210 dealii::update_JxW_values);
1212 const unsigned int dofs_per_cell =
1213 discretization_->finite_element().dofs_per_cell;
1215 const auto support_points =
1216 discretization_->finite_element().get_unit_support_points();
1218 for (
auto cell = begin; cell != end; ++cell) {
1226 if ((cell->is_active() && cell->is_artificial()) ||
1227 (!cell->is_active() && cell->is_artificial_on_level()))
1230 local_dof_indices.resize(dofs_per_cell);
1231 cell->get_active_or_mg_dof_indices(local_dof_indices);
1233 for (
auto f : cell->face_indices()) {
1234 const auto face = cell->face(f);
1235 const auto id = face->boundary_id();
1237 if (!face->at_boundary())
1248 fe_face_values.reinit(cell, f);
1250 for (
unsigned int j : fe_face_values.dof_indices()) {
1251 if (!discretization_->finite_element().has_support_on_face(j, f))
1254 Number boundary_mass = 0.;
1255 dealii::Tensor<1, dim, Number> normal;
1257 for (
unsigned int q : fe_face_values.quadrature_point_indices()) {
1258 const auto JxW = fe_face_values.JxW(q);
1259 const auto phi_i = fe_face_values.shape_value(j, q);
1261 boundary_mass += phi_i * JxW;
1262 normal += phi_i * fe_face_values.normal_vector(q) * JxW;
1265 const auto global_index = local_dof_indices[j];
1266 const auto index = partitioner.global_to_local(global_index);
1269 if (index >= n_locally_owned_)
1273 const unsigned int row_length =
1274 sparsity_pattern_simd_.row_length(index);
1275 if (row_length == 1)
1278 Point<dim> position =
1279 discretization_->mapping().transform_unit_to_real_cell(
1280 cell, support_points[j]);
1286 preliminary_map.insert(
1287 {index, {normal, boundary_mass, boundary_mass, id, position}});
1303 std::multimap<unsigned int, BoundaryData> filtered_map;
1304 std::set<dealii::types::global_dof_index> boundary_dofs;
1305 for (
auto entry : preliminary_map) {
1306 bool inserted =
false;
1307 const auto range = filtered_map.equal_range(entry.first);
1308 for (
auto it = range.first; it != range.second; ++it) {
1313 new_point] = entry.second;
1314 auto &[normal, normal_mass, boundary_mass, id, point] = it->second;
1319 Assert(point.distance(new_point) < 1.0e-14, dealii::ExcInternalError());
1321 if (normal * new_normal / normal.norm() / new_normal.norm() > 0.50) {
1326 normal += new_normal;
1327 boundary_mass += new_boundary_mass;
1331 }
else if constexpr (dim == 2) {
1347 filtered_map.insert(entry);
1354 BoundaryMap boundary_map;
1356 std::begin(filtered_map),
1357 std::end(filtered_map),
1358 std::back_inserter(boundary_map),
1360 auto index = it.first;
1361 const auto &[normal, normal_mass, boundary_mass, id, point] =
1364 const auto new_normal_mass =
1365 normal.norm() + std::numeric_limits<Number>::epsilon();
1366 const auto new_normal = normal / new_normal_mass;
1368 return {index, new_normal, new_normal_mass, boundary_mass, id, point};
1371 return boundary_map;
1375 template <
int dim,
typename Number>
1376 template <
typename ITERATOR1,
typename ITERATOR2>
1378 const ITERATOR1 &begin,
1379 const ITERATOR2 &end,
1380 const Utilities::MPI::Partitioner &partitioner)
const
1381 -> CouplingBoundaryPairs
1384 std::cout <<
"OfflineData<dim, Number>::collect_coupling_boundary_pairs()"
1394 std::set<unsigned int> locally_relevant_boundary_indices;
1396 std::vector<dealii::types::global_dof_index> local_dof_indices;
1398 const auto &finite_element = discretization_->finite_element();
1399 const unsigned int dofs_per_cell = finite_element.dofs_per_cell;
1400 local_dof_indices.resize(dofs_per_cell);
1402 for (
auto cell = begin; cell != end; ++cell) {
1405 if (cell->is_artificial())
1408 cell->get_active_or_mg_dof_indices(local_dof_indices);
1410 for (
auto f : cell->face_indices()) {
1411 const auto face = cell->face(f);
1412 const auto id = face->boundary_id();
1414 if (!face->at_boundary())
1421 for (
unsigned int j = 0; j < dofs_per_cell; ++j) {
1423 if (!finite_element.has_support_on_face(j, f))
1426 const auto global_index = local_dof_indices[j];
1427 const auto index = partitioner.global_to_local(global_index);
1430 if (index >= n_locally_relevant_)
1433 locally_relevant_boundary_indices.insert(index);
1442 CouplingBoundaryPairs result;
1444 for (
const auto i : locally_relevant_boundary_indices) {
1447 if (i >= n_locally_owned_)
1450 const unsigned int row_length = sparsity_pattern_simd_.row_length(i);
1453 if (row_length == 1)
1456 const unsigned int *js = sparsity_pattern_simd_.columns(i);
1457 constexpr auto simd_length = VectorizedArray<Number>::size();
1459 for (
unsigned int col_idx = 1; col_idx < row_length; ++col_idx) {
1460 const auto j = *(i < n_locally_internal_ ? js + col_idx * simd_length
1463 if (locally_relevant_boundary_indices.count(j) != 0) {
1464 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)