Reference documentation for deal.II version 9.2.0
\(\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\}}\)
mesh_loop.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_mesh_worker_mesh_loop_h
18 #define dealii_mesh_worker_mesh_loop_h
19 
20 #include <deal.II/base/config.h>
21 
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
41 template <typename>
42 class TriaActiveIterator;
43 #endif
44 
45 namespace MeshWorker
46 {
47  namespace internal
48  {
53  template <class CellIteratorType>
55  {
59  using type = CellIteratorType;
60  };
61 
69  template <class CellIteratorType>
70  struct CellIteratorBaseType<IteratorOverIterators<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 <class 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 
231  template <class CellIteratorType,
232  class ScratchData,
233  class CopyData,
234  class CellIteratorBaseType =
236  void
238  const CellIteratorType & begin,
239  const typename identity<CellIteratorType>::type &end,
240 
241  const typename identity<std::function<
242  void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type
243  &cell_worker,
244  const typename identity<std::function<void(const CopyData &)>>::type
245  &copier,
246 
247  const ScratchData &sample_scratch_data,
248  const CopyData & sample_copy_data,
249 
250  const AssembleFlags flags = assemble_own_cells,
251 
252  const typename identity<std::function<void(const CellIteratorBaseType &,
253  const unsigned int,
254  ScratchData &,
255  CopyData &)>>::type
256  &boundary_worker = std::function<void(const CellIteratorBaseType &,
257  const unsigned int,
258  ScratchData &,
259  CopyData &)>(),
260 
261  const typename identity<std::function<void(const CellIteratorBaseType &,
262  const unsigned int,
263  const unsigned int,
264  const CellIteratorBaseType &,
265  const unsigned int,
266  const unsigned int,
267  ScratchData &,
268  CopyData &)>>::type
269  &face_worker = std::function<void(const CellIteratorBaseType &,
270  const unsigned int,
271  const unsigned int,
272  const CellIteratorBaseType &,
273  const unsigned int,
274  const unsigned int,
275  ScratchData &,
276  CopyData &)>(),
277 
278  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
279  const unsigned int chunk_size = 8)
280  {
281  Assert(
282  (!cell_worker) == !(flags & work_on_cells),
283  ExcMessage(
284  "If you specify a cell_worker, you need to set assemble_own_cells or assemble_ghost_cells."));
285 
286  Assert(
287  (flags &
290  ExcMessage(
291  "You can only specify assemble_own_interior_faces_once OR assemble_own_interior_faces_both."));
292 
293  Assert(
296  ExcMessage(
297  "You can only specify assemble_ghost_faces_once OR assemble_ghost_faces_both."));
298 
299  Assert(
300  !(flags & cells_after_faces) ||
302  ExcMessage(
303  "The option cells_after_faces only makes sense if you assemble on cells."));
304 
305  Assert((!face_worker) == !(flags & work_on_faces),
306  ExcMessage(
307  "If you specify a face_worker, assemble_face_* needs to be set."));
308 
309  Assert(
310  (!boundary_worker) == !(flags & assemble_boundary_faces),
311  ExcMessage(
312  "If you specify a boundary_worker, assemble_boundary_faces needs to be set."));
313 
314  auto cell_action = [&](const CellIteratorBaseType &cell,
315  ScratchData & scratch,
316  CopyData & copy) {
317  // First reset the CopyData class to the empty copy_data given by the
318  // user.
319  copy = sample_copy_data;
320 
321  const bool ignore_subdomain =
322  (cell->get_triangulation().locally_owned_subdomain() ==
324 
325  types::subdomain_id current_subdomain_id =
326  (cell->is_level_cell() ? cell->level_subdomain_id() :
327  cell->subdomain_id());
328 
329  const bool own_cell =
330  ignore_subdomain ||
331  (current_subdomain_id ==
332  cell->get_triangulation().locally_owned_subdomain());
333 
334  if ((!ignore_subdomain) &&
335  (current_subdomain_id == numbers::artificial_subdomain_id))
336  return;
337 
338  if (!(flags & (cells_after_faces)) &&
339  (((flags & (assemble_own_cells)) && own_cell) ||
340  ((flags & assemble_ghost_cells) && !own_cell)))
341  cell_worker(cell, scratch, copy);
342 
343  if (flags & (work_on_faces | work_on_boundary))
344  for (const unsigned int face_no :
345  GeometryInfo<CellIteratorBaseType::AccessorType::Container::
346  dimension>::face_indices())
347  {
348  if (cell->at_boundary(face_no) &&
349  !cell->has_periodic_neighbor(face_no))
350  {
351  // only integrate boundary faces of own cells
352  if ((flags & assemble_boundary_faces) && own_cell)
353  boundary_worker(cell, face_no, scratch, copy);
354  }
355  else
356  {
357  // interior face, potentially assemble
359  neighbor = cell->neighbor_or_periodic_neighbor(face_no);
360 
361  types::subdomain_id neighbor_subdomain_id =
363  if (neighbor->is_level_cell())
364  neighbor_subdomain_id = neighbor->level_subdomain_id();
365  // subdomain id is only valid for active cells
366  else if (neighbor->is_active())
367  neighbor_subdomain_id = neighbor->subdomain_id();
368 
369  const bool own_neighbor =
370  ignore_subdomain ||
371  (neighbor_subdomain_id ==
372  cell->get_triangulation().locally_owned_subdomain());
373 
374  // skip all faces between two ghost cells
375  if (!own_cell && !own_neighbor)
376  continue;
377 
378  // skip if the user doesn't want faces between own cells
379  if (own_cell && own_neighbor &&
382  continue;
383 
384  // skip face to ghost
385  if (own_cell != own_neighbor &&
386  !(flags &
388  continue;
389 
390  // Deal with refinement edges from the refined side. Assuming
391  // one-irregular meshes, this situation should only occur if
392  // both cells are active.
393  const bool periodic_neighbor =
394  cell->has_periodic_neighbor(face_no);
395 
396  if ((!periodic_neighbor &&
397  cell->neighbor_is_coarser(face_no)) ||
398  (periodic_neighbor &&
399  cell->periodic_neighbor_is_coarser(face_no)))
400  {
401  Assert(cell->is_active(), ExcInternalError());
402  Assert(neighbor->is_active(), ExcInternalError());
403 
404  // skip if only one processor needs to assemble the face
405  // to a ghost cell and the fine cell is not ours.
406  if (!own_cell && (flags & assemble_ghost_faces_once))
407  continue;
408 
409  const std::pair<unsigned int, unsigned int>
410  neighbor_face_no =
411  periodic_neighbor ?
412  cell->periodic_neighbor_of_coarser_periodic_neighbor(
413  face_no) :
414  cell->neighbor_of_coarser_neighbor(face_no);
415 
416  face_worker(cell,
417  face_no,
419  neighbor,
420  neighbor_face_no.first,
421  neighbor_face_no.second,
422  scratch,
423  copy);
424 
426  {
427  // If own faces are to be assembled from both sides,
428  // call the faceworker again with swapped arguments.
429  // This is because we won't be looking at an adaptively
430  // refined edge coming from the other side.
431  face_worker(neighbor,
432  neighbor_face_no.first,
433  neighbor_face_no.second,
434  cell,
435  face_no,
437  scratch,
438  copy);
439  }
440  }
441  else
442  {
443  // If iterator is active and neighbor is refined, skip
444  // internal face.
445  if (::internal::is_active_iterator(cell) &&
446  neighbor->has_children())
447  continue;
448 
449  // Now neighbor is on same level, double-check this:
450  Assert(cell->level() == neighbor->level(),
451  ExcInternalError());
452 
453  // If we own both cells only do faces from one side (unless
454  // AssembleFlags says otherwise). Here, we rely on cell
455  // comparison that will look at cell->index().
456  if (own_cell && own_neighbor &&
458  (neighbor < cell))
459  continue;
460 
461  // We only look at faces to ghost on the same level once
462  // (only where own_cell=true and own_neighbor=false)
463  if (!own_cell)
464  continue;
465 
466  // now only one processor assembles faces_to_ghost. We let
467  // the processor with the smaller (level-)subdomain id
468  // assemble the face.
469  if (own_cell && !own_neighbor &&
470  (flags & assemble_ghost_faces_once) &&
471  (neighbor_subdomain_id < current_subdomain_id))
472  continue;
473 
474  const unsigned int neighbor_face_no =
475  periodic_neighbor ?
476  cell->periodic_neighbor_face_no(face_no) :
477  cell->neighbor_face_no(face_no);
478  Assert(periodic_neighbor ||
479  neighbor->face(neighbor_face_no) ==
480  cell->face(face_no),
481  ExcInternalError());
482 
483  face_worker(cell,
484  face_no,
486  neighbor,
487  neighbor_face_no,
489  scratch,
490  copy);
491  }
492  }
493  } // faces
494 
495  // Execute the cell_worker if faces are handled before cells
496  if ((flags & cells_after_faces) &&
497  (((flags & assemble_own_cells) && own_cell) ||
498  ((flags & assemble_ghost_cells) && !own_cell)))
499  cell_worker(cell, scratch, copy);
500  };
501 
502  // Submit to workstream
504  end,
505  cell_action,
506  copier,
507  sample_scratch_data,
508  sample_copy_data,
509  queue_length,
510  chunk_size);
511  }
512 
582  template <class CellIteratorType,
583  class ScratchData,
584  class CopyData,
585  class CellIteratorBaseType =
587  void
589  IteratorRange<CellIteratorType> iterator_range,
590  const typename identity<std::function<
591  void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type
592  &cell_worker,
593  const typename identity<std::function<void(const CopyData &)>>::type
594  &copier,
595 
596  const ScratchData &sample_scratch_data,
597  const CopyData & sample_copy_data,
598 
599  const AssembleFlags flags = assemble_own_cells,
600 
601  const typename identity<std::function<void(const CellIteratorBaseType &,
602  const unsigned int,
603  ScratchData &,
604  CopyData &)>>::type
605  &boundary_worker = std::function<void(const CellIteratorBaseType &,
606  const unsigned int,
607  ScratchData &,
608  CopyData &)>(),
609 
610  const typename identity<std::function<void(const CellIteratorBaseType &,
611  const unsigned int,
612  const unsigned int,
613  const CellIteratorBaseType &,
614  const unsigned int,
615  const unsigned int,
616  ScratchData &,
617  CopyData &)>>::type
618  &face_worker = std::function<void(const CellIteratorBaseType &,
619  const unsigned int,
620  const unsigned int,
621  const CellIteratorBaseType &,
622  const unsigned int,
623  const unsigned int,
624  ScratchData &,
625  CopyData &)>(),
626 
627  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
628  const unsigned int chunk_size = 8)
629  {
630  // Call the function above
631  mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
632  ScratchData,
633  CopyData,
634  CellIteratorBaseType>(iterator_range.begin(),
635  iterator_range.end(),
636  cell_worker,
637  copier,
638  sample_scratch_data,
639  sample_copy_data,
640  flags,
641  boundary_worker,
642  face_worker,
643  queue_length,
644  chunk_size);
645  }
646 
707  template <class CellIteratorType,
708  class ScratchData,
709  class CopyData,
710  class MainClass>
711  void
712  mesh_loop(const CellIteratorType & begin,
713  const typename identity<CellIteratorType>::type &end,
714  MainClass & main_class,
715  void (MainClass::*cell_worker)(const CellIteratorType &,
716  ScratchData &,
717  CopyData &),
718  void (MainClass::*copier)(const CopyData &),
719  const ScratchData & sample_scratch_data,
720  const CopyData & sample_copy_data,
721  const AssembleFlags flags = assemble_own_cells,
722  void (MainClass::*boundary_worker)(const CellIteratorType &,
723  const unsigned int,
724  ScratchData &,
725  CopyData &) = nullptr,
726  void (MainClass::*face_worker)(const CellIteratorType &,
727  const unsigned int,
728  const unsigned int,
729  const CellIteratorType &,
730  const unsigned int,
731  const unsigned int,
732  ScratchData &,
733  CopyData &) = nullptr,
734  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
735  const unsigned int chunk_size = 8)
736  {
737  std::function<void(const CellIteratorType &, ScratchData &, CopyData &)>
738  f_cell_worker;
739 
740  std::function<void(
741  const CellIteratorType &, const unsigned int, ScratchData &, CopyData &)>
742  f_boundary_worker;
743 
744  std::function<void(const CellIteratorType &,
745  const unsigned int,
746  const unsigned int,
747  const CellIteratorType &,
748  const unsigned int,
749  const unsigned int,
750  ScratchData &,
751  CopyData &)>
752  f_face_worker;
753 
754  if (cell_worker != nullptr)
755  f_cell_worker = [&main_class,
756  cell_worker](const CellIteratorType &cell_iterator,
757  ScratchData & scratch_data,
758  CopyData & copy_data) {
759  (main_class.*cell_worker)(cell_iterator, scratch_data, copy_data);
760  };
761 
762  if (boundary_worker != nullptr)
763  f_boundary_worker =
764  [&main_class, boundary_worker](const CellIteratorType &cell_iterator,
765  const unsigned int face_no,
766  ScratchData & scratch_data,
767  CopyData & copy_data) {
768  (main_class.*
769  boundary_worker)(cell_iterator, face_no, scratch_data, copy_data);
770  };
771 
772  if (face_worker != nullptr)
773  f_face_worker = [&main_class,
774  face_worker](const CellIteratorType &cell_iterator_1,
775  const unsigned int face_index_1,
776  const unsigned int subface_index_1,
777  const CellIteratorType &cell_iterator_2,
778  const unsigned int face_index_2,
779  const unsigned int subface_index_2,
780  ScratchData & scratch_data,
781  CopyData & copy_data) {
782  (main_class.*face_worker)(cell_iterator_1,
783  face_index_1,
784  subface_index_1,
785  cell_iterator_2,
786  face_index_2,
787  subface_index_2,
788  scratch_data,
789  copy_data);
790  };
791 
793  end,
794  f_cell_worker,
795  [&main_class, copier](const CopyData &copy_data) {
796  (main_class.*copier)(copy_data);
797  },
798  sample_scratch_data,
799  sample_copy_data,
800  flags,
801  f_boundary_worker,
802  f_face_worker,
803  queue_length,
804  chunk_size);
805  }
806 
886  template <class CellIteratorType,
887  class ScratchData,
888  class CopyData,
889  class MainClass,
890  class CellIteratorBaseType =
892  void
894  MainClass & main_class,
895  void (MainClass::*cell_worker)(const CellIteratorBaseType &,
896  ScratchData &,
897  CopyData &),
898  void (MainClass::*copier)(const CopyData &),
899  const ScratchData & sample_scratch_data,
900  const CopyData & sample_copy_data,
901  const AssembleFlags flags = assemble_own_cells,
902  void (MainClass::*boundary_worker)(const CellIteratorBaseType &,
903  const unsigned int,
904  ScratchData &,
905  CopyData &) = nullptr,
906  void (MainClass::*face_worker)(const CellIteratorBaseType &,
907  const unsigned int,
908  const unsigned int,
909  const CellIteratorBaseType &,
910  const unsigned int,
911  const unsigned int,
912  ScratchData &,
913  CopyData &) = nullptr,
914  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
915  const unsigned int chunk_size = 8)
916  {
917  // Call the function above
918  mesh_loop<typename IteratorRange<CellIteratorType>::IteratorOverIterators,
919  ScratchData,
920  CopyData,
921  MainClass,
922  CellIteratorBaseType>(iterator_range.begin(),
923  iterator_range.end(),
924  main_class,
925  cell_worker,
926  copier,
927  sample_scratch_data,
928  sample_copy_data,
929  flags,
930  boundary_worker,
931  face_worker,
932  queue_length,
933  chunk_size);
934  }
935 } // namespace MeshWorker
936 
938 
939 #endif
identity::type
T type
Definition: template_constraints.h:270
MeshWorker::assemble_ghost_faces_once
@ assemble_ghost_faces_once
Definition: assemble_flags.h:71
assemble_flags.h
MeshWorker::CopyData
Definition: copy_data.h:52
MultithreadInfo::n_threads
static unsigned int n_threads()
Definition: multithread_info.cc:169
MeshWorker::cell_action
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)
Definition: loop.h:193
MeshWorker::assemble_own_cells
@ assemble_own_cells
Definition: assemble_flags.h:50
tria.h
MeshWorker
Definition: assemble_flags.h:30
MeshWorker::mesh_loop
void mesh_loop(const CellIteratorType &begin, const typename identity< CellIteratorType >::type &end, const typename identity< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type &cell_worker, const typename identity< std::function< void(const CopyData &)>>::type &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>>::type &boundary_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>(), const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>>::type &face_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition: mesh_loop.h:237
GeometryInfo
Definition: geometry_info.h:1224
MeshWorker::work_on_cells
@ work_on_cells
Definition: assemble_flags.h:92
MeshWorker::internal::CellIteratorBaseType::type
CellIteratorType type
Definition: mesh_loop.h:59
IteratorOverIterators
Definition: iterator_range.h:200
MeshWorker::work_on_boundary
@ work_on_boundary
Definition: assemble_flags.h:105
identity
Definition: template_constraints.h:268
MeshWorker::assemble_ghost_faces_both
@ assemble_ghost_faces_both
Definition: assemble_flags.h:77
integration_info.h
internal::is_active_iterator
bool is_active_iterator(const DI &)
Definition: loop.h:49
MeshWorker::cells_after_faces
@ cells_after_faces
Definition: assemble_flags.h:87
MeshWorker::internal::CellIteratorBaseType< IteratorOverIterators< CellIteratorType > >::type
typename CellIteratorBaseType< CellIteratorType >::type type
Definition: mesh_loop.h:78
WorkStream::run
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)
Definition: work_stream.h:1185
work_stream.h
MeshWorker::ScratchData
Definition: scratch_data.h:232
FilteredIterator
Definition: filtered_iterator.h:529
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
IteratorRange
Definition: iterator_range.h:129
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
MeshWorker::assemble_boundary_faces
@ assemble_boundary_faces
Definition: assemble_flags.h:81
dof_info.h
TriaActiveIterator
Definition: tria_iterator.h:759
MeshWorker::work_on_faces
@ work_on_faces
Definition: assemble_flags.h:97
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
MeshWorker::assemble_own_interior_faces_both
@ assemble_own_interior_faces_both
Definition: assemble_flags.h:65
loop.h
CUDAWrappers::chunk_size
constexpr int chunk_size
Definition: cuda_size.h:35
MeshWorker::internal::CellIteratorBaseType< FilteredIterator< CellIteratorType > >::type
typename CellIteratorBaseType< CellIteratorType >::type type
Definition: mesh_loop.h:97
MeshWorker::internal::CellIteratorBaseType
Definition: mesh_loop.h:54
IteratorRange::begin
IteratorOverIterators begin()
Definition: iterator_range.h:390
MeshWorker::assemble_own_interior_faces_once
@ assemble_own_interior_faces_once
Definition: assemble_flags.h:59
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
local_integrator.h
template_constraints.h
config.h
numbers::invalid_subdomain_id
const types::subdomain_id invalid_subdomain_id
Definition: types.h:285
numbers::artificial_subdomain_id
const types::subdomain_id artificial_subdomain_id
Definition: types.h:302
internal
Definition: aligned_vector.h:369
MeshWorker::assemble_ghost_cells
@ assemble_ghost_cells
Definition: assemble_flags.h:54
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TriaIterator
Definition: tria_iterator.h:578
IteratorRange::end
IteratorOverIterators end() const
Definition: iterator_range.h:406
filtered_iterator.h
MeshWorker::AssembleFlags
AssembleFlags
Definition: assemble_flags.h:41
internal::VectorOperations::copy
void copy(const T *begin, const T *end, U *dest)
Definition: vector_operations_internal.h:67