16#ifndef dealii_face_setup_internal_h
17#define dealii_face_setup_internal_h
40 namespace MatrixFreeFunctions
82 const unsigned int mg_level,
83 const bool hold_all_faces_to_owned_cells,
84 const bool build_inner_faces,
85 std::vector<std::pair<unsigned int, unsigned int>> &cell_levels);
96 const std::vector<std::pair<unsigned int, unsigned int>> &cell_levels,
107 const unsigned int face_no,
108 const typename ::Triangulation<dim>::cell_iterator &cell,
109 const unsigned int number_cell_interior,
110 const typename ::Triangulation<dim>::cell_iterator &neighbor,
111 const unsigned int number_cell_exterior,
112 const bool is_mixed_mesh);
142 template <
int vectorization_w
idth>
146 const std::vector<bool> &hard_vectorization_boundary,
147 std::vector<unsigned int> &face_partition_data,
158 : use_active_cells(true)
165 FaceSetup<dim>::initialize(
167 const unsigned int mg_level,
168 const bool hold_all_faces_to_owned_cells,
169 const bool build_inner_faces,
170 std::vector<std::pair<unsigned int, unsigned int>> &cell_levels)
176 if (use_active_cells)
177 for (
const auto &cell_level : cell_levels)
179 typename ::Triangulation<dim>::cell_iterator dcell(
188 at_processor_boundary.resize(cell_levels.size(),
false);
191 FaceCategory::locally_active_done_elsewhere);
195 std::map<types::subdomain_id, FaceIdentifier>
196 inner_faces_at_proc_boundary;
202 for (
unsigned int i = 0; i < cell_levels.size(); ++i)
204 if (i > 0 && cell_levels[i] == cell_levels[i - 1])
206 typename ::Triangulation<dim>::cell_iterator dcell(
208 for (
const unsigned int f : dcell->face_indices())
210 if (dcell->at_boundary(f) && !dcell->has_periodic_neighbor(f))
212 typename ::Triangulation<dim>::cell_iterator neighbor =
213 dcell->neighbor_or_periodic_neighbor(f);
219 const CellId id_mine = dcell->id();
220 if (use_active_cells && neighbor->has_children())
221 for (
unsigned int c = 0;
222 c < (dcell->has_periodic_neighbor(f) ?
223 dcell->periodic_neighbor(f)
224 ->face(dcell->periodic_neighbor_face_no(f))
226 dcell->face(f)->n_children());
229 typename ::Triangulation<dim>::cell_iterator
231 dcell->at_boundary(f) ?
232 dcell->periodic_neighbor_child_on_subface(f, c) :
233 dcell->neighbor_child_on_subface(f, c);
235 neighbor_c->subdomain_id();
236 if (my_domain < neigh_domain)
237 inner_faces_at_proc_boundary[neigh_domain]
238 .n_hanging_faces_larger_subdomain++;
239 else if (my_domain > neigh_domain)
240 inner_faces_at_proc_boundary[neigh_domain]
241 .n_hanging_faces_smaller_subdomain++;
246 use_active_cells ? neighbor->subdomain_id() :
247 neighbor->level_subdomain_id();
248 if (neighbor->level() < dcell->level() &&
251 if (my_domain < neigh_domain)
252 inner_faces_at_proc_boundary[neigh_domain]
253 .n_hanging_faces_smaller_subdomain++;
254 else if (my_domain > neigh_domain)
255 inner_faces_at_proc_boundary[neigh_domain]
256 .n_hanging_faces_larger_subdomain++;
258 else if (neighbor->level() == dcell->level() &&
259 my_domain != neigh_domain)
265 const CellId id_neigh = neighbor->id();
266 if (my_domain < neigh_domain)
267 inner_faces_at_proc_boundary[neigh_domain]
268 .shared_faces.emplace_back(id_mine, id_neigh);
270 inner_faces_at_proc_boundary[neigh_domain]
271 .shared_faces.emplace_back(id_neigh, id_mine);
280 for (
auto &inner_face : inner_faces_at_proc_boundary)
282 Assert(inner_face.first != my_domain,
284 std::sort(inner_face.second.shared_faces.begin(),
285 inner_face.second.shared_faces.end());
286 inner_face.second.shared_faces.erase(
287 std::unique(inner_face.second.shared_faces.begin(),
288 inner_face.second.shared_faces.end()),
289 inner_face.second.shared_faces.end());
294# if defined(DEAL_II_WITH_MPI) && defined(DEBUG)
296 if (const ::parallel::TriangulationBase<dim> *ptria =
297 dynamic_cast<const ::parallel::TriangulationBase<dim>
299 comm = ptria->get_communicator();
302 unsigned int mysize = inner_face.second.shared_faces.size();
305 int ierr = MPI_Sendrecv(&mysize,
314 600 + inner_face.first,
319 mysize = inner_face.second.n_hanging_faces_smaller_subdomain;
320 ierr = MPI_Sendrecv(&mysize,
329 700 + inner_face.first,
334 mysize = inner_face.second.n_hanging_faces_larger_subdomain;
335 ierr = MPI_Sendrecv(&mysize,
344 800 + inner_face.first,
366 std::vector<std::tuple<CellId, CellId, unsigned int>> other_range(
367 inner_face.second.shared_faces.size());
368 for (
unsigned int i = 0; i < other_range.size(); ++i)
370 std::make_tuple(inner_face.second.shared_faces[i].second,
371 inner_face.second.shared_faces[i].first,
373 std::sort(other_range.begin(), other_range.end());
381 unsigned int n_faces_lower_proc = 0, n_faces_higher_proc = 0;
382 std::vector<signed char> assignment(other_range.size(), 0);
383 if (inner_face.second.shared_faces.size() > 0)
387 unsigned int count = 0;
388 for (
unsigned int i = 1;
389 i < inner_face.second.shared_faces.size();
391 if (inner_face.second.shared_faces[i].first ==
392 inner_face.second.shared_faces[i - 1 - count].first)
399 for (
unsigned int k = 0; k <= count; ++k)
400 assignment[i - 1 - k] = 1;
401 n_faces_higher_proc += count + 1;
411 for (
unsigned int i = 1; i < other_range.size(); ++i)
412 if (std::get<0>(other_range[i]) ==
413 std::get<0>(other_range[i - 1 - count]))
420 for (
unsigned int k = 0; k <= count; ++k)
423 .shared_faces[std::get<2>(
427 .shared_faces[std::get<2>(
428 other_range[i - 1 - k])]
433 if (assignment[std::get<2>(
434 other_range[i - 1 - k])] == 0)
436 assignment[std::get<2>(
437 other_range[i - 1 - k])] = -1;
438 ++n_faces_lower_proc;
452 n_faces_lower_proc +=
453 inner_face.second.n_hanging_faces_smaller_subdomain;
454 n_faces_higher_proc +=
455 inner_face.second.n_hanging_faces_larger_subdomain;
456 const unsigned int n_total_faces_at_proc_boundary =
457 (inner_face.second.shared_faces.size() +
458 inner_face.second.n_hanging_faces_smaller_subdomain +
459 inner_face.second.n_hanging_faces_larger_subdomain);
460 unsigned int split_index = n_total_faces_at_proc_boundary / 2;
461 if (split_index < n_faces_lower_proc)
463 else if (split_index <
464 n_total_faces_at_proc_boundary - n_faces_higher_proc)
465 split_index -= n_faces_lower_proc;
467 split_index = n_total_faces_at_proc_boundary -
468 n_faces_higher_proc - n_faces_lower_proc;
471# if defined(DEAL_II_WITH_MPI) && defined(DEBUG)
472 ierr = MPI_Sendrecv(&split_index,
481 900 + inner_face.first,
486 ierr = MPI_Sendrecv(&n_faces_lower_proc,
495 1000 + inner_face.first,
500 ierr = MPI_Sendrecv(&n_faces_higher_proc,
509 1100 + inner_face.first,
517 std::vector<std::pair<CellId, CellId>> owned_faces_lower,
519 for (
unsigned int i = 0; i < assignment.size(); ++i)
520 if (assignment[i] < 0)
521 owned_faces_lower.push_back(
522 inner_face.second.shared_faces[i]);
523 else if (assignment[i] > 0)
524 owned_faces_higher.push_back(
525 inner_face.second.shared_faces[i]);
527 inner_face.second.shared_faces.size() + 1 -
528 owned_faces_lower.size() -
529 owned_faces_higher.size());
531 unsigned int i = 0, c = 0;
532 for (; i < assignment.size() && c < split_index; ++i)
533 if (assignment[i] == 0)
535 owned_faces_lower.push_back(
536 inner_face.second.shared_faces[i]);
539 for (; i < assignment.size(); ++i)
540 if (assignment[i] == 0)
542 owned_faces_higher.push_back(
543 inner_face.second.shared_faces[i]);
548 std::vector<std::pair<CellId, CellId>> check_faces;
549 check_faces.insert(check_faces.end(),
550 owned_faces_lower.begin(),
551 owned_faces_lower.end());
552 check_faces.insert(check_faces.end(),
553 owned_faces_higher.begin(),
554 owned_faces_higher.end());
555 std::sort(check_faces.begin(), check_faces.end());
557 inner_face.second.shared_faces.size());
558 for (
unsigned int i = 0; i < check_faces.size(); ++i)
559 Assert(check_faces[i] == inner_face.second.shared_faces[i],
564 if (my_domain < inner_face.first)
565 inner_face.second.shared_faces.swap(owned_faces_lower);
567 inner_face.second.shared_faces.swap(owned_faces_higher);
569 std::sort(inner_face.second.shared_faces.begin(),
570 inner_face.second.shared_faces.end());
576 std::set<std::pair<unsigned int, unsigned int>> ghost_cells;
577 for (
unsigned int i = 0; i < cell_levels.size(); ++i)
579 typename ::Triangulation<dim>::cell_iterator dcell(
581 if (use_active_cells)
583 for (
const auto f : dcell->face_indices())
585 if (dcell->at_boundary(f) && !dcell->has_periodic_neighbor(f))
586 face_is_owned[dcell->face(f)->index()] =
587 FaceCategory::locally_active_at_boundary;
588 else if (!build_inner_faces)
593 else if ((dcell->at_boundary(f) ==
false ||
594 dcell->has_periodic_neighbor(f)) &&
596 dcell->neighbor_or_periodic_neighbor(f)->level() <
599 face_is_owned[dcell->face(f)->index()] =
600 FaceCategory::multigrid_refinement_edge;
604 typename ::Triangulation<dim>::cell_iterator neighbor =
605 dcell->neighbor_or_periodic_neighbor(f);
608 if (use_active_cells && neighbor->has_children() &&
609 hold_all_faces_to_owned_cells ==
false)
612 bool add_to_ghost =
false;
614 id1 = use_active_cells ? dcell->subdomain_id() :
615 dcell->level_subdomain_id(),
616 id2 = use_active_cells ?
617 (neighbor->has_children() ?
618 dcell->neighbor_child_on_subface(f, 0)
620 neighbor->subdomain_id()) :
621 neighbor->level_subdomain_id();
629 (use_active_cells ==
false || neighbor->is_active())) ||
630 dcell->level() > neighbor->level() ||
632 inner_faces_at_proc_boundary[id2].shared_faces.begin(),
633 inner_faces_at_proc_boundary[id2].shared_faces.end(),
634 std::make_pair(id1 < id2 ? dcell->
id() : neighbor->id(),
635 id1 < id2 ? neighbor->id() :
638 face_is_owned[dcell->face(f)->index()] =
639 FaceCategory::locally_active_done_here;
640 if (dcell->level() == neighbor->level())
643 ->face(dcell->has_periodic_neighbor(f) ?
644 dcell->periodic_neighbor_face_no(f) :
645 dcell->neighbor_face_no(f))
647 FaceCategory::locally_active_done_here;
653 if (use_active_cells)
655 (dcell->subdomain_id() != neighbor->subdomain_id());
657 add_to_ghost = (dcell->level_subdomain_id() !=
658 neighbor->level_subdomain_id());
660 else if (hold_all_faces_to_owned_cells ==
true)
663 face_is_owned[dcell->face(f)->index()] =
664 FaceCategory::ghosted;
665 if (use_active_cells)
667 if (neighbor->has_children())
668 for (
unsigned int s = 0;
669 s < dcell->face(f)->n_children();
671 if (dcell->at_boundary(f))
674 ->periodic_neighbor_child_on_subface(f,
677 dcell->subdomain_id())
682 if (dcell->neighbor_child_on_subface(f, s)
684 dcell->subdomain_id())
688 add_to_ghost = (dcell->subdomain_id() !=
689 neighbor->subdomain_id());
692 add_to_ghost = (dcell->level_subdomain_id() !=
693 neighbor->level_subdomain_id());
698 if (use_active_cells && neighbor->has_children())
699 for (
unsigned int s = 0;
700 s < dcell->face(f)->n_children();
703 typename ::Triangulation<dim>::cell_iterator
705 dcell->at_boundary(f) ?
706 dcell->periodic_neighbor_child_on_subface(f,
708 dcell->neighbor_child_on_subface(f, s);
709 if (neighbor_child->subdomain_id() !=
710 dcell->subdomain_id())
712 std::pair<unsigned int, unsigned int>(
713 neighbor_child->level(),
714 neighbor_child->index()));
718 std::pair<unsigned int, unsigned int>(
719 neighbor->level(), neighbor->index()));
720 at_processor_boundary[i] =
true;
728 for (
const auto &ghost_cell : ghost_cells)
729 cell_levels.push_back(ghost_cell);
736 FaceSetup<dim>::generate_faces(
738 const std::vector<std::pair<unsigned int, unsigned int>> &cell_levels,
745 std::map<std::pair<unsigned int, unsigned int>,
unsigned int>
747 for (
unsigned int cell = 0; cell < cell_levels.size(); ++cell)
748 if (cell == 0 || cell_levels[cell] != cell_levels[cell - 1])
750 typename ::Triangulation<dim>::cell_iterator dcell(
752 cell_levels[cell].
first,
753 cell_levels[cell].
second);
754 std::pair<unsigned int, unsigned int> level_index(dcell->level(),
756 map_to_vectorized[level_index] = cell;
760 const unsigned int vectorization_length = task_info.vectorization_length;
761 task_info.face_partition_data.resize(
762 task_info.cell_partition_data.size() - 1, 0);
763 task_info.boundary_partition_data.resize(
764 task_info.cell_partition_data.size() - 1, 0);
765 std::vector<unsigned char> face_visited(face_is_owned.size(), 0);
766 for (
unsigned int partition = 0;
767 partition < task_info.cell_partition_data.size() - 2;
770 unsigned int boundary_counter = 0;
771 unsigned int inner_counter = 0;
772 for (
unsigned int cell = task_info.cell_partition_data[partition] *
773 vectorization_length;
774 cell < task_info.cell_partition_data[
partition + 1] *
775 vectorization_length;
777 if (cell == 0 || cell_levels[cell] != cell_levels[cell - 1])
779 typename ::Triangulation<dim>::cell_iterator dcell(
781 cell_levels[cell].
first,
782 cell_levels[cell].
second);
783 for (
const auto f : dcell->face_indices())
786 if (face_is_owned[dcell->face(f)->index()] ==
787 FaceCategory::locally_active_at_boundary)
791 FaceToCellTopology<1> info;
792 info.cells_interior[0] = cell;
794 info.interior_face_no = f;
795 info.exterior_face_no = dcell->face(f)->boundary_id();
798 (dcell->face(f)->reference_cell() !=
803 info.face_orientation = 0;
804 boundary_faces.push_back(info);
806 face_visited[dcell->face(f)->index()]++;
811 typename ::Triangulation<dim>::cell_iterator
812 neighbor = dcell->neighbor_or_periodic_neighbor(f);
813 if (use_active_cells && neighbor->has_children())
815 for (
unsigned int c = 0;
816 c < dcell->face(f)->n_children();
819 typename ::Triangulation<
820 dim>::cell_iterator neighbor_c =
821 dcell->at_boundary(f) ?
822 dcell->periodic_neighbor_child_on_subface(
824 dcell->neighbor_child_on_subface(f, c);
826 neighbor_c->subdomain_id();
827 const unsigned int neighbor_face_no =
828 dcell->has_periodic_neighbor(f) ?
829 dcell->periodic_neighbor_face_no(f) :
830 dcell->neighbor_face_no(f);
831 if (neigh_domain != dcell->subdomain_id() ||
833 [dcell->face(f)->child(c)->index()] ==
836 std::pair<unsigned int, unsigned int>
837 level_index(neighbor_c->level(),
838 neighbor_c->index());
840 [dcell->face(f)->child(c)->index()] ==
841 FaceCategory::locally_active_done_here)
844 inner_faces.push_back(create_face(
847 map_to_vectorized[level_index],
852 else if (face_is_owned[dcell->face(f)
855 FaceCategory::ghosted)
857 inner_ghost_faces.push_back(create_face(
860 map_to_vectorized[level_index],
867 face_is_owned[dcell->face(f)
870 locally_active_done_elsewhere ||
871 face_is_owned[dcell->face(f)
873 FaceCategory::ghosted,
879 [dcell->face(f)->child(c)->index()] = 1;
886 use_active_cells ? dcell->subdomain_id() :
887 dcell->level_subdomain_id();
889 use_active_cells ? neighbor->subdomain_id() :
890 neighbor->level_subdomain_id();
891 if (neigh_domain != my_domain ||
892 face_visited[dcell->face(f)->index()] == 1)
894 std::pair<unsigned int, unsigned int>
895 level_index(neighbor->level(),
897 if (face_is_owned[dcell->face(f)->index()] ==
898 FaceCategory::locally_active_done_here)
900 Assert(use_active_cells ||
905 inner_faces.push_back(create_face(
910 map_to_vectorized[level_index],
913 else if (face_is_owned[dcell->face(f)
915 FaceCategory::ghosted)
917 inner_ghost_faces.push_back(create_face(
922 map_to_vectorized[level_index],
928 face_visited[dcell->face(f)->index()] = 1;
929 if (dcell->has_periodic_neighbor(f))
933 dcell->periodic_neighbor_face_no(f))
936 if (face_is_owned[dcell->face(f)->index()] ==
937 FaceCategory::multigrid_refinement_edge)
939 refinement_edge_faces.push_back(
944 refinement_edge_faces.size(),
951 task_info.face_partition_data[
partition + 1] =
952 task_info.face_partition_data[
partition] + inner_counter;
953 task_info.boundary_partition_data[
partition + 1] =
954 task_info.boundary_partition_data[
partition] + boundary_counter;
956 task_info.ghost_face_partition_data.resize(2);
957 task_info.ghost_face_partition_data[0] = 0;
958 task_info.ghost_face_partition_data[1] = inner_ghost_faces.size();
959 task_info.refinement_edge_face_partition_data.resize(2);
960 task_info.refinement_edge_face_partition_data[0] = 0;
961 task_info.refinement_edge_face_partition_data[1] =
962 refinement_edge_faces.size();
968 FaceToCellTopology<1>
969 FaceSetup<dim>::create_face(
970 const unsigned int face_no,
971 const typename ::Triangulation<dim>::cell_iterator &cell,
972 const unsigned int number_cell_interior,
973 const typename ::Triangulation<dim>::cell_iterator &neighbor,
974 const unsigned int number_cell_exterior,
975 const bool is_mixed_mesh)
977 FaceToCellTopology<1> info;
978 info.cells_interior[0] = number_cell_interior;
979 info.cells_exterior[0] = number_cell_exterior;
980 info.interior_face_no = face_no;
981 if (cell->has_periodic_neighbor(face_no))
982 info.exterior_face_no = cell->periodic_neighbor_face_no(face_no);
984 info.exterior_face_no = cell->neighbor_face_no(face_no);
986 info.face_type = is_mixed_mesh ?
987 (cell->face(face_no)->reference_cell() !=
993 if (cell->level() > neighbor->level())
995 if (cell->has_periodic_neighbor(face_no))
997 cell->periodic_neighbor_of_coarser_periodic_neighbor(face_no)
1000 info.subface_index =
1001 cell->neighbor_of_coarser_neighbor(face_no).second;
1005 if (dim == 3 && cell->has_periodic_neighbor(face_no))
1007 unsigned char exterior_face_orientation = cell->get_triangulation()
1008 .get_periodic_face_map()
1009 .at({cell, face_no})
1011 const auto [orientation, rotation, flip] =
1013 exterior_face_orientation);
1015 info.face_orientation =
1016 (orientation ? 0u : 1u) + 2 * flip + 4 * rotation;
1021 info.face_orientation = 0;
1022 const unsigned int interior_face_orientation =
1023 !cell->face_orientation(face_no) + 2 * cell->face_flip(face_no) +
1024 4 * cell->face_rotation(face_no);
1025 const unsigned int exterior_face_orientation =
1026 !neighbor->face_orientation(info.exterior_face_no) +
1027 2 * neighbor->face_flip(info.exterior_face_no) +
1028 4 * neighbor->face_rotation(info.exterior_face_no);
1029 if (interior_face_orientation != 0)
1031 info.face_orientation = 8 + interior_face_orientation;
1032 Assert(exterior_face_orientation == 0,
1034 "Face seems to be wrongly oriented from both sides"));
1037 info.face_orientation = exterior_face_orientation;
1041 if (cell->level() > neighbor->level() && exterior_face_orientation > 0)
1044 ShapeInfo<double>::compute_orientation_table(2);
1045 const std::array<unsigned int, 8> inverted_orientations{
1046 {0, 1, 2, 3, 6, 5, 4, 7}};
1047 info.subface_index =
1048 orientation[inverted_orientations[exterior_face_orientation]]
1049 [info.subface_index];
1064 compare_faces_for_vectorization(
1065 const FaceToCellTopology<1> &face1,
1066 const FaceToCellTopology<1> &face2,
1067 const std::vector<unsigned int> &active_fe_indices,
1068 const unsigned int length)
1070 if (face1.interior_face_no != face2.interior_face_no)
1072 if (face1.exterior_face_no != face2.exterior_face_no)
1074 if (face1.subface_index != face2.subface_index)
1076 if (face1.face_orientation != face2.face_orientation)
1078 if (face1.face_type != face2.face_type)
1081 if (active_fe_indices.size() > 0)
1083 if (active_fe_indices[face1.cells_interior[0] / length] !=
1084 active_fe_indices[face2.cells_interior[0] / length])
1088 if (active_fe_indices[face1.cells_exterior[0] / length] !=
1089 active_fe_indices[face2.cells_exterior[0] / length])
1104 template <
int length>
1105 struct FaceComparator
1107 FaceComparator(
const std::vector<unsigned int> &active_fe_indices)
1108 : active_fe_indices(active_fe_indices)
1112 operator()(
const FaceToCellTopology<length> &face1,
1113 const FaceToCellTopology<length> &face2)
const
1116 if (face1.face_type < face2.face_type)
1118 else if (face1.face_type > face2.face_type)
1122 if (active_fe_indices.size() > 0)
1125 if (active_fe_indices[face1.cells_interior[0] / length] <
1126 active_fe_indices[face2.cells_interior[0] / length])
1128 else if (active_fe_indices[face1.cells_interior[0] / length] >
1129 active_fe_indices[face2.cells_interior[0] / length])
1135 if (active_fe_indices[face1.cells_exterior[0] / length] <
1136 active_fe_indices[face2.cells_exterior[0] / length])
1138 else if (active_fe_indices[face1.cells_exterior[0] / length] >
1139 active_fe_indices[face2.cells_exterior[0] / length])
1144 for (
unsigned int i = 0; i < length; ++i)
1145 if (face1.cells_interior[i] < face2.cells_interior[i])
1147 else if (face1.cells_interior[i] > face2.cells_interior[i])
1149 for (
unsigned int i = 0; i < length; ++i)
1150 if (face1.cells_exterior[i] < face2.cells_exterior[i])
1152 else if (face1.cells_exterior[i] > face2.cells_exterior[i])
1154 if (face1.interior_face_no < face2.interior_face_no)
1156 else if (face1.interior_face_no > face2.interior_face_no)
1158 if (face1.exterior_face_no < face2.exterior_face_no)
1160 else if (face1.exterior_face_no > face2.exterior_face_no)
1173 const std::vector<unsigned int> &active_fe_indices;
1178 template <
int vectorization_w
idth>
1181 const std::vector<FaceToCellTopology<1>> &faces_in,
1182 const std::vector<bool> &hard_vectorization_boundary,
1183 std::vector<unsigned int> &face_partition_data,
1184 std::vector<FaceToCellTopology<vectorization_width>> &faces_out,
1185 const std::vector<unsigned int> &active_fe_indices)
1187 FaceToCellTopology<vectorization_width> face_batch;
1188 std::vector<std::vector<unsigned int>> faces_type;
1190 unsigned int face_start = face_partition_data[0],
1191 face_end = face_partition_data[0];
1193 face_partition_data[0] = faces_out.size();
1194 for (
unsigned int partition = 0;
1195 partition < face_partition_data.size() - 1;
1198 std::vector<std::vector<unsigned int>> new_faces_type;
1201 face_start = face_end;
1202 face_end = face_partition_data[
partition + 1];
1209 for (
unsigned int face = face_start; face < face_end; ++face)
1211 for (
auto &face_type : faces_type)
1214 if (compare_faces_for_vectorization(faces_in[face],
1215 faces_in[face_type[0]],
1217 vectorization_width))
1219 face_type.push_back(face);
1223 faces_type.emplace_back(1, face);
1229 FaceComparator<vectorization_width> face_comparator(
1231 std::set<FaceToCellTopology<vectorization_width>,
1232 FaceComparator<vectorization_width>>
1233 new_faces(face_comparator);
1234 for (
const auto &face_type : faces_type)
1236 face_batch.face_type = faces_in[face_type[0]].face_type;
1237 face_batch.interior_face_no =
1238 faces_in[face_type[0]].interior_face_no;
1239 face_batch.exterior_face_no =
1240 faces_in[face_type[0]].exterior_face_no;
1241 face_batch.subface_index = faces_in[face_type[0]].subface_index;
1242 face_batch.face_orientation =
1243 faces_in[face_type[0]].face_orientation;
1244 unsigned int no_faces = face_type.size();
1245 std::vector<unsigned char> touched(no_faces, 0);
1251 unsigned int n_vectorized = 0;
1252 for (
unsigned int f = 0; f < no_faces; ++f)
1253 if (faces_in[face_type[f]].cells_interior[0] %
1254 vectorization_width ==
1257 bool is_contiguous =
true;
1258 if (f + vectorization_width > no_faces)
1259 is_contiguous =
false;
1261 for (
unsigned int v = 1; v < vectorization_width; ++v)
1262 if (faces_in[face_type[f + v]].cells_interior[0] !=
1263 faces_in[face_type[f]].cells_interior[0] + v)
1264 is_contiguous =
false;
1269 vectorization_width + 1);
1270 for (
unsigned int v = 0; v < vectorization_width; ++v)
1272 face_batch.cells_interior[v] =
1273 faces_in[face_type[f + v]].cells_interior[0];
1274 face_batch.cells_exterior[v] =
1275 faces_in[face_type[f + v]].cells_exterior[0];
1278 new_faces.insert(face_batch);
1279 f += vectorization_width - 1;
1280 n_vectorized += vectorization_width;
1284 std::vector<unsigned int> untouched;
1285 untouched.reserve(no_faces - n_vectorized);
1286 for (
unsigned int f = 0; f < no_faces; ++f)
1287 if (touched[f] == 0)
1288 untouched.push_back(f);
1290 for (
const auto f : untouched)
1292 face_batch.cells_interior[v] =
1293 faces_in[face_type[f]].cells_interior[0];
1294 face_batch.cells_exterior[v] =
1295 faces_in[face_type[f]].cells_exterior[0];
1297 if (v == vectorization_width)
1299 new_faces.insert(face_batch);
1303 if (v > 0 && v < vectorization_width)
1306 if (hard_vectorization_boundary[partition + 1] ||
1307 partition == face_partition_data.size() - 2)
1309 for (; v < vectorization_width; ++v)
1312 face_batch.cells_interior[v] =
1314 face_batch.cells_exterior[v] =
1317 new_faces.insert(face_batch);
1322 std::vector<unsigned int> untreated(v);
1323 for (
unsigned int f = 0; f < v; ++f)
1324 untreated[f] = face_type[*(untouched.end() - 1 - f)];
1325 new_faces_type.push_back(untreated);
1331 for (
auto it = new_faces.begin(); it != new_faces.end(); ++it)
1332 faces_out.push_back(*it);
1333 face_partition_data[
partition + 1] += new_faces.size();
1336 faces_type = std::move(new_faces_type);
1341 for (
const auto &face_type : faces_type)
1345 unsigned int nfaces = 0;
1346 for (
unsigned int i = face_partition_data[0];
1347 i < face_partition_data.back();
1349 for (
unsigned int v = 0; v < vectorization_width; ++v)
1354 std::vector<std::pair<unsigned int, unsigned int>> in_faces, out_faces;
1355 for (
const auto &face_in : faces_in)
1356 in_faces.emplace_back(face_in.cells_interior[0],
1357 face_in.cells_exterior[0]);
1358 for (
unsigned int i = face_partition_data[0];
1359 i < face_partition_data.back();
1361 for (
unsigned int v = 0;
1362 v < vectorization_width &&
1365 out_faces.emplace_back(faces_out[i].cells_interior[v],
1366 faces_out[i].cells_exterior[v]);
1367 std::sort(in_faces.begin(), in_faces.end());
1368 std::sort(out_faces.begin(), out_faces.end());
1370 for (
unsigned int i = 0; i < in_faces.size(); ++i)
unsigned int n_raw_faces() const
bool is_mixed_mesh() const
types::subdomain_id locally_owned_subdomain() const override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
constexpr const ReferenceCell & get_hypercube()
void collect_faces_vectorization(const std::vector< FaceToCellTopology< 1 > > &faces_in, const std::vector< bool > &hard_vectorization_boundary, std::vector< unsigned int > &face_partition_data, std::vector< FaceToCellTopology< vectorization_width > > &faces_out)
std::tuple< bool, bool, bool > split_face_orientation(const unsigned char combined_face_orientation)
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
*braid_SplitCommworld & comm
unsigned int n_hanging_faces_smaller_subdomain
unsigned int n_hanging_faces_larger_subdomain
std::vector< std::pair< CellId, CellId > > shared_faces
std::vector< FaceToCellTopology< 1 > > inner_faces
std::vector< bool > at_processor_boundary
std::vector< FaceToCellTopology< 1 > > boundary_faces
std::vector< FaceToCellTopology< 1 > > refinement_edge_faces
FaceToCellTopology< 1 > create_face(const unsigned int face_no, const typename ::Triangulation< dim >::cell_iterator &cell, const unsigned int number_cell_interior, const typename ::Triangulation< dim >::cell_iterator &neighbor, const unsigned int number_cell_exterior, const bool is_mixed_mesh)
std::vector< FaceCategory > face_is_owned
@ locally_active_done_here
@ multigrid_refinement_edge
@ locally_active_at_boundary
@ locally_active_done_elsewhere
void initialize(const ::Triangulation< dim > &triangulation, const unsigned int mg_level, const bool hold_all_faces_to_owned_cells, const bool build_inner_faces, std::vector< std::pair< unsigned int, unsigned int > > &cell_levels)
std::vector< FaceToCellTopology< 1 > > inner_ghost_faces
void generate_faces(const ::Triangulation< dim > &triangulation, const std::vector< std::pair< unsigned int, unsigned int > > &cell_levels, TaskInfo &task_info)