deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
hanging_nodes_internal.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2018 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_hanging_nodes_internal_h
16#define dealii_hanging_nodes_internal_h
17
18#include <deal.II/base/config.h>
19
22
24
25#include <deal.II/fe/fe_q.h>
27#include <deal.II/fe/fe_tools.h>
28
30
31#include <boost/container/small_vector.hpp>
32
33
35
36namespace internal
37{
38 namespace MatrixFreeFunctions
39 {
53 enum class ConstraintKinds : std::uint16_t
54 {
55 // default: unconstrained cell
57
58 // subcell
59 subcell_x = 1 << 0,
60 subcell_y = 1 << 1,
61 subcell_z = 1 << 2,
62
63 // face is constrained
64 face_x = 1 << 3,
65 face_y = 1 << 4,
66 face_z = 1 << 5,
67
68 // edge is constrained
69 edge_x = 1 << 6,
70 edge_y = 1 << 7,
71 edge_z = 1 << 8
72 };
73
74
75
79 using compressed_constraint_kind = std::uint8_t;
80
81
82
86 constexpr compressed_constraint_kind
88
89
90
94 inline bool
95 check(const ConstraintKinds kind_in, const unsigned int dim)
96 {
97 const std::uint16_t kind = static_cast<std::uint16_t>(kind_in);
98 const std::uint16_t subcell = (kind >> 0) & 7;
99 const std::uint16_t face = (kind >> 3) & 7;
100 const std::uint16_t edge = (kind >> 6) & 7;
101
102 if ((kind >> 9) > 0)
103 return false;
104
105 if (dim == 2)
106 {
107 if (edge > 0)
108 return false; // in 2d there are no edge constraints
109
110 if (subcell == 0 && face == 0)
111 return true; // no constraints
112 else if (0 < face)
113 return true; // at least one face is constrained
114 }
115 else if (dim == 3)
116 {
117 if (subcell == 0 && face == 0 && edge == 0)
118 return true; // no constraints
119 else if (0 < face && edge == 0)
120 return true; // at least one face is constrained
121 else if (0 == face && 0 < edge)
122 return true; // at least one edge is constrained
123 else if ((face == edge) && (face == 1 || face == 2 || face == 4))
124 return true; // one face and its orthogonal edge is constrained
125 }
126
127 return false;
128 }
129
130
131
136 inline compressed_constraint_kind
137 compress(const ConstraintKinds kind_in, const unsigned int dim)
138 {
139 Assert(check(kind_in, dim), ExcInternalError());
140
141 if (dim == 2)
142 return static_cast<compressed_constraint_kind>(kind_in);
143
144 if (kind_in == ConstraintKinds::unconstrained)
146
147 const std::uint16_t kind = static_cast<std::uint16_t>(kind_in);
148 const std::uint16_t subcell = (kind >> 0) & 7;
149 const std::uint16_t face = (kind >> 3) & 7;
150 const std::uint16_t edge = (kind >> 6) & 7;
151
152 return subcell + ((face > 0) << 3) + ((edge > 0) << 4) +
153 (std::max(face, edge) << 5);
154 }
155
156
157
162 inline ConstraintKinds
163 decompress(const compressed_constraint_kind kind_in, const unsigned int dim)
164 {
165 if (dim == 2)
166 return static_cast<ConstraintKinds>(kind_in);
167
170
171 const std::uint16_t subcell = (kind_in >> 0) & 7;
172 const std::uint16_t flag_0 = (kind_in >> 3) & 3;
173 const std::uint16_t flag_1 = (kind_in >> 5) & 7;
174
175 const auto result = static_cast<ConstraintKinds>(
176 subcell + (((flag_0 & 0b01) ? flag_1 : 0) << 3) +
177 (((flag_0 & 0b10) ? flag_1 : 0) << 6));
178
179 Assert(check(result, dim), ExcInternalError());
180
181 return result;
182 }
183
184
185
189 inline std::size_t
191 {
192 return sizeof(ConstraintKinds);
193 }
194
195
196
204 DEAL_II_HOST_DEVICE inline ConstraintKinds
206 {
207 return static_cast<ConstraintKinds>(static_cast<std::uint16_t>(f1) |
208 static_cast<std::uint16_t>(f2));
209 }
210
211
212
217 DEAL_II_HOST_DEVICE inline ConstraintKinds &
219 {
220 f1 = f1 | f2;
221 return f1;
222 }
223
224
225
229 DEAL_II_HOST_DEVICE inline bool
231 {
232 return static_cast<std::uint16_t>(f1) != static_cast<std::uint16_t>(f2);
233 }
234
235
236
242 operator<(const ConstraintKinds f1, const ConstraintKinds f2)
243 {
244 return static_cast<std::uint16_t>(f1) < static_cast<std::uint16_t>(f2);
245 }
246
247
248
252 DEAL_II_HOST_DEVICE inline ConstraintKinds
254 {
255 return static_cast<ConstraintKinds>(static_cast<std::uint16_t>(f1) &
256 static_cast<std::uint16_t>(f2));
257 }
258
259
260
267 template <int dim>
269 {
270 public:
274 HangingNodes(const Triangulation<dim> &triangualtion);
275
279 std::size_t
280 memory_consumption() const;
281
285 template <typename CellIterator>
286 bool
288 const CellIterator &cell,
289 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
290 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
291 std::vector<types::global_dof_index> &dof_indices,
292 const ArrayView<ConstraintKinds> &mask) const;
293
298 static std::vector<std::vector<bool>>
299 compute_supported_components(const ::hp::FECollection<dim> &fe);
300
304 template <typename CellIterator>
306 compute_refinement_configuration(const CellIterator &cell) const;
307
312 template <typename CellIterator>
313 void
315 const CellIterator &cell,
316 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
317 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
318 const std::vector<std::vector<bool>> &component_mask,
319 const ConstraintKinds &refinement_configuration,
320 std::vector<types::global_dof_index> &dof_indices) const;
321
322 private:
326 void
327 setup_line_to_cell(const Triangulation<dim> &triangulation);
328
329 void
330 rotate_subface_index(int times, unsigned int &subface_index) const;
331
332 void
333 orient_face(const types::geometric_orientation combined_orientation,
334 const unsigned int n_dofs_1d,
335 std::vector<types::global_dof_index> &dofs) const;
336
337 unsigned int
338 line_dof_idx(int local_line,
339 unsigned int dof,
340 unsigned int n_dofs_1d) const;
341
342 void
343 transpose_subface_index(unsigned int &subface) const;
344
345 std::vector<
346 boost::container::small_vector<std::array<unsigned int, 3>, 6>>
348
349 const ::ndarray<unsigned int, 3, 2, 2> local_lines = {
350 {{{{{7, 3}}, {{6, 2}}}},
351 {{{{5, 1}}, {{4, 0}}}},
352 {{{{11, 9}}, {{10, 8}}}}}};
353 };
354
355
356
357 template <int dim>
359 const Triangulation<dim> &triangulation)
360 {
361 // Set up line-to-cell mapping for edge constraints (only if dim = 3 and
362 // for pure hex meshes)
363 if (triangulation.all_reference_cells_are_hyper_cube())
364 setup_line_to_cell(triangulation);
365 }
366
367
368
369 template <int dim>
370 inline std::size_t
372 {
373 std::size_t size = 0;
374 for (const auto &a : line_to_cells)
375 size +=
376 (a.capacity() > 6 ? a.capacity() : 0) * sizeof(a[0]) + sizeof(a);
377 return size;
378 }
379
380
381
382 template <int dim>
383 inline void
385 const Triangulation<dim> & /*triangulation*/)
386 {}
387
388
389
390 template <>
391 inline void
393 {
394 // Check if we there are no hanging nodes on the current MPI process,
395 // which we do by checking if the second finest level holds no active
396 // non-artificial cell
397 if (triangulation.n_levels() <= 1 ||
398 std::none_of(triangulation.begin_active(triangulation.n_levels() - 2),
399 triangulation.end_active(triangulation.n_levels() - 2),
400 [](const CellAccessor<3, 3> &cell) {
401 return !cell.is_artificial();
402 }))
403 return;
404
405 const unsigned int n_raw_lines = triangulation.n_raw_lines();
406 this->line_to_cells.resize(n_raw_lines);
407
408 // In 3d, we can have DoFs on only an edge being constrained (e.g. in a
409 // cartesian 2x2x2 grid, where only the upper left 2 cells are refined).
410 // This sets up a helper data structure in the form of a mapping from
411 // edges (i.e. lines) to neighboring cells.
412
413 // Mapping from an edge to which children that share that edge.
414 const unsigned int line_to_children[12][2] = {{0, 2},
415 {1, 3},
416 {0, 1},
417 {2, 3},
418 {4, 6},
419 {5, 7},
420 {4, 5},
421 {6, 7},
422 {0, 4},
423 {1, 5},
424 {2, 6},
425 {3, 7}};
426
427 std::vector<
428 boost::container::small_vector<std::array<unsigned int, 3>, 6>>
429 line_to_inactive_cells(n_raw_lines);
430
431 // First add active and inactive cells to their lines:
432 for (const auto &cell : triangulation.cell_iterators())
433 {
434 const unsigned int cell_level = cell->level();
435 const unsigned int cell_index = cell->index();
436 for (unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
437 ++line)
438 {
439 const unsigned int line_idx = cell->line_index(line);
440 if (cell->is_active())
441 line_to_cells[line_idx].push_back(
442 {{cell_level, cell_index, line}});
443 else
444 line_to_inactive_cells[line_idx].push_back(
445 {{cell_level, cell_index, line}});
446 }
447 }
448
449 // Now, we can access edge-neighboring active cells on same level to also
450 // access of an edge to the edges "children". These are found from looking
451 // at the corresponding edge of children of inactive edge neighbors.
452 for (unsigned int line_idx = 0; line_idx < n_raw_lines; ++line_idx)
453 {
454 if ((line_to_cells[line_idx].size() > 0) &&
455 line_to_inactive_cells[line_idx].size() > 0)
456 {
457 // We now have cells to add (active ones) and edges to which they
458 // should be added (inactive cells).
459 const Triangulation<3>::cell_iterator inactive_cell(
460 &triangulation,
461 line_to_inactive_cells[line_idx][0][0],
462 line_to_inactive_cells[line_idx][0][1]);
463 const unsigned int neighbor_line =
464 line_to_inactive_cells[line_idx][0][2];
465
466 for (unsigned int c = 0; c < 2; ++c)
467 {
468 const auto &child =
469 inactive_cell->child(line_to_children[neighbor_line][c]);
470 const unsigned int child_line_idx =
471 child->line_index(neighbor_line);
472
473 Assert(child->is_active(), ExcInternalError());
474
475 // Now add all active cells
476 for (const auto &cl : line_to_cells[line_idx])
477 line_to_cells[child_line_idx].push_back(cl);
478 }
479 }
480 }
481 }
482
483
484
485 template <int dim>
486 inline std::vector<std::vector<bool>>
488 const ::hp::FECollection<dim> &fe_collection)
489 {
490 std::vector<std::vector<bool>> supported_components(
491 fe_collection.size(),
492 std::vector<bool>(fe_collection.n_components(), false));
493
494 for (unsigned int i = 0; i < fe_collection.size(); ++i)
495 {
496 for (unsigned int base_element_index = 0, comp = 0;
497 base_element_index < fe_collection[i].n_base_elements();
498 ++base_element_index)
499 for (unsigned int c = 0;
500 c < fe_collection[i].element_multiplicity(base_element_index);
501 ++c, ++comp)
502 if (dim == 1 ||
503 (dynamic_cast<const FE_Q<dim> *>(
504 &fe_collection[i].base_element(base_element_index)) ==
505 nullptr &&
506 dynamic_cast<const FE_Q_iso_Q1<dim> *>(
507 &fe_collection[i].base_element(base_element_index)) ==
508 nullptr))
509 supported_components[i][comp] = false;
510 else
511 supported_components[i][comp] = true;
512 }
513
514 return supported_components;
515 }
516
517
518
519 template <int dim>
520 template <typename CellIterator>
521 inline ConstraintKinds
523 const CellIterator &cell) const
524 {
525 // TODO: for simplex or mixed meshes: nothing to do
526 if ((dim == 3 && line_to_cells.empty()) ||
527 (cell->reference_cell().is_hyper_cube() == false))
529
530 if (cell->level() == 0)
532
533 const std::uint16_t subcell =
534 cell->parent()->child_iterator_to_index(cell);
535 const std::uint16_t subcell_x = (subcell >> 0) & 1;
536 const std::uint16_t subcell_y = (subcell >> 1) & 1;
537 const std::uint16_t subcell_z = (subcell >> 2) & 1;
538
539 std::uint16_t face = 0;
540 std::uint16_t edge = 0;
541
542 for (unsigned int direction = 0; direction < dim; ++direction)
543 {
544 const auto side = (subcell >> direction) & 1U;
545 const auto face_no = direction * 2 + side;
546
547 // ignore if at boundary
548 if (cell->at_boundary(face_no))
549 continue;
550
551 const auto &neighbor = cell->neighbor(face_no);
552
553 // ignore neighbors that are artificial or have the same level or
554 // have children
555 if (neighbor->has_children() || neighbor->is_artificial() ||
556 neighbor->level() == cell->level())
557 continue;
558
559 // Ignore if the neighbors are FE_Nothing
560 if (neighbor->get_fe().n_dofs_per_cell() == 0)
561 continue;
562
563 face |= 1 << direction;
564 }
565
566 if (dim == 3)
567 for (unsigned int direction = 0; direction < dim; ++direction)
568 if (face == 0 || face == (1 << direction))
569 {
570 const unsigned int line_no =
571 direction == 0 ?
572 (local_lines[0][subcell_y == 0][subcell_z == 0]) :
573 (direction == 1 ?
574 (local_lines[1][subcell_x == 0][subcell_z == 0]) :
575 (local_lines[2][subcell_x == 0][subcell_y == 0]));
576
577 const unsigned int line_index = cell->line_index(line_no);
578
579 const auto edge_neighbor =
580 std::find_if(line_to_cells[line_index].begin(),
581 line_to_cells[line_index].end(),
582 [&cell](const auto &edge_neighbor) {
584 &cell->get_triangulation(),
585 edge_neighbor[0],
586 edge_neighbor[1],
587 &cell->get_dof_handler());
588 return dof_cell.is_artificial() == false &&
589 dof_cell.level() < cell->level() &&
590 dof_cell.get_fe().n_dofs_per_cell() > 0;
591 });
592
593 if (edge_neighbor == line_to_cells[line_index].end())
594 continue;
595
596 edge |= 1 << direction;
597 }
598
599 if ((face == 0) && (edge == 0))
601
602 const std::uint16_t inverted_subcell = (subcell ^ (dim == 2 ? 3 : 7));
603
604 const auto refinement_configuration = static_cast<ConstraintKinds>(
605 inverted_subcell + (face << 3) + (edge << 6));
606 Assert(check(refinement_configuration, dim), ExcInternalError());
607 return refinement_configuration;
608 }
609
610
611
612 template <int dim>
613 template <typename CellIterator>
614 inline void
616 const CellIterator &cell,
617 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
618 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
619 const std::vector<std::vector<bool>> &supported_components,
620 const ConstraintKinds &refinement_configuration,
621 std::vector<types::global_dof_index> &dof_indices) const
622 {
623 if (std::find(supported_components[cell->active_fe_index()].begin(),
624 supported_components[cell->active_fe_index()].end(),
625 true) ==
626 supported_components[cell->active_fe_index()].end())
627 return;
628
629 const auto &fe = cell->get_fe();
630 AssertDimension(fe.n_unique_faces(), 1);
631
632 std::vector<std::vector<unsigned int>>
633 component_to_system_index_face_array(fe.n_components());
634
635 for (unsigned int i = 0; i < fe.n_dofs_per_face(0); ++i)
636 component_to_system_index_face_array
637 [fe.face_system_to_component_index(i, /*face_no=*/0).first]
638 .push_back(i);
639
640 std::vector<unsigned int> idx_offset = {0};
641
642 for (unsigned int base_element_index = 0;
643 base_element_index < fe.n_base_elements();
644 ++base_element_index)
645 for (unsigned int c = 0;
646 c < fe.element_multiplicity(base_element_index);
647 ++c)
648 idx_offset.push_back(
649 idx_offset.back() +
650 fe.base_element(base_element_index).n_dofs_per_cell());
651
652 std::vector<types::global_dof_index> neighbor_dofs_all(idx_offset.back());
653 std::vector<types::global_dof_index> neighbor_dofs_all_temp(
654 idx_offset.back());
655 std::vector<types::global_dof_index> neighbor_dofs_face(
656 fe.n_dofs_per_face(/*face_no=*/0));
657
658
659 const auto get_face_idx = [](const auto n_dofs_1d,
660 const auto face_no,
661 const auto i,
662 const auto j) -> unsigned int {
663 const auto direction = face_no / 2;
664 const auto side = face_no % 2;
665 const auto offset = (side == 1) ? (n_dofs_1d - 1) : 0;
666
667 if (dim == 2)
668 return (direction == 0) ? (n_dofs_1d * i + offset) :
669 (n_dofs_1d * offset + i);
670 else if (dim == 3)
671 switch (direction)
672 {
673 case 0:
674 return n_dofs_1d * n_dofs_1d * i + n_dofs_1d * j + offset;
675 case 1:
676 return n_dofs_1d * n_dofs_1d * j + n_dofs_1d * offset + i;
677 case 2:
678 return n_dofs_1d * n_dofs_1d * offset + n_dofs_1d * i + j;
679 default:
681 }
682
684
685 return 0;
686 };
687
688 const std::uint16_t kind =
689 static_cast<std::uint16_t>(refinement_configuration);
690 const std::uint16_t subcell = (kind >> 0) & 7;
691 const std::uint16_t subcell_x = (subcell >> 0) & 1;
692 const std::uint16_t subcell_y = (subcell >> 1) & 1;
693 const std::uint16_t subcell_z = (subcell >> 2) & 1;
694 const std::uint16_t face = (kind >> 3) & 7;
695 const std::uint16_t edge = (kind >> 6) & 7;
696
697 for (unsigned int direction = 0; direction < dim; ++direction)
698 if ((face >> direction) & 1U)
699 {
700 const auto side = ((subcell >> direction) & 1U) == 0;
701 const auto face_no = direction * 2 + side;
702
703 // read DoFs of parent of face, ...
704 cell->neighbor(face_no)
705 ->face(cell->neighbor_face_no(face_no))
706 ->get_dof_indices(neighbor_dofs_face,
707 cell->neighbor(face_no)->active_fe_index());
708
709 // ... convert the global DoFs to serial ones, and ...
710 if (partitioner)
711 for (auto &index : neighbor_dofs_face)
712 index = partitioner->global_to_local(index);
713
714 for (unsigned int base_element_index = 0, comp = 0;
715 base_element_index < fe.n_base_elements();
716 ++base_element_index)
717 for (unsigned int c = 0;
718 c < fe.element_multiplicity(base_element_index);
719 ++c, ++comp)
720 {
721 if (supported_components[cell->active_fe_index()][comp] ==
722 false)
723 continue;
724
725 const unsigned int n_dofs_1d =
726 cell->get_fe()
727 .base_element(base_element_index)
728 .tensor_degree() +
729 1;
730 const unsigned int dofs_per_face =
731 Utilities::pow(n_dofs_1d, dim - 1);
732 std::vector<types::global_dof_index> neighbor_dofs(
733 dofs_per_face);
734 const auto lex_face_mapping =
736 n_dofs_1d - 1);
737
738 // ... extract the DoFs of the current component
739 for (unsigned int i = 0; i < dofs_per_face; ++i)
740 neighbor_dofs[i] = neighbor_dofs_face
741 [component_to_system_index_face_array[comp][i]];
742
743 // fix DoFs depending on orientation, rotation, and flip
744 if (dim == 2)
745 {
746 // TODO: this needs to be implemented for simplices but
747 // all-quad meshes are OK
748 Assert(cell->combined_face_orientation(face_no) ==
751 }
752 else if (dim == 3)
753 {
754 orient_face(cell->combined_face_orientation(face_no),
755 n_dofs_1d,
756 neighbor_dofs);
757 }
758 else
759 {
761 }
762
763 // update DoF map
764 for (unsigned int i = 0, k = 0; i < n_dofs_1d; ++i)
765 for (unsigned int j = 0; j < (dim == 2 ? 1 : n_dofs_1d);
766 ++j, ++k)
767 dof_indices[get_face_idx(n_dofs_1d, face_no, i, j) +
768 idx_offset[comp]] =
769 neighbor_dofs[lex_face_mapping[k]];
770 }
771 }
772
773 if (dim == 3)
774 for (unsigned int direction = 0; direction < dim; ++direction)
775 if ((edge >> direction) & 1U)
776 {
777 const unsigned int line_no =
778 direction == 0 ?
780 (direction == 1 ? (local_lines[1][subcell_x][subcell_z]) :
782
783 const unsigned int line_index = cell->line(line_no)->index();
784
785 const auto edge_neighbor =
786 std::find_if(line_to_cells[line_index].begin(),
787 line_to_cells[line_index].end(),
788 [&cell](const auto &edge_array) {
790 edge_neighbor(&cell->get_triangulation(),
791 edge_array[0],
792 edge_array[1]);
793 return edge_neighbor->is_artificial() == false &&
794 edge_neighbor->level() < cell->level();
795 });
796
797 if (edge_neighbor == line_to_cells[line_index].end())
798 continue;
799
800 const DoFCellAccessor<dim, dim, false> neighbor_cell(
801 &cell->get_triangulation(),
802 (*edge_neighbor)[0],
803 (*edge_neighbor)[1],
804 &cell->get_dof_handler());
805 const auto local_line_neighbor = (*edge_neighbor)[2];
806
807 neighbor_cell.get_dof_indices(neighbor_dofs_all);
808
809 if (partitioner)
810 for (auto &index : neighbor_dofs_all)
811 index = partitioner->global_to_local(index);
812
813 for (unsigned int i = 0; i < neighbor_dofs_all_temp.size(); ++i)
814 neighbor_dofs_all_temp[i] = neighbor_dofs_all
815 [lexicographic_mapping[cell->active_fe_index()][i]];
816
817 const bool flipped =
818 cell->line_orientation(line_no) !=
819 neighbor_cell.line_orientation(local_line_neighbor);
820
821 for (unsigned int base_element_index = 0, comp = 0;
822 base_element_index < fe.n_base_elements();
823 ++base_element_index)
824 for (unsigned int c = 0;
825 c < fe.element_multiplicity(base_element_index);
826 ++c, ++comp)
827 {
828 if (supported_components[cell->active_fe_index()][comp] ==
829 false)
830 continue;
831
832 const unsigned int n_dofs_1d =
833 cell->get_fe()
834 .base_element(base_element_index)
835 .tensor_degree() +
836 1;
837
838 for (unsigned int i = 0; i < n_dofs_1d; ++i)
839 dof_indices[line_dof_idx(line_no, i, n_dofs_1d) +
840 idx_offset[comp]] = neighbor_dofs_all_temp
841 [line_dof_idx(local_line_neighbor,
842 flipped ? (n_dofs_1d - 1 - i) : i,
843 n_dofs_1d) +
844 idx_offset[comp]];
845 }
846 }
847 }
848
849
850
851 template <int dim>
852 template <typename CellIterator>
853 inline bool
855 const CellIterator &cell,
856 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
857 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
858 std::vector<types::global_dof_index> &dof_indices,
859 const ArrayView<ConstraintKinds> &masks) const
860 {
861 // 1) check if finite elements support fast hanging-node algorithm
862 const auto supported_components = compute_supported_components(
863 cell->get_dof_handler().get_fe_collection());
864
865 if (std::none_of(supported_components.begin(),
866 supported_components.end(),
867 [](const auto &a) {
868 return *std::max_element(a.begin(), a.end());
869 }))
870 return false;
871
872 // 2) determine the refinement configuration of the cell
873 const auto refinement_configuration =
875
876 if (refinement_configuration == ConstraintKinds::unconstrained)
877 return false;
878
879 // 3) update DoF indices of cell for specified components
881 partitioner,
882 lexicographic_mapping,
883 supported_components,
884 refinement_configuration,
885 dof_indices);
886
887 // 4) TODO: copy refinement configuration to all components
888 for (unsigned int c = 0; c < supported_components[0].size(); ++c)
889 if (supported_components[cell->active_fe_index()][c])
890 masks[c] = refinement_configuration;
891
892 return true;
893 }
894
895
896
897 template <int dim>
898 inline void
900 unsigned int &subface_index) const
901 {
902 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
903
904 times = times % 4;
905 times = times < 0 ? times + 4 : times;
906 for (int t = 0; t < times; ++t)
907 subface_index = rot_mapping[subface_index];
908 }
909
910
911
912 template <int dim>
913 inline void
915 const types::geometric_orientation combined_orientation,
916 const unsigned int n_dofs_1d,
917 std::vector<types::global_dof_index> &dofs) const
918 {
919 const auto [orientation, rotation, flip] =
920 ::internal::split_face_orientation(combined_orientation);
921 const int n_rotations =
922 rotation || flip ? 4 - int(rotation) - 2 * int(flip) : 0;
923
924 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
925 Assert(n_dofs_1d > 1, ExcInternalError());
926 // 'per line' has the same meaning here as in FiniteElementData, i.e., the
927 // number of dofs assigned to a line (not including vertices).
928 const unsigned int dofs_per_line = n_dofs_1d - 2;
929
930 // rotate:
931 std::vector<types::global_dof_index> copy(dofs.size());
932 for (int t = 0; t < n_rotations; ++t)
933 {
934 std::swap(copy, dofs);
935
936 // Vertices
937 for (unsigned int i = 0; i < 4; ++i)
938 dofs[rot_mapping[i]] = copy[i];
939
940 // Edges
941 unsigned int offset = 4;
942 for (unsigned int i = 0; i < dofs_per_line; ++i)
943 {
944 // Left edge
945 dofs[offset + i] =
946 copy[offset + 2 * dofs_per_line + (dofs_per_line - 1 - i)];
947 // Right edge
948 dofs[offset + dofs_per_line + i] =
949 copy[offset + 3 * dofs_per_line + (dofs_per_line - 1 - i)];
950 // Bottom edge
951 dofs[offset + 2 * dofs_per_line + i] =
952 copy[offset + dofs_per_line + i];
953 // Top edge
954 dofs[offset + 3 * dofs_per_line + i] = copy[offset + i];
955 }
956
957 // Interior points
958 offset += 4 * dofs_per_line;
959
960 for (unsigned int i = 0; i < dofs_per_line; ++i)
961 for (unsigned int j = 0; j < dofs_per_line; ++j)
962 dofs[offset + i * dofs_per_line + j] =
963 copy[offset + j * dofs_per_line + (dofs_per_line - 1 - i)];
964 }
965
966 // transpose (note that we are using the standard geometric orientation
967 // here so orientation = true is the default):
968 if (!orientation)
969 {
970 copy = dofs;
971
972 // Vertices
973 dofs[1] = copy[2];
974 dofs[2] = copy[1];
975
976 // Edges
977 unsigned int offset = 4;
978 for (unsigned int i = 0; i < dofs_per_line; ++i)
979 {
980 // Right edge
981 dofs[offset + i] = copy[offset + 2 * dofs_per_line + i];
982 // Left edge
983 dofs[offset + dofs_per_line + i] =
984 copy[offset + 3 * dofs_per_line + i];
985 // Bottom edge
986 dofs[offset + 2 * dofs_per_line + i] = copy[offset + i];
987 // Top edge
988 dofs[offset + 3 * dofs_per_line + i] =
989 copy[offset + dofs_per_line + i];
990 }
991
992 // Interior
993 offset += 4 * dofs_per_line;
994 for (unsigned int i = 0; i < dofs_per_line; ++i)
995 for (unsigned int j = 0; j < dofs_per_line; ++j)
996 dofs[offset + i * dofs_per_line + j] =
997 copy[offset + j * dofs_per_line + i];
998 }
999 }
1000
1001
1002
1003 template <int dim>
1004 inline unsigned int
1006 unsigned int dof,
1007 unsigned int n_dofs_1d) const
1008 {
1009 unsigned int x, y, z;
1010
1011 const unsigned int fe_degree = n_dofs_1d - 1;
1012
1013 if (local_line < 8)
1014 {
1015 x = (local_line % 4 == 0) ? 0 :
1016 (local_line % 4 == 1) ? fe_degree :
1017 dof;
1018 y = (local_line % 4 == 2) ? 0 :
1019 (local_line % 4 == 3) ? fe_degree :
1020 dof;
1021 z = (local_line / 4) * fe_degree;
1022 }
1023 else
1024 {
1025 x = ((local_line - 8) % 2) * fe_degree;
1026 y = ((local_line - 8) / 2) * fe_degree;
1027 z = dof;
1028 }
1029
1030 return n_dofs_1d * n_dofs_1d * z + n_dofs_1d * y + x;
1031 }
1032
1033
1034
1035 template <int dim>
1036 void
1038 {
1039 if (subface == 1)
1040 subface = 2;
1041 else if (subface == 2)
1042 subface = 1;
1043 }
1044 } // namespace MatrixFreeFunctions
1045} // namespace internal
1046
1047
1049
1050#endif
bool is_active() const
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
bool is_artificial() const
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
int index() const
int level() const
unsigned int line_index(const unsigned int i) const
bool all_reference_cells_are_hyper_cube() const
unsigned int n_raw_lines() const
unsigned int n_levels() const
active_cell_iterator end_active(const unsigned int level) const
active_cell_iterator begin_active(const unsigned int level=0) const
ConstraintKinds compute_refinement_configuration(const CellIterator &cell) const
void setup_line_to_cell(const Triangulation< dim > &triangulation)
void orient_face(const types::geometric_orientation combined_orientation, const unsigned int n_dofs_1d, std::vector< types::global_dof_index > &dofs) const
std::vector< boost::container::small_vector< std::array< unsigned int, 3 >, 6 > > line_to_cells
static std::vector< std::vector< bool > > compute_supported_components(const ::hp::FECollection< dim > &fe)
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
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
Definition config.h:40
#define DEAL_II_HOST_DEVICE
Definition config.h:171
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_NOT_IMPLEMENTED()
IteratorRange< cell_iterator > cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
std::vector< unsigned int > lexicographic_to_hierarchic_numbering(unsigned int degree)
constexpr T pow(const T base, const int iexp)
Definition utilities.h:967
std::uint8_t compressed_constraint_kind
Definition dof_info.h:86
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 &)
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)
std::tuple< bool, bool, bool > split_face_orientation(const types::geometric_orientation combined_face_orientation)
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned char geometric_orientation
Definition types.h:40