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,
305 const typename identity<std::function<
316 const typename identity<std::function<
void(
const CellIteratorBaseType &,
320 &boundary_worker = std::function<
void(
const CellIteratorBaseType &,
325 const typename identity<std::function<
void(
const CellIteratorBaseType &,
328 const CellIteratorBaseType &,
333 &face_worker = std::function<
void(
const CellIteratorBaseType &,
336 const CellIteratorBaseType &,
343 const unsigned int chunk_size = 8
350 "If you provide a cell worker function, you also need to request "
351 "that work should be done on cells by setting the 'work_on_cells' flag. "
352 "Conversely, if you don't provide a cell worker function, you "
353 "cannot set the 'work_on_cells' flag. One of these two "
354 "conditions is not satisfied."));
361 "If you provide a face worker function, you also need to request "
362 "that work should be done on interior faces by setting either the "
363 "'assemble_own_interior_faces_once' flag or the "
364 "'assemble_own_interior_faces_both' flag. "
365 "Conversely, if you don't provide a face worker function, you "
366 "cannot set either of these two flags. One of these two "
367 "conditions is not satisfied."));
372 "You can only 'specify assemble_ghost_faces_once' "
373 "OR 'assemble_ghost_faces_both', but not both of these flags."));
379 "The option 'cells_after_faces' only makes sense if you assemble on cells."));
384 "If you provide a face worker function, you also need to request "
385 "that work should be done on faces by setting the 'work_on_faces' flag. "
386 "Conversely, if you don't provide a face worker function, you "
387 "cannot set the 'work_on_faces' flag. One of these two "
388 "conditions is not satisfied."));
393 "If you provide a boundary face worker function, you also need to request "
394 "that work should be done on boundary faces by setting the 'assemble_boundary_faces' flag. "
395 "Conversely, if you don't provide a boundary face worker function, you "
396 "cannot set the 'assemble_boundary_faces' flag. One of these two "
397 "conditions is not satisfied."));
399 auto cell_action = [&](
const CellIteratorBaseType &cell,
404 copy = sample_copy_data;
407 const auto dim = cell->get_triangulation().dimension;
409 const bool ignore_subdomain =
410 (cell->get_triangulation().locally_owned_subdomain() ==
414 (cell->is_level_cell() ? cell->level_subdomain_id() :
415 cell->subdomain_id());
417 const bool own_cell =
419 (current_subdomain_id ==
420 cell->get_triangulation().locally_owned_subdomain());
422 if ((!ignore_subdomain) &&
429 cell_worker(cell, scratch, copy);
432 for (
const unsigned int face_no : cell->face_indices())
434 if (cell->at_boundary(face_no) &&
435 !cell->has_periodic_neighbor(face_no))
439 boundary_worker(cell, face_no, scratch, copy);
445 neighbor = cell->neighbor_or_periodic_neighbor(face_no);
449 if (neighbor->is_level_cell())
450 neighbor_subdomain_id = neighbor->level_subdomain_id();
452 else if (neighbor->is_active())
453 neighbor_subdomain_id = neighbor->subdomain_id();
455 const bool own_neighbor =
457 (neighbor_subdomain_id ==
458 cell->get_triangulation().locally_owned_subdomain());
461 if (!own_cell && !own_neighbor)
465 if (own_cell && own_neighbor &&
471 if (own_cell != own_neighbor &&
479 const bool periodic_neighbor =
480 cell->has_periodic_neighbor(face_no);
482 if (dim > 1 && ((!periodic_neighbor &&
483 cell->neighbor_is_coarser(face_no) &&
484 neighbor->is_active()) ||
485 (periodic_neighbor &&
486 cell->periodic_neighbor_is_coarser(face_no) &&
487 neighbor->is_active())))
496 const std::pair<unsigned int, unsigned int>
499 cell->periodic_neighbor_of_coarser_periodic_neighbor(
501 cell->neighbor_of_coarser_neighbor(face_no);
507 neighbor_face_no.first,
508 neighbor_face_no.second,
518 face_worker(neighbor,
519 neighbor_face_no.first,
520 neighbor_face_no.second,
528 else if (dim == 1 && cell->level() > neighbor->level())
531 const unsigned int neighbor_face_no =
533 cell->periodic_neighbor_face_no(face_no) :
534 cell->neighbor_face_no(face_no);
535 Assert(periodic_neighbor ||
536 neighbor->face(neighbor_face_no) ==
553 face_worker(neighbor,
568 neighbor->has_children())
573 Assert(!cell->neighbor_is_coarser(face_no),
579 if (own_cell && own_neighbor &&
592 if (own_cell && !own_neighbor &&
594 (neighbor_subdomain_id < current_subdomain_id))
597 const unsigned int neighbor_face_no =
599 cell->periodic_neighbor_face_no(face_no) :
600 cell->neighbor_face_no(face_no);
601 Assert(periodic_neighbor ||
602 neighbor->face(neighbor_face_no) ==
622 cell_worker(cell, scratch, copy);
705 template <
class CellIteratorType,
708 class CellIteratorBaseType =
713 const typename identity<std::function<
724 const typename identity<std::function<
void(
const CellIteratorBaseType &,
728 &boundary_worker = std::function<
void(
const CellIteratorBaseType &,
733 const typename identity<std::function<
void(
const CellIteratorBaseType &,
736 const CellIteratorBaseType &,
741 &face_worker = std::function<
void(
const CellIteratorBaseType &,
744 const CellIteratorBaseType &,
751 const unsigned int chunk_size = 8)
754 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
757 CellIteratorBaseType>(iterator_range.
begin(),
758 iterator_range.
end(),
829 template <
class CellIteratorType,
836 MainClass & main_class,
837 void (MainClass::*cell_worker)(
const CellIteratorType &,
840 void (MainClass::*copier)(
const CopyData &),
844 void (MainClass::*boundary_worker)(
const CellIteratorType &,
848 void (MainClass::*face_worker)(
const CellIteratorType &,
851 const CellIteratorType &,
857 const unsigned int chunk_size = 8)
866 std::function<void(
const CellIteratorType &,
869 const CellIteratorType &,
876 if (cell_worker !=
nullptr)
877 f_cell_worker = [&main_class,
878 cell_worker](
const CellIteratorType &cell_iterator,
881 (main_class.*cell_worker)(cell_iterator, scratch_data, copy_data);
884 if (boundary_worker !=
nullptr)
886 [&main_class, boundary_worker](
const CellIteratorType &cell_iterator,
887 const unsigned int face_no,
891 boundary_worker)(cell_iterator, face_no, scratch_data, copy_data);
894 if (face_worker !=
nullptr)
895 f_face_worker = [&main_class,
896 face_worker](
const CellIteratorType &cell_iterator_1,
897 const unsigned int face_index_1,
898 const unsigned int subface_index_1,
899 const CellIteratorType &cell_iterator_2,
900 const unsigned int face_index_2,
901 const unsigned int subface_index_2,
904 (main_class.*face_worker)(cell_iterator_1,
918 [&main_class, copier](
const CopyData ©_data) {
919 (main_class.*copier)(copy_data);
1009 template <
class CellIteratorType,
1013 class CellIteratorBaseType =
1017 MainClass & main_class,
1018 void (MainClass::*cell_worker)(
const CellIteratorBaseType &,
1021 void (MainClass::*copier)(
const CopyData &),
1025 void (MainClass::*boundary_worker)(
const CellIteratorBaseType &,
1029 void (MainClass::*face_worker)(
const CellIteratorBaseType &,
1032 const CellIteratorBaseType &,
1038 const unsigned int chunk_size = 8)
1041 mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
1045 CellIteratorBaseType>(iterator_range.
begin(),
1046 iterator_range.
end(),
1050 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 CellIteratorBaseType< CellIteratorType >::type type
typename CellIteratorBaseType< CellIteratorType >::type type