Loading [MathJax]/extensions/TeX/AMSsymbols.js
 deal.II version GIT relicensing-3083-g7b89508ac7 2025-04-18 12:50:00+00:00
\(\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\}}\)
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
face_setup_internal.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2018 - 2024 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
16#ifndef dealii_face_setup_internal_h
17#define dealii_face_setup_internal_h
18
19#include <deal.II/base/config.h>
20
22
24
25#include <deal.II/grid/tria.h>
27
31
32#include <fstream>
33#include <set>
34
35
37
38
39namespace internal
40{
41 namespace MatrixFreeFunctions
42 {
58
59
60
69 template <int dim>
70 struct FaceSetup
71 {
73
80 void
82 const ::Triangulation<dim> &triangulation,
83 const unsigned int mg_level,
84 const bool hold_all_faces_to_owned_cells,
85 const bool build_inner_faces,
86 std::vector<std::pair<unsigned int, unsigned int>> &cell_levels);
87
94 void
96 const ::Triangulation<dim> &triangulation,
97 const std::vector<std::pair<unsigned int, unsigned int>> &cell_levels,
98 TaskInfo &task_info);
99
108 const unsigned int face_no,
109 const typename ::Triangulation<dim>::cell_iterator &cell,
110 const unsigned int number_cell_interior,
111 const typename ::Triangulation<dim>::cell_iterator &neighbor,
112 const unsigned int number_cell_exterior,
113 const bool is_mixed_mesh);
114
116
129
130 std::vector<FaceCategory> face_is_owned;
131 std::vector<bool> at_processor_boundary;
132 std::vector<FaceToCellTopology<1>> inner_faces;
133 std::vector<FaceToCellTopology<1>> boundary_faces;
134 std::vector<FaceToCellTopology<1>> inner_ghost_faces;
135 std::vector<FaceToCellTopology<1>> refinement_edge_faces;
136 };
137
138
139
143 template <int vectorization_width>
144 void
146 const std::vector<FaceToCellTopology<1>> &faces_in,
147 const std::vector<bool> &hard_vectorization_boundary,
148 std::vector<unsigned int> &face_partition_data,
149 std::vector<FaceToCellTopology<vectorization_width>> &faces_out);
150
151
152
153 /* -------------------------------------------------------------------- */
154
155#ifndef DOXYGEN
156
157 template <int dim>
159 : use_active_cells(true)
160 {}
161
162
163
164 template <int dim>
165 void
166 FaceSetup<dim>::initialize(
167 const ::Triangulation<dim> &triangulation,
168 const unsigned int mg_level,
169 const bool hold_all_faces_to_owned_cells,
170 const bool build_inner_faces,
171 std::vector<std::pair<unsigned int, unsigned int>> &cell_levels)
172 {
173 use_active_cells = mg_level == numbers::invalid_unsigned_int;
174
175 if constexpr (running_in_debug_mode())
176 {
177 // safety check
178 if (use_active_cells)
179 for (const auto &cell_level : cell_levels)
180 {
181 typename ::Triangulation<dim>::cell_iterator dcell(
182 &triangulation, cell_level.first, cell_level.second);
183 Assert(dcell->is_active(), ExcInternalError());
184 }
185 }
186
187 // step 1: add ghost cells for those cells that we identify as
188 // interesting
189
190 at_processor_boundary.resize(cell_levels.size(), false);
191 face_is_owned.resize(dim > 1 ? triangulation.n_raw_faces() :
192 triangulation.n_vertices(),
193 FaceCategory::locally_active_done_elsewhere);
194
195 // go through the mesh and divide the faces on the processor
196 // boundaries as evenly as possible between the processors
197 std::map<types::subdomain_id, FaceIdentifier>
198 inner_faces_at_proc_boundary;
199 if (triangulation.locally_owned_subdomain() !=
201 {
202 const types::subdomain_id my_domain =
203 triangulation.locally_owned_subdomain();
204 for (unsigned int i = 0; i < cell_levels.size(); ++i)
205 {
206 if (i > 0 && cell_levels[i] == cell_levels[i - 1])
207 continue;
208 typename ::Triangulation<dim>::cell_iterator dcell(
209 &triangulation, cell_levels[i].first, cell_levels[i].second);
210 for (const unsigned int f : dcell->face_indices())
211 {
212 if (dcell->at_boundary(f) && !dcell->has_periodic_neighbor(f))
213 continue;
214 typename ::Triangulation<dim>::cell_iterator neighbor =
215 dcell->neighbor_or_periodic_neighbor(f);
216
217 // faces at hanging nodes are always treated by the processor
218 // who owns the element on the fine side. but we need to count
219 // the number of inner faces in order to balance the remaining
220 // faces properly
221 const CellId id_mine = dcell->id();
222 if (use_active_cells && neighbor->has_children())
223 for (unsigned int c = 0;
224 c < (dcell->has_periodic_neighbor(f) ?
225 dcell->periodic_neighbor(f)
226 ->face(dcell->periodic_neighbor_face_no(f))
227 ->n_children() :
228 dcell->face(f)->n_children());
229 ++c)
230 {
231 typename ::Triangulation<dim>::cell_iterator
232 neighbor_c =
233 dcell->at_boundary(f) ?
234 dcell->periodic_neighbor_child_on_subface(f, c) :
235 dcell->neighbor_child_on_subface(f, c);
236 const types::subdomain_id neigh_domain =
237 neighbor_c->subdomain_id();
238 if (my_domain < neigh_domain)
239 inner_faces_at_proc_boundary[neigh_domain]
240 .n_hanging_faces_larger_subdomain++;
241 else if (my_domain > neigh_domain)
242 inner_faces_at_proc_boundary[neigh_domain]
243 .n_hanging_faces_smaller_subdomain++;
244 }
245 else
246 {
247 const types::subdomain_id neigh_domain =
248 use_active_cells ? neighbor->subdomain_id() :
249 neighbor->level_subdomain_id();
250 if (neighbor->level() < dcell->level() &&
251 use_active_cells)
252 {
253 if (my_domain < neigh_domain)
254 inner_faces_at_proc_boundary[neigh_domain]
255 .n_hanging_faces_smaller_subdomain++;
256 else if (my_domain > neigh_domain)
257 inner_faces_at_proc_boundary[neigh_domain]
258 .n_hanging_faces_larger_subdomain++;
259 }
260 else if (neighbor->level() == dcell->level() &&
261 my_domain != neigh_domain)
262 {
263 // always list the cell whose owner has the lower
264 // subdomain id first. this applies to both processors
265 // involved, so both processors will generate the same
266 // list that we will later order
267 const CellId id_neigh = neighbor->id();
268 if (my_domain < neigh_domain)
269 inner_faces_at_proc_boundary[neigh_domain]
270 .shared_faces.emplace_back(id_mine, id_neigh);
271 else
272 inner_faces_at_proc_boundary[neigh_domain]
273 .shared_faces.emplace_back(id_neigh, id_mine);
274 }
275 }
276 }
277 }
278
279 // sort the cell ids related to each neighboring processor. This
280 // algorithm is symmetric so every processor combination should
281 // arrive here and no deadlock should be possible
282 for (auto &inner_face : inner_faces_at_proc_boundary)
283 {
284 Assert(inner_face.first != my_domain,
285 ExcInternalError("Should not send info to myself"));
286 std::sort(inner_face.second.shared_faces.begin(),
287 inner_face.second.shared_faces.end());
288 inner_face.second.shared_faces.erase(
289 std::unique(inner_face.second.shared_faces.begin(),
290 inner_face.second.shared_faces.end()),
291 inner_face.second.shared_faces.end());
292
293 // safety check: both involved processors should see the same list
294 // because the pattern of ghosting is symmetric. We test this by
295 // looking at the length of the lists of faces
296# if defined(DEAL_II_WITH_MPI) && defined(DEBUG)
297 MPI_Comm comm = MPI_COMM_SELF;
298 if (const ::parallel::TriangulationBase<dim> *ptria =
299 dynamic_cast<const ::parallel::TriangulationBase<dim>
300 *>(&triangulation))
301 comm = ptria->get_mpi_communicator();
302
303 MPI_Status status;
304 unsigned int mysize = inner_face.second.shared_faces.size();
305 unsigned int othersize = numbers::invalid_unsigned_int;
306
307 int ierr = MPI_Sendrecv(&mysize,
308 1,
309 MPI_UNSIGNED,
310 inner_face.first,
311 600 + my_domain,
312 &othersize,
313 1,
314 MPI_UNSIGNED,
315 inner_face.first,
316 600 + inner_face.first,
317 comm,
318 &status);
319 AssertThrowMPI(ierr);
320 AssertDimension(mysize, othersize);
321 mysize = inner_face.second.n_hanging_faces_smaller_subdomain;
322 ierr = MPI_Sendrecv(&mysize,
323 1,
324 MPI_UNSIGNED,
325 inner_face.first,
326 700 + my_domain,
327 &othersize,
328 1,
329 MPI_UNSIGNED,
330 inner_face.first,
331 700 + inner_face.first,
332 comm,
333 &status);
334 AssertThrowMPI(ierr);
335 AssertDimension(mysize, othersize);
336 mysize = inner_face.second.n_hanging_faces_larger_subdomain;
337 ierr = MPI_Sendrecv(&mysize,
338 1,
339 MPI_UNSIGNED,
340 inner_face.first,
341 800 + my_domain,
342 &othersize,
343 1,
344 MPI_UNSIGNED,
345 inner_face.first,
346 800 + inner_face.first,
347 comm,
348 &status);
349 AssertThrowMPI(ierr);
350 AssertDimension(mysize, othersize);
351# endif
352
353 // Arrange the face "ownership" such that cells that are access
354 // by more than one face (think of a cell in a corner) get
355 // ghosted. This arrangement has the advantage that we need to
356 // send less data because the same data is used twice. The
357 // strategy applied here is to ensure the same order of face
358 // pairs on both processors that share some faces, and make the
359 // same decision on both sides.
360
361 // Create a vector with cell ids sorted over the processor with
362 // the larger rank. In the code below we need to be able to
363 // identify the same cell once for the processor with higher
364 // rank and once for the processor with the lower rank. The
365 // format for the processor with the higher rank is already
366 // contained in `shared_faces`, whereas we need a copy that we
367 // sort differently for the other way around.
368 std::vector<std::tuple<CellId, CellId, unsigned int>> other_range(
369 inner_face.second.shared_faces.size());
370 for (unsigned int i = 0; i < other_range.size(); ++i)
371 other_range[i] =
372 std::make_tuple(inner_face.second.shared_faces[i].second,
373 inner_face.second.shared_faces[i].first,
374 i);
375 std::sort(other_range.begin(), other_range.end());
376
377 // the vector 'assignment' sets whether a particular cell
378 // appears more often and acts as a pre-selection of the rank. A
379 // value of 1 means that the process with the higher rank gets
380 // those faces, a value -1 means that the process with the lower
381 // rank gets it, whereas a value 0 means that the decision can
382 // be made in an arbitrary way.
383 unsigned int n_faces_lower_proc = 0, n_faces_higher_proc = 0;
384 std::vector<signed char> assignment(other_range.size(), 0);
385 if (inner_face.second.shared_faces.size() > 0)
386 {
387 // identify faces that go to the processor with the higher
388 // rank
389 unsigned int count = 0;
390 for (unsigned int i = 1;
391 i < inner_face.second.shared_faces.size();
392 ++i)
393 if (inner_face.second.shared_faces[i].first ==
394 inner_face.second.shared_faces[i - 1 - count].first)
395 ++count;
396 else
397 {
398 AssertThrow(count < 2 * dim, ExcInternalError());
399 if (count > 0)
400 {
401 for (unsigned int k = 0; k <= count; ++k)
402 assignment[i - 1 - k] = 1;
403 n_faces_higher_proc += count + 1;
404 }
405 count = 0;
406 }
407
408 // identify faces that definitely go to the processor with
409 // the lower rank - this must use the sorting of CellId
410 // variables from the processor with the higher rank, i.e.,
411 // other_range rather than `shared_faces`.
412 count = 0;
413 for (unsigned int i = 1; i < other_range.size(); ++i)
414 if (std::get<0>(other_range[i]) ==
415 std::get<0>(other_range[i - 1 - count]))
416 ++count;
417 else
418 {
419 AssertThrow(count < 2 * dim, ExcInternalError());
420 if (count > 0)
421 {
422 for (unsigned int k = 0; k <= count; ++k)
423 {
424 Assert(inner_face.second
425 .shared_faces[std::get<2>(
426 other_range[i - 1])]
427 .second ==
428 inner_face.second
429 .shared_faces[std::get<2>(
430 other_range[i - 1 - k])]
431 .second,
433 // only assign to -1 if higher rank was not
434 // yet set
435 if (assignment[std::get<2>(
436 other_range[i - 1 - k])] == 0)
437 {
438 assignment[std::get<2>(
439 other_range[i - 1 - k])] = -1;
440 ++n_faces_lower_proc;
441 }
442 }
443 }
444 count = 0;
445 }
446 }
447
448
449 // divide the faces evenly between the two processors. the
450 // processor with small rank takes the first half, the processor
451 // with larger rank the second half. Adjust for the hanging
452 // faces that always get assigned to one side, and the faces we
453 // have already assigned due to the criterion above
454 n_faces_lower_proc +=
455 inner_face.second.n_hanging_faces_smaller_subdomain;
456 n_faces_higher_proc +=
457 inner_face.second.n_hanging_faces_larger_subdomain;
458 const unsigned int n_total_faces_at_proc_boundary =
459 (inner_face.second.shared_faces.size() +
460 inner_face.second.n_hanging_faces_smaller_subdomain +
461 inner_face.second.n_hanging_faces_larger_subdomain);
462 unsigned int split_index = n_total_faces_at_proc_boundary / 2;
463 if (split_index < n_faces_lower_proc)
464 split_index = 0;
465 else if (split_index <
466 n_total_faces_at_proc_boundary - n_faces_higher_proc)
467 split_index -= n_faces_lower_proc;
468 else
469 split_index = n_total_faces_at_proc_boundary -
470 n_faces_higher_proc - n_faces_lower_proc;
471
472 // make sure the splitting is consistent between both sides
473# if defined(DEAL_II_WITH_MPI) && defined(DEBUG)
474 ierr = MPI_Sendrecv(&split_index,
475 1,
476 MPI_UNSIGNED,
477 inner_face.first,
478 900 + my_domain,
479 &othersize,
480 1,
481 MPI_UNSIGNED,
482 inner_face.first,
483 900 + inner_face.first,
484 comm,
485 &status);
486 AssertThrowMPI(ierr);
487 AssertDimension(split_index, othersize);
488 ierr = MPI_Sendrecv(&n_faces_lower_proc,
489 1,
490 MPI_UNSIGNED,
491 inner_face.first,
492 1000 + my_domain,
493 &othersize,
494 1,
495 MPI_UNSIGNED,
496 inner_face.first,
497 1000 + inner_face.first,
498 comm,
499 &status);
500 AssertThrowMPI(ierr);
501 AssertDimension(n_faces_lower_proc, othersize);
502 ierr = MPI_Sendrecv(&n_faces_higher_proc,
503 1,
504 MPI_UNSIGNED,
505 inner_face.first,
506 1100 + my_domain,
507 &othersize,
508 1,
509 MPI_UNSIGNED,
510 inner_face.first,
511 1100 + inner_face.first,
512 comm,
513 &status);
514 AssertThrowMPI(ierr);
515 AssertDimension(n_faces_higher_proc, othersize);
516# endif
517
518 // collect the faces on both sides
519 std::vector<std::pair<CellId, CellId>> owned_faces_lower,
520 owned_faces_higher;
521 for (unsigned int i = 0; i < assignment.size(); ++i)
522 if (assignment[i] < 0)
523 owned_faces_lower.push_back(
524 inner_face.second.shared_faces[i]);
525 else if (assignment[i] > 0)
526 owned_faces_higher.push_back(
527 inner_face.second.shared_faces[i]);
528 AssertIndexRange(split_index,
529 inner_face.second.shared_faces.size() + 1 -
530 owned_faces_lower.size() -
531 owned_faces_higher.size());
532
533 unsigned int i = 0, c = 0;
534 for (; i < assignment.size() && c < split_index; ++i)
535 if (assignment[i] == 0)
536 {
537 owned_faces_lower.push_back(
538 inner_face.second.shared_faces[i]);
539 ++c;
540 }
541 for (; i < assignment.size(); ++i)
542 if (assignment[i] == 0)
543 {
544 owned_faces_higher.push_back(
545 inner_face.second.shared_faces[i]);
546 }
547
548 if constexpr (running_in_debug_mode())
549 {
550 // check consistency of faces on both sides
551 std::vector<std::pair<CellId, CellId>> check_faces;
552 check_faces.insert(check_faces.end(),
553 owned_faces_lower.begin(),
554 owned_faces_lower.end());
555 check_faces.insert(check_faces.end(),
556 owned_faces_higher.begin(),
557 owned_faces_higher.end());
558 std::sort(check_faces.begin(), check_faces.end());
559 AssertDimension(check_faces.size(),
560 inner_face.second.shared_faces.size());
561 for (unsigned int i = 0; i < check_faces.size(); ++i)
562 Assert(check_faces[i] == inner_face.second.shared_faces[i],
564 }
565
566 // now only set half of the faces as the ones to keep
567 if (my_domain < inner_face.first)
568 inner_face.second.shared_faces.swap(owned_faces_lower);
569 else
570 inner_face.second.shared_faces.swap(owned_faces_higher);
571
572 std::sort(inner_face.second.shared_faces.begin(),
573 inner_face.second.shared_faces.end());
574 }
575 }
576
577 // fill in the additional cells that we need access to via ghosting to
578 // cell_levels
579 std::set<std::pair<unsigned int, unsigned int>> ghost_cells;
580 for (unsigned int i = 0; i < cell_levels.size(); ++i)
581 {
582 typename ::Triangulation<dim>::cell_iterator dcell(
583 &triangulation, cell_levels[i].first, cell_levels[i].second);
584 if (use_active_cells)
585 Assert(dcell->is_active(), ExcNotImplemented());
586 for (const auto f : dcell->face_indices())
587 {
588 if (dcell->at_boundary(f) && !dcell->has_periodic_neighbor(f))
589 face_is_owned[dcell->face(f)->index()] =
590 FaceCategory::locally_active_at_boundary;
591 else if (!build_inner_faces)
592 continue;
593
594 // treat boundaries of cells of different refinement level
595 // inside the domain in case of multigrid separately
596 else if ((dcell->at_boundary(f) == false ||
597 dcell->has_periodic_neighbor(f)) &&
598 mg_level != numbers::invalid_unsigned_int &&
599 dcell->neighbor_or_periodic_neighbor(f)->level() <
600 dcell->level())
601 {
602 face_is_owned[dcell->face(f)->index()] =
603 FaceCategory::multigrid_refinement_edge;
604 }
605 else
606 {
607 typename ::Triangulation<dim>::cell_iterator neighbor =
608 dcell->neighbor_or_periodic_neighbor(f);
609
610 // neighbor is refined -> face will be treated by neighbor
611 if (use_active_cells && neighbor->has_children() &&
612 hold_all_faces_to_owned_cells == false)
613 continue;
614
615 bool add_to_ghost = false;
617 id1 = use_active_cells ? dcell->subdomain_id() :
618 dcell->level_subdomain_id(),
619 id2 = use_active_cells ?
620 (neighbor->has_children() ?
621 dcell->neighbor_child_on_subface(f, 0)
622 ->subdomain_id() :
623 neighbor->subdomain_id()) :
624 neighbor->level_subdomain_id();
625
626 // Check whether the current face should be processed
627 // locally (instead of being processed from the other
628 // side). We process a face locally when we are more refined
629 // (in the active cell case) or when the face is listed in
630 // the `shared_faces` data structure that we built above.
631 if ((id1 == id2 &&
632 (use_active_cells == false || neighbor->is_active())) ||
633 dcell->level() > neighbor->level() ||
634 std::binary_search(
635 inner_faces_at_proc_boundary[id2].shared_faces.begin(),
636 inner_faces_at_proc_boundary[id2].shared_faces.end(),
637 std::make_pair(id1 < id2 ? dcell->id() : neighbor->id(),
638 id1 < id2 ? neighbor->id() :
639 dcell->id())))
640 {
641 face_is_owned[dcell->face(f)->index()] =
642 FaceCategory::locally_active_done_here;
643 if (dcell->level() == neighbor->level())
644 face_is_owned
645 [neighbor
646 ->face(dcell->has_periodic_neighbor(f) ?
647 dcell->periodic_neighbor_face_no(f) :
648 dcell->neighbor_face_no(f))
649 ->index()] =
650 FaceCategory::locally_active_done_here;
651
652 // If neighbor is a ghost element (i.e.
653 // dcell->subdomain_id !
654 // dcell->neighbor(f)->subdomain_id()), we need to add its
655 // index into cell level list.
656 if (use_active_cells)
657 add_to_ghost =
658 (dcell->subdomain_id() != neighbor->subdomain_id());
659 else
660 add_to_ghost = (dcell->level_subdomain_id() !=
661 neighbor->level_subdomain_id());
662 }
663 else if (hold_all_faces_to_owned_cells == true)
664 {
665 // add all cells to ghost layer...
666 face_is_owned[dcell->face(f)->index()] =
667 FaceCategory::ghosted;
668 if (use_active_cells)
669 {
670 if (neighbor->has_children())
671 for (unsigned int s = 0;
672 s < dcell->face(f)->n_children();
673 ++s)
674 if (dcell->at_boundary(f))
675 {
676 if (dcell
677 ->periodic_neighbor_child_on_subface(f,
678 s)
679 ->subdomain_id() !=
680 dcell->subdomain_id())
681 add_to_ghost = true;
682 }
683 else
684 {
685 if (dcell->neighbor_child_on_subface(f, s)
686 ->subdomain_id() !=
687 dcell->subdomain_id())
688 add_to_ghost = true;
689 }
690 else
691 add_to_ghost = (dcell->subdomain_id() !=
692 neighbor->subdomain_id());
693 }
694 else
695 add_to_ghost = (dcell->level_subdomain_id() !=
696 neighbor->level_subdomain_id());
697 }
698
699 if (add_to_ghost)
700 {
701 if (use_active_cells && neighbor->has_children())
702 for (unsigned int s = 0;
703 s < dcell->face(f)->n_children();
704 ++s)
705 {
706 typename ::Triangulation<dim>::cell_iterator
707 neighbor_child =
708 dcell->at_boundary(f) ?
709 dcell->periodic_neighbor_child_on_subface(f,
710 s) :
711 dcell->neighbor_child_on_subface(f, s);
712 if (neighbor_child->subdomain_id() !=
713 dcell->subdomain_id())
714 ghost_cells.insert(
715 std::pair<unsigned int, unsigned int>(
716 neighbor_child->level(),
717 neighbor_child->index()));
718 }
719 else
720 ghost_cells.insert(
721 std::pair<unsigned int, unsigned int>(
722 neighbor->level(), neighbor->index()));
723 at_processor_boundary[i] = true;
724 }
725 }
726 }
727 }
728
729 // step 2: append the ghost cells at the end of the locally owned
730 // cells
731 for (const auto &ghost_cell : ghost_cells)
732 cell_levels.push_back(ghost_cell);
733 }
734
735
736
737 template <int dim>
738 void
739 FaceSetup<dim>::generate_faces(
740 const ::Triangulation<dim> &triangulation,
741 const std::vector<std::pair<unsigned int, unsigned int>> &cell_levels,
742 TaskInfo &task_info)
743 {
744 const bool is_mixed_mesh = triangulation.is_mixed_mesh();
745
746 // step 1: create the inverse map between cell iterators and the
747 // cell_level_index field
748 std::map<std::pair<unsigned int, unsigned int>, unsigned int>
749 map_to_vectorized;
750 for (unsigned int cell = 0; cell < cell_levels.size(); ++cell)
751 if (cell == 0 || cell_levels[cell] != cell_levels[cell - 1])
752 {
753 typename ::Triangulation<dim>::cell_iterator dcell(
755 cell_levels[cell].first,
756 cell_levels[cell].second);
757 std::pair<unsigned int, unsigned int> level_index(dcell->level(),
758 dcell->index());
759 map_to_vectorized[level_index] = cell;
760 }
761
762 // step 2: fill the information about inner faces and boundary faces
763 const unsigned int vectorization_length = task_info.vectorization_length;
764 task_info.face_partition_data.resize(
765 task_info.cell_partition_data.size() - 1, 0);
766 task_info.boundary_partition_data.resize(
767 task_info.cell_partition_data.size() - 1, 0);
768 std::vector<unsigned char> face_visited(face_is_owned.size(), 0);
769 for (unsigned int partition = 0;
770 partition < task_info.cell_partition_data.size() - 2;
771 ++partition)
772 {
773 unsigned int boundary_counter = 0;
774 unsigned int inner_counter = 0;
775 for (unsigned int cell = task_info.cell_partition_data[partition] *
776 vectorization_length;
777 cell < task_info.cell_partition_data[partition + 1] *
778 vectorization_length;
779 ++cell)
780 if (cell == 0 || cell_levels[cell] != cell_levels[cell - 1])
781 {
782 typename ::Triangulation<dim>::cell_iterator dcell(
784 cell_levels[cell].first,
785 cell_levels[cell].second);
786 for (const auto f : dcell->face_indices())
787 {
788 // boundary face
789 if (face_is_owned[dcell->face(f)->index()] ==
790 FaceCategory::locally_active_at_boundary)
791 {
792 Assert(dcell->at_boundary(f), ExcInternalError());
793 ++boundary_counter;
794 FaceToCellTopology<1> info;
795 info.cells_interior[0] = cell;
796 info.cells_exterior[0] = numbers::invalid_unsigned_int;
797 info.interior_face_no = f;
798 info.exterior_face_no = dcell->face(f)->boundary_id();
799 info.face_type =
800 is_mixed_mesh ?
801 (dcell->face(f)->reference_cell() !=
803 0;
804 info.subface_index =
806 info.face_orientation = 0;
807 boundary_faces.push_back(info);
808
809 face_visited[dcell->face(f)->index()]++;
810 }
811 // interior face, including faces over periodic boundaries
812 else
813 {
814 typename ::Triangulation<dim>::cell_iterator
815 neighbor = dcell->neighbor_or_periodic_neighbor(f);
816 if (use_active_cells && neighbor->has_children() &&
817 dim > 1)
818 {
819 // dim > 1 because face()->n_children() = 0 when dim
820 // == 1
821 for (unsigned int c = 0;
822 c < dcell->face(f)->n_children();
823 ++c)
824 {
825 typename ::Triangulation<
826 dim>::cell_iterator neighbor_c =
827 dcell->at_boundary(f) ?
828 dcell->periodic_neighbor_child_on_subface(
829 f, c) :
830 dcell->neighbor_child_on_subface(f, c);
831 const types::subdomain_id neigh_domain =
832 neighbor_c->subdomain_id();
833 const unsigned int neighbor_face_no =
834 dcell->has_periodic_neighbor(f) ?
835 dcell->periodic_neighbor_face_no(f) :
836 dcell->neighbor_face_no(f);
837 if (neigh_domain != dcell->subdomain_id() ||
838 face_visited
839 [dcell->face(f)->child(c)->index()] ==
840 1)
841 {
842 std::pair<unsigned int, unsigned int>
843 level_index(neighbor_c->level(),
844 neighbor_c->index());
845 if (face_is_owned
846 [dcell->face(f)->child(c)->index()] ==
847 FaceCategory::locally_active_done_here)
848 {
849 ++inner_counter;
850 inner_faces.push_back(create_face(
851 neighbor_face_no,
852 neighbor_c,
853 map_to_vectorized[level_index],
854 dcell,
855 cell,
856 is_mixed_mesh));
857 }
858 else if (face_is_owned[dcell->face(f)
859 ->child(c)
860 ->index()] ==
861 FaceCategory::ghosted ||
862 face_is_owned[dcell->face(f)
863 ->index()] ==
864 FaceCategory::ghosted)
865 {
866 inner_ghost_faces.push_back(create_face(
867 neighbor_face_no,
868 neighbor_c,
869 map_to_vectorized[level_index],
870 dcell,
871 cell,
872 is_mixed_mesh));
873 }
874 else
875 Assert(face_is_owned[dcell->face(f)
876 ->index()] ==
877 FaceCategory::
878 locally_active_done_elsewhere,
880 }
881 else
882 {
883 face_visited
884 [dcell->face(f)->child(c)->index()] = 1;
885 }
886 }
887 }
888 else if (dim == 1)
889 {
890 // Follow much the same procedure of dim > 1 with
891 // one large exception: Face is created on first
892 // visitation as long as neighbor has no children
893 // face_visited is used as a flag that a face has
894 // already been created
895 if (face_visited[dcell->face(f)->index()] == 0 &&
896 !(neighbor->has_children()))
897 {
898 std::pair<unsigned int, unsigned int>
899 level_index(neighbor->level(),
900 neighbor->index());
901 if (face_is_owned[dcell->face(f)->index()] ==
902 FaceCategory::locally_active_done_here)
903 {
904 Assert(use_active_cells ||
905 dcell->level() ==
906 neighbor->level(),
908 ++inner_counter;
909 inner_faces.push_back(create_face(
910 f,
911 dcell,
912 cell,
913 neighbor,
914 map_to_vectorized[level_index],
915 is_mixed_mesh));
916 face_visited[dcell->face(f)->index()] = 1;
917 }
918 else if (face_is_owned[dcell->face(f)
919 ->index()] ==
920 FaceCategory::ghosted)
921 {
922 inner_ghost_faces.push_back(create_face(
923 f,
924 dcell,
925 cell,
926 neighbor,
927 map_to_vectorized[level_index],
928 is_mixed_mesh));
929 face_visited[dcell->face(f)->index()] = 1;
930 }
931 }
932 if (face_is_owned[dcell->face(f)->index()] ==
933 FaceCategory::multigrid_refinement_edge)
934 {
935 refinement_edge_faces.push_back(
936 create_face(f,
937 dcell,
938 cell,
939 neighbor,
940 refinement_edge_faces.size(),
941 is_mixed_mesh));
942 }
943 }
944 else
945 {
946 const types::subdomain_id my_domain =
947 use_active_cells ? dcell->subdomain_id() :
948 dcell->level_subdomain_id();
949 const types::subdomain_id neigh_domain =
950 use_active_cells ? neighbor->subdomain_id() :
951 neighbor->level_subdomain_id();
952 if (neigh_domain != my_domain ||
953 face_visited[dcell->face(f)->index()] == 1)
954 {
955 std::pair<unsigned int, unsigned int>
956 level_index(neighbor->level(),
957 neighbor->index());
958 if (face_is_owned[dcell->face(f)->index()] ==
959 FaceCategory::locally_active_done_here)
960 {
961 Assert(use_active_cells ||
962 dcell->level() ==
963 neighbor->level(),
965 ++inner_counter;
966 inner_faces.push_back(create_face(
967 f,
968 dcell,
969 cell,
970 neighbor,
971 map_to_vectorized[level_index],
972 is_mixed_mesh));
973 }
974 else if (face_is_owned[dcell->face(f)
975 ->index()] ==
976 FaceCategory::ghosted)
977 {
978 inner_ghost_faces.push_back(create_face(
979 f,
980 dcell,
981 cell,
982 neighbor,
983 map_to_vectorized[level_index],
984 is_mixed_mesh));
985 }
986 }
987 else
988 {
989 face_visited[dcell->face(f)->index()] = 1;
990 if (dcell->has_periodic_neighbor(f))
991 face_visited
992 [neighbor
993 ->face(
994 dcell->periodic_neighbor_face_no(f))
995 ->index()] = 1;
996 }
997 if (face_is_owned[dcell->face(f)->index()] ==
998 FaceCategory::multigrid_refinement_edge)
999 {
1000 refinement_edge_faces.push_back(
1001 create_face(f,
1002 dcell,
1003 cell,
1004 neighbor,
1005 refinement_edge_faces.size(),
1006 is_mixed_mesh));
1007 }
1008 }
1009 }
1010 }
1011 }
1012 task_info.face_partition_data[partition + 1] =
1013 task_info.face_partition_data[partition] + inner_counter;
1014 task_info.boundary_partition_data[partition + 1] =
1015 task_info.boundary_partition_data[partition] + boundary_counter;
1016 }
1017 task_info.ghost_face_partition_data.resize(2);
1018 task_info.ghost_face_partition_data[0] = 0;
1019 task_info.ghost_face_partition_data[1] = inner_ghost_faces.size();
1020 task_info.refinement_edge_face_partition_data.resize(2);
1021 task_info.refinement_edge_face_partition_data[0] = 0;
1022 task_info.refinement_edge_face_partition_data[1] =
1023 refinement_edge_faces.size();
1024 }
1025
1026
1027
1028 template <int dim>
1029 FaceToCellTopology<1>
1030 FaceSetup<dim>::create_face(
1031 const unsigned int face_no,
1032 const typename ::Triangulation<dim>::cell_iterator &cell,
1033 const unsigned int number_cell_interior,
1034 const typename ::Triangulation<dim>::cell_iterator &neighbor,
1035 const unsigned int number_cell_exterior,
1036 const bool is_mixed_mesh)
1037 {
1038 FaceToCellTopology<1> info;
1039 info.cells_interior[0] = number_cell_interior;
1040 info.cells_exterior[0] = number_cell_exterior;
1041 info.interior_face_no = face_no;
1042 if (cell->has_periodic_neighbor(face_no))
1043 info.exterior_face_no = cell->periodic_neighbor_face_no(face_no);
1044 else
1045 info.exterior_face_no = cell->neighbor_face_no(face_no);
1046
1047 info.face_type = is_mixed_mesh ?
1048 (cell->face(face_no)->reference_cell() !=
1049 ReferenceCells::get_hypercube<dim - 1>()) :
1050 0;
1051
1052 info.subface_index = GeometryInfo<dim>::max_children_per_cell;
1053 Assert(neighbor->level() <= cell->level(), ExcInternalError());
1054
1055 // for dim > 1 and hanging faces we must set a subface index
1056 if (dim > 1 && cell->level() > neighbor->level())
1057 {
1058 if (cell->has_periodic_neighbor(face_no))
1059 info.subface_index =
1060 cell->periodic_neighbor_of_coarser_periodic_neighbor(face_no)
1061 .second;
1062 else
1063 info.subface_index =
1064 cell->neighbor_of_coarser_neighbor(face_no).second;
1065 }
1066
1067 // special treatment of periodic boundaries
1068 if (cell->has_periodic_neighbor(face_no))
1069 {
1070 info.face_orientation = cell->get_triangulation()
1071 .get_periodic_face_map()
1072 .at({cell, face_no})
1073 .second;
1074 }
1075 else
1076 {
1077 const auto interior_face_orientation =
1078 cell->combined_face_orientation(face_no);
1079 const auto exterior_face_orientation =
1080 neighbor->combined_face_orientation(info.exterior_face_no);
1081 if (interior_face_orientation !=
1083 {
1084 info.face_orientation = 8 + interior_face_orientation;
1085 Assert(exterior_face_orientation ==
1087 ExcMessage(
1088 "Face seems to be wrongly oriented from both sides"));
1089 }
1090 else
1091 info.face_orientation = exterior_face_orientation;
1092
1093 // make sure to select correct subface index in case of non-standard
1094 // orientation of the coarser neighbor face
1095 if (cell->level() > neighbor->level() &&
1096 exterior_face_orientation > 0)
1097 {
1098 const Table<2, unsigned int> orientation =
1099 ShapeInfo<double>::compute_orientation_table(2);
1100 const auto face_reference_cell =
1101 cell->face(face_no)->reference_cell();
1102 info.subface_index = orientation(
1103 face_reference_cell.get_inverse_combined_orientation(
1104 exterior_face_orientation),
1105 info.subface_index);
1106 }
1107 }
1108
1109 return info;
1110 }
1111
1112
1113
1120 inline bool
1121 compare_faces_for_vectorization(
1122 const FaceToCellTopology<1> &face1,
1123 const FaceToCellTopology<1> &face2,
1124 const std::vector<unsigned int> &active_fe_indices,
1125 const unsigned int length)
1126 {
1127 if (face1.interior_face_no != face2.interior_face_no)
1128 return false;
1129 if (face1.exterior_face_no != face2.exterior_face_no)
1130 return false;
1131 if (face1.subface_index != face2.subface_index)
1132 return false;
1133 if (face1.face_orientation != face2.face_orientation)
1134 return false;
1135 if (face1.face_type != face2.face_type)
1136 return false;
1137
1138 if (active_fe_indices.size() > 0)
1139 {
1140 if (active_fe_indices[face1.cells_interior[0] / length] !=
1141 active_fe_indices[face2.cells_interior[0] / length])
1142 return false;
1143
1144 if (face2.cells_exterior[0] != numbers::invalid_unsigned_int)
1145 if (active_fe_indices[face1.cells_exterior[0] / length] !=
1146 active_fe_indices[face2.cells_exterior[0] / length])
1147 return false;
1148 }
1149
1150 return true;
1151 }
1152
1153
1154
1161 template <int length>
1162 struct FaceComparator
1163 {
1164 FaceComparator(const std::vector<unsigned int> &active_fe_indices)
1165 : active_fe_indices(active_fe_indices)
1166 {}
1167
1168 bool
1169 operator()(const FaceToCellTopology<length> &face1,
1170 const FaceToCellTopology<length> &face2) const
1171 {
1172 // check if active FE indices match
1173 if (face1.face_type < face2.face_type)
1174 return true;
1175 else if (face1.face_type > face2.face_type)
1176 return false;
1177
1178 // check if active FE indices match
1179 if (active_fe_indices.size() > 0)
1180 {
1181 // ... for interior faces
1182 if (active_fe_indices[face1.cells_interior[0] / length] <
1183 active_fe_indices[face2.cells_interior[0] / length])
1184 return true;
1185 else if (active_fe_indices[face1.cells_interior[0] / length] >
1186 active_fe_indices[face2.cells_interior[0] / length])
1187 return false;
1188
1189 // ... for exterior faces
1190 if (face2.cells_exterior[0] != numbers::invalid_unsigned_int)
1191 {
1192 if (active_fe_indices[face1.cells_exterior[0] / length] <
1193 active_fe_indices[face2.cells_exterior[0] / length])
1194 return true;
1195 else if (active_fe_indices[face1.cells_exterior[0] / length] >
1196 active_fe_indices[face2.cells_exterior[0] / length])
1197 return false;
1198 }
1199 }
1200
1201 for (unsigned int i = 0; i < length; ++i)
1202 if (face1.cells_interior[i] < face2.cells_interior[i])
1203 return true;
1204 else if (face1.cells_interior[i] > face2.cells_interior[i])
1205 return false;
1206 for (unsigned int i = 0; i < length; ++i)
1207 if (face1.cells_exterior[i] < face2.cells_exterior[i])
1208 return true;
1209 else if (face1.cells_exterior[i] > face2.cells_exterior[i])
1210 return false;
1211 if (face1.interior_face_no < face2.interior_face_no)
1212 return true;
1213 else if (face1.interior_face_no > face2.interior_face_no)
1214 return false;
1215 if (face1.exterior_face_no < face2.exterior_face_no)
1216 return true;
1217 else if (face1.exterior_face_no > face2.exterior_face_no)
1218 return false;
1219
1220 // we do not need to check for subface_index and orientation because
1221 // those cannot be different if when all the other values are the
1222 // same.
1223 AssertDimension(face1.subface_index, face2.subface_index);
1224 AssertDimension(face1.face_orientation, face2.face_orientation);
1225
1226 return false;
1227 }
1228
1229 private:
1230 const std::vector<unsigned int> &active_fe_indices;
1231 };
1232
1233
1234
1235 template <int vectorization_width>
1236 void
1238 const std::vector<FaceToCellTopology<1>> &faces_in,
1239 const std::vector<bool> &hard_vectorization_boundary,
1240 std::vector<unsigned int> &face_partition_data,
1241 std::vector<FaceToCellTopology<vectorization_width>> &faces_out,
1242 const std::vector<unsigned int> &active_fe_indices)
1243 {
1244 FaceToCellTopology<vectorization_width> face_batch;
1245 std::vector<std::vector<unsigned int>> faces_type;
1246
1247 unsigned int face_start = face_partition_data[0],
1248 face_end = face_partition_data[0];
1249
1250 face_partition_data[0] = faces_out.size();
1251 for (unsigned int partition = 0;
1252 partition < face_partition_data.size() - 1;
1253 ++partition)
1254 {
1255 std::vector<std::vector<unsigned int>> new_faces_type;
1256
1257 // start with the end point for the last partition
1258 face_start = face_end;
1259 face_end = face_partition_data[partition + 1];
1260
1261 // set the partitioner to the new vectorized lengths
1262 face_partition_data[partition + 1] = face_partition_data[partition];
1263
1264 // loop over the faces in the current partition and reorder according
1265 // to the face type
1266 for (unsigned int face = face_start; face < face_end; ++face)
1267 {
1268 for (auto &face_type : faces_type)
1269 {
1270 // Compare current face with first face of type type
1271 if (compare_faces_for_vectorization(faces_in[face],
1272 faces_in[face_type[0]],
1273 active_fe_indices,
1274 vectorization_width))
1275 {
1276 face_type.push_back(face);
1277 goto face_found;
1278 }
1279 }
1280 faces_type.emplace_back(1, face);
1281 face_found:
1282 {}
1283 }
1284
1285 // insert new faces in sorted list to get good data locality
1286 FaceComparator<vectorization_width> face_comparator(
1287 active_fe_indices);
1288 std::set<FaceToCellTopology<vectorization_width>,
1289 FaceComparator<vectorization_width>>
1290 new_faces(face_comparator);
1291 for (const auto &face_type : faces_type)
1292 {
1293 face_batch.face_type = faces_in[face_type[0]].face_type;
1294 face_batch.interior_face_no =
1295 faces_in[face_type[0]].interior_face_no;
1296 face_batch.exterior_face_no =
1297 faces_in[face_type[0]].exterior_face_no;
1298 face_batch.subface_index = faces_in[face_type[0]].subface_index;
1299 face_batch.face_orientation =
1300 faces_in[face_type[0]].face_orientation;
1301 unsigned int no_faces = face_type.size();
1302 std::vector<unsigned char> touched(no_faces, 0);
1303
1304 // do two passes through the data. The first is to identify
1305 // similar faces within the same index range as the cells which
1306 // will allow for vectorized read operations, the second picks up
1307 // all the rest
1308 unsigned int n_vectorized = 0;
1309 for (unsigned int f = 0; f < no_faces; ++f)
1310 if (faces_in[face_type[f]].cells_interior[0] %
1311 vectorization_width ==
1312 0)
1313 {
1314 bool is_contiguous = true;
1315 if (f + vectorization_width > no_faces)
1316 is_contiguous = false;
1317 else
1318 for (unsigned int v = 1; v < vectorization_width; ++v)
1319 if (faces_in[face_type[f + v]].cells_interior[0] !=
1320 faces_in[face_type[f]].cells_interior[0] + v)
1321 is_contiguous = false;
1322 if (is_contiguous)
1323 {
1325 face_type.size() -
1326 vectorization_width + 1);
1327 for (unsigned int v = 0; v < vectorization_width; ++v)
1328 {
1329 face_batch.cells_interior[v] =
1330 faces_in[face_type[f + v]].cells_interior[0];
1331 face_batch.cells_exterior[v] =
1332 faces_in[face_type[f + v]].cells_exterior[0];
1333 touched[f + v] = 1;
1334 }
1335 new_faces.insert(face_batch);
1336 f += vectorization_width - 1;
1337 n_vectorized += vectorization_width;
1338 }
1339 }
1340
1341 std::vector<unsigned int> untouched;
1342 untouched.reserve(no_faces - n_vectorized);
1343 for (unsigned int f = 0; f < no_faces; ++f)
1344 if (touched[f] == 0)
1345 untouched.push_back(f);
1346 unsigned int v = 0;
1347 for (const auto f : untouched)
1348 {
1349 face_batch.cells_interior[v] =
1350 faces_in[face_type[f]].cells_interior[0];
1351 face_batch.cells_exterior[v] =
1352 faces_in[face_type[f]].cells_exterior[0];
1353 ++v;
1354 if (v == vectorization_width)
1355 {
1356 new_faces.insert(face_batch);
1357 v = 0;
1358 }
1359 }
1360 if (v > 0 && v < vectorization_width)
1361 {
1362 // must add non-filled face
1363 if (hard_vectorization_boundary[partition + 1] ||
1364 partition == face_partition_data.size() - 2)
1365 {
1366 for (; v < vectorization_width; ++v)
1367 {
1368 // Dummy cell, not used
1369 face_batch.cells_interior[v] =
1371 face_batch.cells_exterior[v] =
1373 }
1374 new_faces.insert(face_batch);
1375 }
1376 else
1377 {
1378 // postpone to the next partition
1379 std::vector<unsigned int> untreated(v);
1380 for (unsigned int f = 0; f < v; ++f)
1381 untreated[f] = face_type[*(untouched.end() - 1 - f)];
1382 new_faces_type.push_back(untreated);
1383 }
1384 }
1385 }
1386
1387 // insert sorted list to vector of faces
1388 for (auto it = new_faces.begin(); it != new_faces.end(); ++it)
1389 faces_out.push_back(*it);
1390 face_partition_data[partition + 1] += new_faces.size();
1391
1392 // set the faces that were left over to faces_type for the next round
1393 faces_type = std::move(new_faces_type);
1394 }
1395
1396 if constexpr (running_in_debug_mode())
1397 {
1398 // final safety checks
1399 for (const auto &face_type : faces_type)
1400 AssertDimension(face_type.size(), 0U);
1401
1402 AssertDimension(faces_out.size(), face_partition_data.back());
1403 unsigned int nfaces = 0;
1404 for (unsigned int i = face_partition_data[0];
1405 i < face_partition_data.back();
1406 ++i)
1407 for (unsigned int v = 0; v < vectorization_width; ++v)
1408 nfaces += (faces_out[i].cells_interior[v] !=
1410 AssertDimension(nfaces, faces_in.size());
1411
1412 std::vector<std::pair<unsigned int, unsigned int>> in_faces,
1413 out_faces;
1414 for (const auto &face_in : faces_in)
1415 in_faces.emplace_back(face_in.cells_interior[0],
1416 face_in.cells_exterior[0]);
1417 for (unsigned int i = face_partition_data[0];
1418 i < face_partition_data.back();
1419 ++i)
1420 for (unsigned int v = 0;
1421 v < vectorization_width && faces_out[i].cells_interior[v] !=
1423 ++v)
1424 out_faces.emplace_back(faces_out[i].cells_interior[v],
1425 faces_out[i].cells_exterior[v]);
1426 std::sort(in_faces.begin(), in_faces.end());
1427 std::sort(out_faces.begin(), out_faces.end());
1428 AssertDimension(in_faces.size(), out_faces.size());
1429 for (unsigned int i = 0; i < in_faces.size(); ++i)
1430 {
1431 AssertDimension(in_faces[i].first, out_faces[i].first);
1432 AssertDimension(in_faces[i].second, out_faces[i].second);
1433 }
1434 }
1435 }
1436
1437#endif // ifndef DOXYGEN
1438
1439 } // namespace MatrixFreeFunctions
1440} // namespace internal
1441
1442
1444
1445#endif
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:35
constexpr bool running_in_debug_mode()
Definition config.h:73
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:36
Point< 2 > second
Definition grid_out.cc:4633
Point< 2 > first
Definition grid_out.cc:4632
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)
std::size_t size
Definition mpi.cc:745
const MPI_Comm comm
Definition mpi.cc:924
constexpr char U
constexpr const ReferenceCell & get_hypercube()
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
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)
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
constexpr types::subdomain_id invalid_subdomain_id
Definition types.h:381
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector< std::pair< CellId, CellId > > shared_faces
std::vector< FaceToCellTopology< 1 > > inner_faces
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)
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)