16#ifndef dealii_hanging_nodes_internal_h
17#define dealii_hanging_nodes_internal_h
35 namespace MatrixFreeFunctions
94 const std::uint16_t kind =
static_cast<std::uint16_t
>(kind_in);
95 const std::uint16_t subcell = (kind >> 0) & 7;
96 const std::uint16_t face = (kind >> 3) & 7;
97 const std::uint16_t edge = (kind >> 6) & 7;
107 if (subcell == 0 && face == 0)
114 if (subcell == 0 && face == 0 && edge == 0)
116 else if (0 < face && edge == 0)
118 else if (0 == face && 0 < edge)
120 else if ((face == edge) && (face == 1 || face == 2 || face == 4))
144 const std::uint16_t kind =
static_cast<std::uint16_t
>(kind_in);
145 const std::uint16_t subcell = (kind >> 0) & 7;
146 const std::uint16_t face = (kind >> 3) & 7;
147 const std::uint16_t edge = (kind >> 6) & 7;
149 return subcell + ((face > 0) << 3) + ((edge > 0) << 4) +
168 const std::uint16_t subcell = (kind_in >> 0) & 7;
169 const std::uint16_t flag_0 = (kind_in >> 3) & 3;
170 const std::uint16_t flag_1 = (kind_in >> 5) & 7;
173 subcell + (((flag_0 & 0b01) ? flag_1 : 0) << 3) +
174 (((flag_0 & 0b10) ? flag_1 : 0) << 6));
205 static_cast<std::uint16_t
>(f2));
229 return static_cast<std::uint16_t
>(f1) !=
static_cast<std::uint16_t
>(f2);
241 return static_cast<std::uint16_t
>(f1) <
static_cast<std::uint16_t
>(f2);
253 static_cast<std::uint16_t
>(f2));
277 template <
typename CellIterator>
280 const CellIterator & cell,
281 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
282 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
283 std::vector<types::global_dof_index> & dof_indices,
290 static std::vector<std::vector<bool>>
296 template <
typename CellIterator>
304 template <
typename CellIterator>
307 const CellIterator & cell,
308 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
309 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
310 const std::vector<std::vector<bool>> & component_mask,
312 std::vector<types::global_dof_index> & dof_indices)
const;
326 unsigned int n_dofs_1d,
327 std::vector<types::global_dof_index> &dofs)
const;
332 unsigned int n_dofs_1d)
const;
336 std::vector<types::global_dof_index> &dofs)
const;
341 std::vector<std::vector<
342 std::pair<typename Triangulation<dim>::cell_iterator,
unsigned int>>>
346 {{{{{7, 3}}, {{6, 2}}}},
347 {{{{5, 1}}, {{4, 0}}}},
348 {{{{11, 9}}, {{10, 8}}}}}};
386 return !cell.is_artificial();
390 const unsigned int n_raw_lines =
triangulation.n_raw_lines();
391 this->line_to_cells.resize(n_raw_lines);
399 const unsigned int line_to_children[12][2] = {{0, 2},
412 std::vector<std::vector<
413 std::pair<typename Triangulation<3>::cell_iterator,
unsigned int>>>
414 line_to_inactive_cells(n_raw_lines);
419 for (
unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
422 const unsigned int line_idx = cell->
line(line)->index();
424 line_to_cells[line_idx].push_back(std::make_pair(cell, line));
426 line_to_inactive_cells[line_idx].push_back(
427 std::make_pair(cell, line));
434 for (
unsigned int line_idx = 0; line_idx < n_raw_lines; ++line_idx)
436 if ((line_to_cells[line_idx].size() > 0) &&
437 line_to_inactive_cells[line_idx].size() > 0)
441 const auto &inactive_cell =
442 line_to_inactive_cells[line_idx][0].first;
443 const unsigned int neighbor_line =
444 line_to_inactive_cells[line_idx][0].second;
446 for (
unsigned int c = 0; c < 2; ++c)
449 inactive_cell->child(line_to_children[neighbor_line][c]);
450 const unsigned int child_line_idx =
451 child->line(neighbor_line)->index();
454 for (
const auto &cl : line_to_cells[line_idx])
455 line_to_cells[child_line_idx].push_back(cl);
464 inline std::vector<std::vector<bool>>
466 const ::hp::FECollection<dim> &fe_collection)
468 std::vector<std::vector<bool>> supported_components(
469 fe_collection.size(),
470 std::vector<bool>(fe_collection.n_components(),
false));
472 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
474 for (
unsigned int base_element_index = 0, comp = 0;
475 base_element_index < fe_collection[i].n_base_elements();
476 ++base_element_index)
477 for (
unsigned int c = 0;
478 c < fe_collection[i].element_multiplicity(base_element_index);
480 if (dim == 1 ||
dynamic_cast<const FE_Q<dim> *
>(
481 &fe_collection[i].base_element(
482 base_element_index)) ==
nullptr)
483 supported_components[i][comp] =
false;
485 supported_components[i][comp] =
true;
488 return supported_components;
494 template <
typename CellIterator>
497 const CellIterator &cell)
const
500 if ((dim == 3 && line_to_cells.size() == 0) ||
501 (cell->reference_cell().is_hyper_cube() ==
false))
504 if (cell->level() == 0)
507 const std::uint16_t subcell =
508 cell->parent()->child_iterator_to_index(cell);
509 const std::uint16_t
subcell_x = (subcell >> 0) & 1;
510 const std::uint16_t
subcell_y = (subcell >> 1) & 1;
511 const std::uint16_t
subcell_z = (subcell >> 2) & 1;
513 std::uint16_t face = 0;
514 std::uint16_t edge = 0;
516 for (
unsigned int direction = 0; direction < dim; ++direction)
518 const auto side = (subcell >> direction) & 1U;
519 const auto face_no = direction * 2 + side;
522 if (cell->at_boundary(face_no))
525 const auto &neighbor = cell->neighbor(face_no);
529 if (neighbor->has_children() || neighbor->is_artificial() ||
530 neighbor->level() == cell->level())
533 face |= 1 << direction;
537 for (
unsigned int direction = 0; direction < dim; ++direction)
538 if (face == 0 || face == (1 << direction))
540 const unsigned int line_no =
547 const unsigned int line_index = cell->line(line_no)->index();
549 const auto edge_neighbor =
550 std::find_if(line_to_cells[line_index].begin(),
551 line_to_cells[line_index].end(),
552 [&cell](
const auto &edge_neighbor) {
553 return edge_neighbor.first->is_artificial() ==
555 edge_neighbor.first->level() <
559 if (edge_neighbor == line_to_cells[line_index].end())
562 edge |= 1 << direction;
565 if ((face == 0) && (edge == 0))
568 const std::uint16_t inverted_subcell = (subcell ^ (dim == 2 ? 3 : 7));
571 inverted_subcell + (face << 3) + (edge << 6));
573 return refinement_configuration;
579 template <
typename CellIterator>
582 const CellIterator & cell,
583 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
584 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
585 const std::vector<std::vector<bool>> & supported_components,
587 std::vector<types::global_dof_index> & dof_indices)
const
589 if (std::find(supported_components[cell->active_fe_index()].begin(),
590 supported_components[cell->active_fe_index()].end(),
592 supported_components[cell->active_fe_index()].end())
595 const auto &fe = cell->get_fe();
599 std::vector<std::vector<unsigned int>>
600 component_to_system_index_face_array(fe.n_components());
602 for (
unsigned int i = 0; i < fe.n_dofs_per_face(0); ++i)
603 component_to_system_index_face_array
604 [fe.face_system_to_component_index(i, 0).first]
607 std::vector<unsigned int> idx_offset = {0};
609 for (
unsigned int base_element_index = 0;
610 base_element_index < cell->get_fe().n_base_elements();
611 ++base_element_index)
612 for (
unsigned int c = 0;
613 c < cell->get_fe().element_multiplicity(base_element_index);
615 idx_offset.push_back(
617 cell->get_fe().base_element(base_element_index).n_dofs_per_cell());
619 std::vector<types::global_dof_index> neighbor_dofs_all(idx_offset.back());
620 std::vector<types::global_dof_index> neighbor_dofs_all_temp(
622 std::vector<types::global_dof_index> neighbor_dofs_face(
623 fe.n_dofs_per_face(0));
626 const auto get_face_idx = [](
const auto n_dofs_1d,
629 const auto j) ->
unsigned int {
630 const auto direction = face_no / 2;
631 const auto side = face_no % 2;
632 const auto offset = (side == 1) ? (n_dofs_1d - 1) : 0;
635 return (direction == 0) ? (n_dofs_1d * i + offset) :
636 (n_dofs_1d * offset + i);
641 return n_dofs_1d * n_dofs_1d * i + n_dofs_1d * j + offset;
643 return n_dofs_1d * n_dofs_1d * j + n_dofs_1d * offset + i;
645 return n_dofs_1d * n_dofs_1d * offset + n_dofs_1d * i + j;
655 const std::uint16_t kind =
656 static_cast<std::uint16_t
>(refinement_configuration);
657 const std::uint16_t subcell = (kind >> 0) & 7;
658 const std::uint16_t
subcell_x = (subcell >> 0) & 1;
659 const std::uint16_t
subcell_y = (subcell >> 1) & 1;
660 const std::uint16_t
subcell_z = (subcell >> 2) & 1;
661 const std::uint16_t face = (kind >> 3) & 7;
662 const std::uint16_t edge = (kind >> 6) & 7;
664 for (
unsigned int direction = 0; direction < dim; ++direction)
665 if ((face >> direction) & 1U)
667 const auto side = ((subcell >> direction) & 1U) == 0;
668 const auto face_no = direction * 2 + side;
671 cell->neighbor(face_no)
672 ->face(cell->neighbor_face_no(face_no))
673 ->get_dof_indices(neighbor_dofs_face,
674 cell->neighbor(face_no)->active_fe_index());
678 for (
auto &
index : neighbor_dofs_face)
681 for (
unsigned int base_element_index = 0, comp = 0;
682 base_element_index < cell->get_fe().n_base_elements();
683 ++base_element_index)
684 for (
unsigned int c = 0;
685 c < cell->get_fe().element_multiplicity(base_element_index);
688 if (supported_components[cell->active_fe_index()][comp] ==
692 const unsigned int n_dofs_1d =
694 .base_element(base_element_index)
697 const unsigned int dofs_per_face =
699 std::vector<types::global_dof_index> neighbor_dofs(
701 const auto lex_face_mapping =
706 for (
unsigned int i = 0; i < dofs_per_face; ++i)
707 neighbor_dofs[i] = neighbor_dofs_face
708 [component_to_system_index_face_array[comp][i]];
715 Assert(cell->face_orientation(face_no),
721 if (cell->face_rotation(face_no))
723 if (cell->face_flip(face_no))
726 rotate_face(rotate, n_dofs_1d, neighbor_dofs);
728 if (cell->face_orientation(face_no) ==
false)
729 transpose_face(n_dofs_1d - 1, neighbor_dofs);
737 for (
unsigned int i = 0, k = 0; i < n_dofs_1d; ++i)
738 for (
unsigned int j = 0; j < (dim == 2 ? 1 : n_dofs_1d);
740 dof_indices[get_face_idx(n_dofs_1d, face_no, i, j) +
742 neighbor_dofs[lex_face_mapping[k]];
747 for (
unsigned int direction = 0; direction < dim; ++direction)
748 if ((edge >> direction) & 1U)
750 const unsigned int line_no =
756 const unsigned int line_index = cell->line(line_no)->index();
758 const auto edge_neighbor =
759 std::find_if(line_to_cells[line_index].begin(),
760 line_to_cells[line_index].end(),
761 [&cell](
const auto &edge_neighbor) {
762 return edge_neighbor.first->is_artificial() ==
764 edge_neighbor.first->level() <
768 if (edge_neighbor == line_to_cells[line_index].end())
771 const auto neighbor_cell = edge_neighbor->first;
772 const auto local_line_neighbor = edge_neighbor->second;
775 &neighbor_cell->get_triangulation(),
776 neighbor_cell->level(),
777 neighbor_cell->index(),
778 &cell->get_dof_handler())
782 for (
auto &
index : neighbor_dofs_all)
785 for (
unsigned int i = 0; i < neighbor_dofs_all_temp.size(); ++i)
786 neighbor_dofs_all_temp[i] = neighbor_dofs_all
787 [lexicographic_mapping[cell->active_fe_index()][i]];
790 cell->line_orientation(line_no) !=
791 neighbor_cell->line_orientation(local_line_neighbor);
793 for (
unsigned int base_element_index = 0, comp = 0;
794 base_element_index < cell->get_fe().n_base_elements();
795 ++base_element_index)
796 for (
unsigned int c = 0;
798 cell->get_fe().element_multiplicity(base_element_index);
801 if (supported_components[cell->active_fe_index()][comp] ==
805 const unsigned int n_dofs_1d =
807 .base_element(base_element_index)
811 for (
unsigned int i = 0; i < n_dofs_1d; ++i)
812 dof_indices[line_dof_idx(line_no, i, n_dofs_1d) +
813 idx_offset[comp]] = neighbor_dofs_all_temp
814 [line_dof_idx(local_line_neighbor,
815 flipped ? (n_dofs_1d - 1 - i) : i,
825 template <
typename CellIterator>
828 const CellIterator & cell,
829 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
830 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
831 std::vector<types::global_dof_index> & dof_indices,
835 const auto supported_components = compute_supported_components(
836 cell->get_dof_handler().get_fe_collection());
838 if ([](
const auto &supported_components) {
839 return std::none_of(supported_components.begin(),
840 supported_components.end(),
842 return *std::max_element(a.begin(), a.end());
844 }(supported_components))
848 const auto refinement_configuration =
849 compute_refinement_configuration(cell);
855 update_dof_indices(cell,
857 lexicographic_mapping,
858 supported_components,
859 refinement_configuration,
863 for (
unsigned int c = 0; c < supported_components[0].size(); ++c)
864 if (supported_components[cell->active_fe_index()][c])
865 masks[c] = refinement_configuration;
875 unsigned int &subface_index)
const
877 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
880 times = times < 0 ? times + 4 : times;
881 for (
int t = 0; t < times; ++t)
882 subface_index = rot_mapping[subface_index];
891 unsigned int n_dofs_1d,
892 std::vector<types::global_dof_index> &dofs)
const
894 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
897 times = times < 0 ? times + 4 : times;
899 std::vector<types::global_dof_index> copy(dofs.size());
900 for (
int t = 0; t < times; ++t)
902 std::swap(copy, dofs);
905 for (
unsigned int i = 0; i < 4; ++i)
906 dofs[rot_mapping[i]] = copy[i];
909 const unsigned int n_int = n_dofs_1d - 2;
910 unsigned int offset = 4;
911 for (
unsigned int i = 0; i < n_int; ++i)
914 dofs[offset + i] = copy[offset + 2 * n_int + (n_int - 1 - i)];
916 dofs[offset + n_int + i] =
917 copy[offset + 3 * n_int + (n_int - 1 - i)];
919 dofs[offset + 2 * n_int + i] = copy[offset + n_int + i];
921 dofs[offset + 3 * n_int + i] = copy[offset + i];
927 for (
unsigned int i = 0; i < n_int; ++i)
928 for (
unsigned int j = 0; j < n_int; ++j)
929 dofs[offset + i * n_int + j] =
930 copy[offset + j * n_int + (n_int - 1 - i)];
940 unsigned int n_dofs_1d)
const
942 unsigned int x, y, z;
944 const unsigned int fe_degree = n_dofs_1d - 1;
948 x = (local_line % 4 == 0) ? 0 :
949 (local_line % 4 == 1) ? fe_degree :
951 y = (local_line % 4 == 2) ? 0 :
952 (local_line % 4 == 3) ? fe_degree :
954 z = (local_line / 4) * fe_degree;
958 x = ((local_line - 8) % 2) * fe_degree;
959 y = ((local_line - 8) / 2) * fe_degree;
963 return n_dofs_1d * n_dofs_1d * z + n_dofs_1d * y + x;
971 const unsigned int fe_degree,
972 std::vector<types::global_dof_index> &dofs)
const
974 const std::vector<types::global_dof_index> copy(dofs);
981 const unsigned int n_int = fe_degree - 1;
982 unsigned int offset = 4;
983 for (
unsigned int i = 0; i < n_int; ++i)
986 dofs[offset + i] = copy[offset + 2 * n_int + i];
988 dofs[offset + n_int + i] = copy[offset + 3 * n_int + i];
990 dofs[offset + 2 * n_int + i] = copy[offset + i];
992 dofs[offset + 3 * n_int + i] = copy[offset + n_int + i];
997 for (
unsigned int i = 0; i < n_int; ++i)
998 for (
unsigned int j = 0; j < n_int; ++j)
999 dofs[offset + i * n_int + j] = copy[offset + j * n_int + i];
1010 else if (subface == 2)
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
ConstraintKinds compute_refinement_configuration(const CellIterator &cell) const
void setup_line_to_cell(const Triangulation< dim > &triangulation)
void transpose_face(const unsigned int fe_degree, std::vector< types::global_dof_index > &dofs) const
static std::vector< std::vector< bool > > compute_supported_components(const ::hp::FECollection< dim > &fe)
std::vector< std::vector< std::pair< typename Triangulation< dim >::cell_iterator, unsigned int > > > line_to_cells
HangingNodes(const Triangulation< dim > &triangualtion)
bool setup_constraints(const CellIterator &cell, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner, const std::vector< std::vector< unsigned int > > &lexicographic_mapping, std::vector< types::global_dof_index > &dof_indices, const ArrayView< ConstraintKinds > &mask) const
void transpose_subface_index(unsigned int &subface) const
void rotate_face(int times, unsigned int n_dofs_1d, std::vector< types::global_dof_index > &dofs) const
unsigned int line_dof_idx(int local_line, unsigned int dof, unsigned int n_dofs_1d) const
void update_dof_indices(const CellIterator &cell, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner, const std::vector< std::vector< unsigned int > > &lexicographic_mapping, const std::vector< std::vector< bool > > &component_mask, const ConstraintKinds &refinement_configuration, std::vector< types::global_dof_index > &dof_indices) const
const ::ndarray< unsigned int, 3, 2, 2 > local_lines
void rotate_subface_index(int times, unsigned int &subface_index) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
constexpr T pow(const T base, const int iexp)
bool operator<(const ConstraintKinds f1, const ConstraintKinds f2)
ConstraintKinds decompress(const compressed_constraint_kind kind_in, const unsigned int dim)
bool check(const ConstraintKinds kind_in, const unsigned int dim)
compressed_constraint_kind compress(const ConstraintKinds kind_in, const unsigned int dim)
ConstraintKinds & operator|=(ConstraintKinds &f1, const ConstraintKinds f2)
std::size_t memory_consumption(const ConstraintKinds &)
std::uint8_t compressed_constraint_kind
ConstraintKinds operator&(const ConstraintKinds f1, const ConstraintKinds f2)
bool operator!=(const ConstraintKinds f1, const ConstraintKinds f2)
constexpr compressed_constraint_kind unconstrained_compressed_constraint_kind
ConstraintKinds operator|(const ConstraintKinds f1, const ConstraintKinds f2)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
#define DEAL_II_CUDA_HOST_DEV
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation