17#ifndef dealii_mesh_worker_mesh_loop_h
18#define dealii_mesh_worker_mesh_loop_h
54 template <
class CellIteratorType>
60 using type = CellIteratorType;
70 template <
class CellIteratorType>
90 template <
class CellIteratorType>
121 std::function<void(
const CellIteratorBaseType &,
131 std::function<void(
const CellIteratorBaseType &,
134 const CellIteratorBaseType &,
276 template <
class CellIteratorType,
279 class CellIteratorBaseType =
284 const CellIteratorType &begin,
285 const CellIteratorType &end,
288 const CopierType & copier,
300 const unsigned int chunk_size = 8
302 const CellIteratorType & begin,
317 std::function<
void(
const CellIteratorBaseType &,
321 std::function<
void(
const CellIteratorBaseType &,
327 std::function<
void(
const CellIteratorBaseType &,
330 const CellIteratorBaseType &,
335 std::function<
void(
const CellIteratorBaseType &,
338 const CellIteratorBaseType &,
345 const unsigned int chunk_size = 8
352 "If you provide a cell worker function, you also need to request "
353 "that work should be done on cells by setting the 'work_on_cells' flag. "
354 "Conversely, if you don't provide a cell worker function, you "
355 "cannot set the 'work_on_cells' flag. One of these two "
356 "conditions is not satisfied."));
363 "If you provide a face worker function, you also need to request "
364 "that work should be done on interior faces by setting either the "
365 "'assemble_own_interior_faces_once' flag or the "
366 "'assemble_own_interior_faces_both' flag. "
367 "Conversely, if you don't provide a face worker function, you "
368 "cannot set either of these two flags. One of these two "
369 "conditions is not satisfied."));
374 "You can only 'specify assemble_ghost_faces_once' "
375 "OR 'assemble_ghost_faces_both', but not both of these flags."));
381 "The option 'cells_after_faces' only makes sense if you assemble on cells."));
386 "If you provide a face worker function, you also need to request "
387 "that work should be done on faces by setting the 'work_on_faces' flag. "
388 "Conversely, if you don't provide a face worker function, you "
389 "cannot set the 'work_on_faces' flag. One of these two "
390 "conditions is not satisfied."));
395 "If you provide a boundary face worker function, you also need to request "
396 "that work should be done on boundary faces by setting the 'assemble_boundary_faces' flag. "
397 "Conversely, if you don't provide a boundary face worker function, you "
398 "cannot set the 'assemble_boundary_faces' flag. One of these two "
399 "conditions is not satisfied."));
401 auto cell_action = [&](
const CellIteratorBaseType &cell,
406 copy = sample_copy_data;
409 const auto dim = cell->get_triangulation().dimension;
411 const bool ignore_subdomain =
412 (cell->get_triangulation().locally_owned_subdomain() ==
416 (cell->is_level_cell() ? cell->level_subdomain_id() :
417 cell->subdomain_id());
419 const bool own_cell =
421 (current_subdomain_id ==
422 cell->get_triangulation().locally_owned_subdomain());
424 if ((!ignore_subdomain) &&
431 cell_worker(cell, scratch, copy);
434 for (
const unsigned int face_no : cell->face_indices())
436 if (cell->at_boundary(face_no) &&
437 !cell->has_periodic_neighbor(face_no))
441 boundary_worker(cell, face_no, scratch, copy);
447 neighbor = cell->neighbor_or_periodic_neighbor(face_no);
451 if (neighbor->is_level_cell())
452 neighbor_subdomain_id = neighbor->level_subdomain_id();
454 else if (neighbor->is_active())
455 neighbor_subdomain_id = neighbor->subdomain_id();
457 const bool own_neighbor =
459 (neighbor_subdomain_id ==
460 cell->get_triangulation().locally_owned_subdomain());
463 if (!own_cell && !own_neighbor)
467 if (own_cell && own_neighbor &&
473 if (own_cell != own_neighbor &&
481 const bool periodic_neighbor =
482 cell->has_periodic_neighbor(face_no);
484 if (dim > 1 && ((!periodic_neighbor &&
485 cell->neighbor_is_coarser(face_no) &&
486 neighbor->is_active()) ||
487 (periodic_neighbor &&
488 cell->periodic_neighbor_is_coarser(face_no) &&
489 neighbor->is_active())))
498 const std::pair<unsigned int, unsigned int>
501 cell->periodic_neighbor_of_coarser_periodic_neighbor(
503 cell->neighbor_of_coarser_neighbor(face_no);
509 neighbor_face_no.first,
510 neighbor_face_no.second,
520 face_worker(neighbor,
521 neighbor_face_no.first,
522 neighbor_face_no.second,
530 else if (dim == 1 && cell->level() > neighbor->level())
533 const unsigned int neighbor_face_no =
535 cell->periodic_neighbor_face_no(face_no) :
536 cell->neighbor_face_no(face_no);
537 Assert(periodic_neighbor ||
538 neighbor->face(neighbor_face_no) ==
555 face_worker(neighbor,
570 neighbor->has_children())
575 Assert((!periodic_neighbor &&
576 !cell->neighbor_is_coarser(face_no)) ||
577 (periodic_neighbor &&
578 !cell->periodic_neighbor_is_coarser(face_no)),
584 if (own_cell && own_neighbor &&
597 if (own_cell && !own_neighbor &&
599 (neighbor_subdomain_id < current_subdomain_id))
602 const unsigned int neighbor_face_no =
604 cell->periodic_neighbor_face_no(face_no) :
605 cell->neighbor_face_no(face_no);
606 Assert(periodic_neighbor ||
607 neighbor->face(neighbor_face_no) ==
627 cell_worker(cell, scratch, copy);
710 template <
class CellIteratorType,
713 class CellIteratorBaseType =
730 std::function<
void(
const CellIteratorBaseType &,
734 std::function<
void(
const CellIteratorBaseType &,
740 std::function<
void(
const CellIteratorBaseType &,
743 const CellIteratorBaseType &,
748 std::function<
void(
const CellIteratorBaseType &,
751 const CellIteratorBaseType &,
758 const unsigned int chunk_size = 8)
761 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
764 CellIteratorBaseType>(iterator_range.
begin(),
765 iterator_range.
end(),
836 template <
class CellIteratorType,
843 MainClass & main_class,
844 void (MainClass::*cell_worker)(const CellIteratorType &,
847 void (MainClass::*copier)(const
CopyData &),
851 void (MainClass::*boundary_worker)(const CellIteratorType &,
855 void (MainClass::*face_worker)(const CellIteratorType &,
858 const CellIteratorType &,
864 const unsigned
int chunk_size = 8)
873 std::function<void(
const CellIteratorType &,
876 const CellIteratorType &,
883 if (cell_worker !=
nullptr)
884 f_cell_worker = [&main_class,
885 cell_worker](
const CellIteratorType &cell_iterator,
888 (main_class.*cell_worker)(cell_iterator, scratch_data, copy_data);
891 if (boundary_worker !=
nullptr)
893 [&main_class, boundary_worker](
const CellIteratorType &cell_iterator,
894 const unsigned int face_no,
898 boundary_worker)(cell_iterator, face_no, scratch_data, copy_data);
901 if (face_worker !=
nullptr)
902 f_face_worker = [&main_class,
903 face_worker](
const CellIteratorType &cell_iterator_1,
904 const unsigned int face_index_1,
905 const unsigned int subface_index_1,
906 const CellIteratorType &cell_iterator_2,
907 const unsigned int face_index_2,
908 const unsigned int subface_index_2,
911 (main_class.*face_worker)(cell_iterator_1,
925 [&main_class, copier](
const CopyData ©_data) {
926 (main_class.*copier)(copy_data);
1016 template <
class CellIteratorType,
1020 class CellIteratorBaseType =
1024 MainClass & main_class,
1025 void (MainClass::*cell_worker)(const CellIteratorBaseType &,
1028 void (MainClass::*copier)(const
CopyData &),
1032 void (MainClass::*boundary_worker)(const CellIteratorBaseType &,
1036 void (MainClass::*face_worker)(const CellIteratorBaseType &,
1039 const CellIteratorBaseType &,
1045 const unsigned
int chunk_size = 8)
1048 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
1052 CellIteratorBaseType>(iterator_range.
begin(),
1053 iterator_range.
end(),
1057 sample_scratch_data,
IteratorOverIterators end() const
IteratorOverIterators begin()
static unsigned int n_threads()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
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)
void cell_action(ITERATOR 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)
std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)> BoundaryWorkerFunctionType
std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)> FaceWorkerFunctionType
std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)> CellWorkerFunctionType
std::function< void(const CopyData &)> CopierFunctionType
@ 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 &)
const types::subdomain_id artificial_subdomain_id
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
typename type_identity< T >::type type_identity_t
typename CellIteratorBaseType< CellIteratorType >::type type
typename CellIteratorBaseType< CellIteratorType >::type type