deal.II version GIT relicensing-1989-gd7a2c90e4e 2024-10-14 01: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\}}\)
Loading...
Searching...
No Matches
mesh_loop.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 2023 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_mesh_worker_mesh_loop_h
17#define dealii_mesh_worker_mesh_loop_h
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/types.h>
24
26#include <deal.II/grid/tria.h>
27
33
34#include <functional>
35#include <type_traits>
36
38
39// Forward declaration
40#ifndef DOXYGEN
41template <typename>
43#endif
44
45namespace MeshWorker
46{
47 namespace internal
48 {
53 template <typename CellIteratorType>
55 {
59 using type = CellIteratorType;
60 };
61
69 template <typename CellIteratorType>
71 {
75 // Since we can have filtered iterators and the like as template
76 // arguments, we recursivelyremove the template layers to retrieve the
77 // underlying iterator type.
79 };
80
89 template <typename CellIteratorType>
90 struct CellIteratorBaseType<FilteredIterator<CellIteratorType>>
91 {
95 // Since we can have nested filtered iterators, we recursively
96 // remove the template layers to retrieve the underlying iterator type.
98 };
99 } // namespace internal
100
101#ifdef DOXYGEN
106 using CellWorkerFunctionType = std::function<
107 void(const CellIteratorBaseType &, ScratchData &, CopyData &)>;
108
113 using CopierFunctionType = std::function<void(const CopyData &)>;
114
120 std::function<void(const CellIteratorBaseType &,
121 const unsigned int,
122 ScratchData &,
123 CopyData &)>;
124
130 std::function<void(const CellIteratorBaseType &,
131 const unsigned int,
132 const unsigned int,
133 const CellIteratorBaseType &,
134 const unsigned int,
135 const unsigned int,
136 ScratchData &,
137 CopyData &)>;
138#endif
139
275 template <typename CellIteratorType,
276 class ScratchData,
277 class CopyData,
278 typename CellIteratorBaseType =
280 void
282#ifdef DOXYGEN
283 const CellIteratorType &begin,
284 const CellIteratorType &end,
285
286 const CellWorkerFunctionType &cell_worker,
287 const CopierType &copier,
288
289 const ScratchData &sample_scratch_data,
290 const CopyData &sample_copy_data,
291
292 const AssembleFlags flags = assemble_own_cells,
293
294 const BoundaryWorkerFunctionType &boundary_worker =
296
297 const FaceWorkerFunctionType &face_worker = FaceWorkerFunctionType(),
298 const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
299 const unsigned int chunk_size = 8
300#else
301 const CellIteratorType &begin,
303
304 const std_cxx20::type_identity_t<std::function<
305 void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>
306 &cell_worker,
307 const std_cxx20::type_identity_t<std::function<void(const CopyData &)>>
308 &copier,
309
310 const ScratchData &sample_scratch_data,
311 const CopyData &sample_copy_data,
312
313 const AssembleFlags flags = assemble_own_cells,
314
316 std::function<void(const CellIteratorBaseType &,
317 const unsigned int,
318 ScratchData &,
319 CopyData &)>> &boundary_worker =
320 std::function<void(const CellIteratorBaseType &,
321 const unsigned int,
322 ScratchData &,
323 CopyData &)>(),
324
326 std::function<void(const CellIteratorBaseType &,
327 const unsigned int,
328 const unsigned int,
329 const CellIteratorBaseType &,
330 const unsigned int,
331 const unsigned int,
332 ScratchData &,
333 CopyData &)>> &face_worker =
334 std::function<void(const CellIteratorBaseType &,
335 const unsigned int,
336 const unsigned int,
337 const CellIteratorBaseType &,
338 const unsigned int,
339 const unsigned int,
340 ScratchData &,
341 CopyData &)>(),
342
343 const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
344 const unsigned int chunk_size = 8
345#endif
346 )
347 {
348 Assert(
349 (!cell_worker) == !(flags & work_on_cells),
351 "If you provide a cell worker function, you also need to request "
352 "that work should be done on cells by setting the 'work_on_cells' flag. "
353 "Conversely, if you don't provide a cell worker function, you "
354 "cannot set the 'work_on_cells' flag. One of these two "
355 "conditions is not satisfied."));
356
362 "If you provide a face worker function, you also need to request "
363 "that work should be done on interior faces by setting either the "
364 "'assemble_own_interior_faces_once' flag or the "
365 "'assemble_own_interior_faces_both' flag. "
366 "Conversely, if you don't provide a face worker function, you "
367 "cannot set either of these two flags. One of these two "
368 "conditions is not satisfied."));
369
373 "You can only 'specify assemble_ghost_faces_once' "
374 "OR 'assemble_ghost_faces_both', but not both of these flags."));
375
376 Assert(
377 !(flags & cells_after_faces) ||
380 "The option 'cells_after_faces' only makes sense if you assemble on cells."));
381
382 Assert(
383 (!face_worker) == !(flags & work_on_faces),
385 "If you provide a face worker function, you also need to request "
386 "that work should be done on faces by setting the 'work_on_faces' flag. "
387 "Conversely, if you don't provide a face worker function, you "
388 "cannot set the 'work_on_faces' flag. One of these two "
389 "conditions is not satisfied."));
390
391 Assert(
392 (!boundary_worker) == !(flags & assemble_boundary_faces),
394 "If you provide a boundary face worker function, you also need to request "
395 "that work should be done on boundary faces by setting the 'assemble_boundary_faces' flag. "
396 "Conversely, if you don't provide a boundary face worker function, you "
397 "cannot set the 'assemble_boundary_faces' flag. One of these two "
398 "conditions is not satisfied."));
399
400 auto cell_action = [&](const CellIteratorBaseType &cell,
401 ScratchData &scratch,
402 CopyData &copy) {
403 // First reset the CopyData class to the empty copy_data given by the
404 // user.
405 copy = sample_copy_data;
406
407 // Store the dimension in which we are working for later use
408 const auto dim = cell->get_triangulation().dimension;
409
410 const bool ignore_subdomain =
411 (cell->get_triangulation().locally_owned_subdomain() ==
413
414 types::subdomain_id current_subdomain_id =
415 (cell->is_level_cell() ? cell->level_subdomain_id() :
416 cell->subdomain_id());
417
418 const bool own_cell =
419 ignore_subdomain ||
420 (current_subdomain_id ==
421 cell->get_triangulation().locally_owned_subdomain());
422
423 if ((!ignore_subdomain) &&
424 (current_subdomain_id == numbers::artificial_subdomain_id))
425 return;
426
427 if (!(flags & (cells_after_faces)) &&
428 (((flags & (assemble_own_cells)) && own_cell) ||
429 ((flags & assemble_ghost_cells) && !own_cell)))
430 cell_worker(cell, scratch, copy);
431
432 if (flags & (work_on_faces | work_on_boundary))
433 for (const unsigned int face_no : cell->face_indices())
434 {
435 if (cell->at_boundary(face_no) &&
436 !cell->has_periodic_neighbor(face_no))
437 {
438 // only integrate boundary faces of own cells
439 if ((flags & assemble_boundary_faces) && own_cell)
440 boundary_worker(cell, face_no, scratch, copy);
441 }
442 else
443 {
444 // interior face, potentially assemble
446 neighbor = cell->neighbor_or_periodic_neighbor(face_no);
447
448 types::subdomain_id neighbor_subdomain_id =
450 if (neighbor->is_level_cell())
451 neighbor_subdomain_id = neighbor->level_subdomain_id();
452 // subdomain id is only valid for active cells
453 else if (neighbor->is_active())
454 neighbor_subdomain_id = neighbor->subdomain_id();
455
456 const bool own_neighbor =
457 ignore_subdomain ||
458 (neighbor_subdomain_id ==
459 cell->get_triangulation().locally_owned_subdomain());
460
461 // skip all faces between two ghost cells
462 if (!own_cell && !own_neighbor)
463 continue;
464
465 // skip if the user doesn't want faces between own cells
466 if (own_cell && own_neighbor &&
469 continue;
470
471 // skip face to ghost
472 if (own_cell != own_neighbor &&
473 !(flags &
475 continue;
476
477 // Deal with refinement edges from the refined side. Assuming
478 // one-irregular meshes, this situation should only occur if
479 // both cells are active.
480 const bool periodic_neighbor =
481 cell->has_periodic_neighbor(face_no);
482
483 if (dim > 1 && ((!periodic_neighbor &&
484 cell->neighbor_is_coarser(face_no) &&
485 neighbor->is_active()) ||
486 (periodic_neighbor &&
487 cell->periodic_neighbor_is_coarser(face_no) &&
488 neighbor->is_active())))
489 {
490 Assert(cell->is_active(), ExcInternalError());
491
492 // skip if only one processor needs to assemble the face
493 // to a ghost cell and the fine cell is not ours.
494 if (!own_cell && (flags & assemble_ghost_faces_once))
495 continue;
496
497 const std::pair<unsigned int, unsigned int>
498 neighbor_face_no =
499 periodic_neighbor ?
500 cell->periodic_neighbor_of_coarser_periodic_neighbor(
501 face_no) :
502 cell->neighbor_of_coarser_neighbor(face_no);
503
504 face_worker(cell,
505 face_no,
507 neighbor,
508 neighbor_face_no.first,
509 neighbor_face_no.second,
510 scratch,
511 copy);
512
514 {
515 // If own faces are to be assembled from both sides,
516 // call the faceworker again with swapped arguments.
517 // This is because we won't be looking at an adaptively
518 // refined edge coming from the other side.
519 face_worker(neighbor,
520 neighbor_face_no.first,
521 neighbor_face_no.second,
522 cell,
523 face_no,
525 scratch,
526 copy);
527 }
528 }
529 else if (dim == 1 && cell->level() > neighbor->level())
530 {
531 // In one dimension, there is no other check to do
532 const unsigned int neighbor_face_no =
533 periodic_neighbor ?
534 cell->periodic_neighbor_face_no(face_no) :
535 cell->neighbor_face_no(face_no);
536 Assert(periodic_neighbor ||
537 neighbor->face(neighbor_face_no) ==
538 cell->face(face_no),
540
541 face_worker(cell,
542 face_no,
544 neighbor,
545 neighbor_face_no,
547 scratch,
548 copy);
549
551 {
552 // If own faces are to be assembled from both sides,
553 // call the faceworker again with swapped arguments.
554 face_worker(neighbor,
555 neighbor_face_no,
557 cell,
558 face_no,
560 scratch,
561 copy);
562 }
563 }
564 else
565 {
566 // If iterator is active and neighbor is refined, skip
567 // internal face.
569 neighbor->has_children())
570 continue;
571
572 // Now neighbor is on the same refinement level.
573 // Double check.
574 Assert((!periodic_neighbor &&
575 !cell->neighbor_is_coarser(face_no)) ||
576 (periodic_neighbor &&
577 !cell->periodic_neighbor_is_coarser(face_no)),
579
580 // If we own both cells only do faces from one side (unless
581 // AssembleFlags says otherwise). Here, we rely on cell
582 // comparison that will look at cell->index().
583 if (own_cell && own_neighbor &&
585 (neighbor < cell))
586 continue;
587
588 // We only look at faces to ghost on the same level once
589 // (only where own_cell=true and own_neighbor=false)
590 if (!own_cell)
591 continue;
592
593 // now only one processor assembles faces_to_ghost. We let
594 // the processor with the smaller (level-)subdomain id
595 // assemble the face.
596 if (own_cell && !own_neighbor &&
597 (flags & assemble_ghost_faces_once) &&
598 (neighbor_subdomain_id < current_subdomain_id))
599 continue;
600
601 const unsigned int neighbor_face_no =
602 periodic_neighbor ?
603 cell->periodic_neighbor_face_no(face_no) :
604 cell->neighbor_face_no(face_no);
605 Assert(periodic_neighbor ||
606 neighbor->face(neighbor_face_no) ==
607 cell->face(face_no),
609
610 face_worker(cell,
611 face_no,
613 neighbor,
614 neighbor_face_no,
616 scratch,
617 copy);
618 }
619 }
620 } // faces
621
622 // Execute the cell_worker if faces are handled before cells
623 if ((flags & cells_after_faces) &&
624 (((flags & assemble_own_cells) && own_cell) ||
625 ((flags & assemble_ghost_cells) && !own_cell)))
626 cell_worker(cell, scratch, copy);
627 };
628
629 // Submit to workstream
630 WorkStream::run(begin,
631 end,
633 copier,
634 sample_scratch_data,
635 sample_copy_data,
636 queue_length,
637 chunk_size);
638 }
639
709 template <typename CellIteratorType,
710 class ScratchData,
711 class CopyData,
712 typename CellIteratorBaseType =
714 void
716 IteratorRange<CellIteratorType> iterator_range,
717 const std_cxx20::type_identity_t<std::function<
718 void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>
719 &cell_worker,
720 const std_cxx20::type_identity_t<std::function<void(const CopyData &)>>
721 &copier,
722
723 const ScratchData &sample_scratch_data,
724 const CopyData &sample_copy_data,
725
726 const AssembleFlags flags = assemble_own_cells,
727
729 std::function<void(const CellIteratorBaseType &,
730 const unsigned int,
731 ScratchData &,
732 CopyData &)>> &boundary_worker =
733 std::function<void(const CellIteratorBaseType &,
734 const unsigned int,
735 ScratchData &,
736 CopyData &)>(),
737
739 std::function<void(const CellIteratorBaseType &,
740 const unsigned int,
741 const unsigned int,
742 const CellIteratorBaseType &,
743 const unsigned int,
744 const unsigned int,
745 ScratchData &,
746 CopyData &)>> &face_worker =
747 std::function<void(const CellIteratorBaseType &,
748 const unsigned int,
749 const unsigned int,
750 const CellIteratorBaseType &,
751 const unsigned int,
752 const unsigned int,
753 ScratchData &,
754 CopyData &)>(),
755
756 const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
757 const unsigned int chunk_size = 8)
758 {
759 // Call the function above
760 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
762 CopyData,
763 CellIteratorBaseType>(iterator_range.begin(),
764 iterator_range.end(),
765 cell_worker,
766 copier,
767 sample_scratch_data,
768 sample_copy_data,
769 flags,
770 boundary_worker,
771 face_worker,
772 queue_length,
773 chunk_size);
774 }
775
835 template <typename CellIteratorType,
836 class ScratchData,
837 class CopyData,
838 class MainClass>
839 void
840 mesh_loop(const CellIteratorType &begin,
842 MainClass &main_class,
843 void (MainClass::*cell_worker)(const CellIteratorType &,
844 ScratchData &,
845 CopyData &),
846 void (MainClass::*copier)(const CopyData &),
847 const ScratchData &sample_scratch_data,
848 const CopyData &sample_copy_data,
849 const AssembleFlags flags = assemble_own_cells,
850 void (MainClass::*boundary_worker)(const CellIteratorType &,
851 const unsigned int,
852 ScratchData &,
853 CopyData &) = nullptr,
854 void (MainClass::*face_worker)(const CellIteratorType &,
855 const unsigned int,
856 const unsigned int,
857 const CellIteratorType &,
858 const unsigned int,
859 const unsigned int,
860 ScratchData &,
861 CopyData &) = nullptr,
862 const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
863 const unsigned int chunk_size = 8)
864 {
865 std::function<void(const CellIteratorType &, ScratchData &, CopyData &)>
866 f_cell_worker;
867
868 std::function<void(
869 const CellIteratorType &, const unsigned int, ScratchData &, CopyData &)>
870 f_boundary_worker;
871
872 std::function<void(const CellIteratorType &,
873 const unsigned int,
874 const unsigned int,
875 const CellIteratorType &,
876 const unsigned int,
877 const unsigned int,
878 ScratchData &,
879 CopyData &)>
880 f_face_worker;
881
882 if (cell_worker != nullptr)
883 f_cell_worker = [&main_class,
884 cell_worker](const CellIteratorType &cell_iterator,
885 ScratchData &scratch_data,
886 CopyData &copy_data) {
887 (main_class.*cell_worker)(cell_iterator, scratch_data, copy_data);
888 };
889
890 if (boundary_worker != nullptr)
891 f_boundary_worker =
892 [&main_class, boundary_worker](const CellIteratorType &cell_iterator,
893 const unsigned int face_no,
894 ScratchData &scratch_data,
895 CopyData &copy_data) {
896 (main_class.*
897 boundary_worker)(cell_iterator, face_no, scratch_data, copy_data);
898 };
899
900 if (face_worker != nullptr)
901 f_face_worker = [&main_class,
902 face_worker](const CellIteratorType &cell_iterator_1,
903 const unsigned int face_index_1,
904 const unsigned int subface_index_1,
905 const CellIteratorType &cell_iterator_2,
906 const unsigned int face_index_2,
907 const unsigned int subface_index_2,
908 ScratchData &scratch_data,
909 CopyData &copy_data) {
910 (main_class.*face_worker)(cell_iterator_1,
911 face_index_1,
912 subface_index_1,
913 cell_iterator_2,
914 face_index_2,
915 subface_index_2,
916 scratch_data,
917 copy_data);
918 };
919
920 mesh_loop(
921 begin,
922 end,
923 f_cell_worker,
924 [&main_class, copier](const CopyData &copy_data) {
925 (main_class.*copier)(copy_data);
926 },
927 sample_scratch_data,
928 sample_copy_data,
929 flags,
930 f_boundary_worker,
931 f_face_worker,
932 queue_length,
933 chunk_size);
934 }
935
1015 template <typename CellIteratorType,
1016 class ScratchData,
1017 class CopyData,
1018 class MainClass,
1019 typename CellIteratorBaseType =
1021 void
1023 MainClass &main_class,
1024 void (MainClass::*cell_worker)(const CellIteratorBaseType &,
1025 ScratchData &,
1026 CopyData &),
1027 void (MainClass::*copier)(const CopyData &),
1028 const ScratchData &sample_scratch_data,
1029 const CopyData &sample_copy_data,
1030 const AssembleFlags flags = assemble_own_cells,
1031 void (MainClass::*boundary_worker)(const CellIteratorBaseType &,
1032 const unsigned int,
1033 ScratchData &,
1034 CopyData &) = nullptr,
1035 void (MainClass::*face_worker)(const CellIteratorBaseType &,
1036 const unsigned int,
1037 const unsigned int,
1038 const CellIteratorBaseType &,
1039 const unsigned int,
1040 const unsigned int,
1041 ScratchData &,
1042 CopyData &) = nullptr,
1043 const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1044 const unsigned int chunk_size = 8)
1045 {
1046 // Call the function above
1047 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
1049 CopyData,
1050 MainClass,
1051 CellIteratorBaseType>(iterator_range.begin(),
1052 iterator_range.end(),
1053 main_class,
1054 cell_worker,
1055 copier,
1056 sample_scratch_data,
1057 sample_copy_data,
1058 flags,
1059 boundary_worker,
1060 face_worker,
1061 queue_length,
1062 chunk_size);
1063 }
1064} // namespace MeshWorker
1065
1067
1068#endif
IteratorOverIterators end() const
IteratorOverIterators begin()
static unsigned int n_threads()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void cell_action(IteratorType cell, DoFInfoBox< dim, DOFINFO > &dof_info, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, const LoopControl &loop_control)
Definition loop.h:316
void mesh_loop(const CellIteratorType &begin, const CellIteratorType &end, const CellWorkerFunctionType &cell_worker, const CopierType &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const BoundaryWorkerFunctionType &boundary_worker=BoundaryWorkerFunctionType(), const FaceWorkerFunctionType &face_worker=FaceWorkerFunctionType(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition mesh_loop.h:281
std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)> BoundaryWorkerFunctionType
Definition mesh_loop.h:123
std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)> FaceWorkerFunctionType
Definition mesh_loop.h:137
std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)> CellWorkerFunctionType
Definition mesh_loop.h:107
std::function< void(const CopyData &)> CopierFunctionType
Definition mesh_loop.h:113
@ assemble_boundary_faces
@ assemble_own_interior_faces_once
@ assemble_ghost_faces_both
@ assemble_own_interior_faces_both
@ assemble_ghost_faces_once
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
bool is_active_iterator(const DI &)
Definition loop.h:48
const types::subdomain_id artificial_subdomain_id
Definition types.h:362
const types::subdomain_id invalid_subdomain_id
Definition types.h:341
static const unsigned int invalid_unsigned_int
Definition types.h:220
typename type_identity< T >::type type_identity_t
Definition type_traits.h:95
typename CellIteratorBaseType< CellIteratorType >::type type
Definition mesh_loop.h:97