Reference documentation for deal.II version 9.0.0
mesh_loop.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/work_stream.h>
22 #include <deal.II/base/template_constraints.h>
23 #include <deal.II/grid/tria.h>
24 #include <deal.II/grid/filtered_iterator.h>
25 #include <deal.II/meshworker/loop.h>
26 #include <deal.II/meshworker/local_integrator.h>
27 #include <deal.II/meshworker/dof_info.h>
28 #include <deal.II/meshworker/integration_info.h>
29 #include <deal.II/meshworker/assemble_flags.h>
30 
31 #include <functional>
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 template <typename> class TriaActiveIterator;
36 
37 namespace MeshWorker
38 {
114  template <class CellIteratorType,
115  class ScratchData, class CopyData>
116  void mesh_loop(const CellIteratorType &begin,
117  const typename identity<CellIteratorType>::type &end,
118 
119  const typename identity<std::function<void (const CellIteratorType &, ScratchData &, CopyData &)>>::type &cell_worker,
120  const typename identity<std::function<void (const CopyData &)>>::type &copier,
121 
122  const ScratchData &sample_scratch_data,
123  const CopyData &sample_copy_data,
124 
125  const AssembleFlags flags = assemble_own_cells,
126 
127  const typename identity<std::function<void (const CellIteratorType &, const unsigned int &, ScratchData &, CopyData &)>>::type &boundary_worker=
128  std::function<void (const CellIteratorType &, const unsigned int &, ScratchData &, CopyData &)>(),
129 
130  const typename identity<std::function<void (const CellIteratorType &, const unsigned int &, const unsigned int &,
131  const CellIteratorType &, const unsigned int &, const unsigned int &,
132  ScratchData &, CopyData &)>>::type &face_worker=
133  std::function<void (const CellIteratorType &, const unsigned int &, const unsigned int &,
134  const CellIteratorType &, const unsigned int &, const unsigned int &,
135  ScratchData &, CopyData &)>(),
136 
137  const unsigned int queue_length = 2*MultithreadInfo::n_threads(),
138  const unsigned int chunk_size = 8)
139  {
140  Assert((!cell_worker) == !(flags & work_on_cells),
141  ExcMessage("If you specify a cell_worker, you need to set assemble_own_cells or assemble_ghost_cells."));
142 
145  ExcMessage("You can only specify assemble_own_interior_faces_once OR assemble_own_interior_faces_both."));
146 
149  ExcMessage("You can only specify assemble_ghost_faces_once OR assemble_ghost_faces_both."));
150 
151  Assert(!(flags & cells_after_faces) ||
153  ExcMessage("The option cells_after_faces only makes sense if you assemble on cells."));
154 
155  Assert((!face_worker) == !(flags & work_on_faces),
156  ExcMessage("If you specify a face_worker, assemble_face_* needs to be set."));
157 
158  Assert((!boundary_worker) == !(flags & assemble_boundary_faces),
159  ExcMessage("If you specify a boundary_worker, assemble_boundary_faces needs to be set."));
160 
161  auto cell_action = [&] (const CellIteratorType &cell, ScratchData &scratch, CopyData &copy)
162  {
163  // First reset the CopyData class to the empty copy_data given by the user.
164  copy = sample_copy_data;
165 
166  const bool ignore_subdomain = (cell->get_triangulation().locally_owned_subdomain()
168 
169  types::subdomain_id current_subdomain_id = (cell->is_level_cell()
170  ? cell->level_subdomain_id()
171  : cell->subdomain_id());
172 
173  const bool own_cell = ignore_subdomain || (current_subdomain_id == cell->get_triangulation().locally_owned_subdomain());
174 
175  if ((!ignore_subdomain) && (current_subdomain_id == numbers::artificial_subdomain_id))
176  return;
177 
178  if ( !(flags & (cells_after_faces)) &&
179  ( ((flags & (assemble_own_cells)) && own_cell)
180  || ( (flags & assemble_ghost_cells) && !own_cell) ) )
181  cell_worker(cell, scratch, copy);
182 
183  if (flags & (work_on_faces | work_on_boundary))
184  for (unsigned int face_no=0; face_no < GeometryInfo<CellIteratorType::AccessorType::Container::dimension>::faces_per_cell; ++face_no)
185  {
186  typename CellIteratorType::AccessorType::Container::face_iterator face = cell->face(face_no);
187  if (cell->at_boundary(face_no) && !cell->has_periodic_neighbor(face_no))
188  {
189  // only integrate boundary faces of own cells
190  if ( (flags & assemble_boundary_faces) && own_cell)
191  boundary_worker(cell, face_no, scratch, copy);
192  }
193  else
194  {
195  // interior face, potentially assemble
196  TriaIterator<typename CellIteratorType::AccessorType> neighbor = cell->neighbor_or_periodic_neighbor(face_no);
197 
199  if (neighbor->is_level_cell())
200  neighbor_subdomain_id = neighbor->level_subdomain_id();
201  //subdomain id is only valid for active cells
202  else if (neighbor->active())
203  neighbor_subdomain_id = neighbor->subdomain_id();
204 
205  const bool own_neighbor = ignore_subdomain ||
206  (neighbor_subdomain_id == cell->get_triangulation().locally_owned_subdomain());
207 
208  // skip all faces between two ghost cells
209  if (!own_cell && !own_neighbor)
210  continue;
211 
212  // skip if the user doesn't want faces between own cells
213  if (own_cell && own_neighbor && !(flags & (assemble_own_interior_faces_both | assemble_own_interior_faces_once)))
214  continue;
215 
216  // skip face to ghost
217  if (own_cell != own_neighbor && !(flags & (assemble_ghost_faces_both | assemble_ghost_faces_once)))
218  continue;
219 
220  // Deal with refinement edges from the refined side. Assuming one-irregular
221  // meshes, this situation should only occur if both cells are active.
222  const bool periodic_neighbor = cell->has_periodic_neighbor(face_no);
223 
224  if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no))
225  || (periodic_neighbor && cell->periodic_neighbor_is_coarser(face_no)))
226  {
227  Assert(!cell->has_children(), ExcInternalError());
228  Assert(!neighbor->has_children(), ExcInternalError());
229 
230  // skip if only one processor needs to assemble the face
231  // to a ghost cell and the fine cell is not ours.
232  if (!own_cell && (flags & assemble_ghost_faces_once))
233  continue;
234 
235  const std::pair<unsigned int, unsigned int> neighbor_face_no
236  = periodic_neighbor?
237  cell->periodic_neighbor_of_coarser_periodic_neighbor(face_no):
238  cell->neighbor_of_coarser_neighbor(face_no);
239 
240  face_worker(cell, face_no, numbers::invalid_unsigned_int,
241  neighbor, neighbor_face_no.first, neighbor_face_no.second,
242  scratch, copy);
243 
245  {
246  // If own faces are to be assembled from both sides, call the
247  // faceworker again with swapped arguments. This is because
248  // we won't be looking at an adaptively refined edge
249  // coming from the other side.
250  face_worker(neighbor, neighbor_face_no.first, neighbor_face_no.second,
251  cell, face_no, numbers::invalid_unsigned_int,
252  scratch, copy);
253  }
254  }
255  else
256  {
257  // If iterator is active and neighbor is refined, skip
258  // internal face.
259  if (internal::is_active_iterator(cell) && neighbor->has_children())
260  continue;
261 
262  // Now neighbor is on same level, double-check this:
263  Assert(cell->level()==neighbor->level(), ExcInternalError());
264 
265  // If we own both cells only do faces from one side (unless
266  // AssembleFlags says otherwise). Here, we rely on cell comparison
267  // that will look at cell->index().
268  if (own_cell && own_neighbor
270  && (neighbor < cell))
271  continue;
272 
273  // We only look at faces to ghost on the same level once
274  // (only where own_cell=true and own_neighbor=false)
275  if (!own_cell)
276  continue;
277 
278  // now only one processor assembles faces_to_ghost. We let the
279  // processor with the smaller (level-)subdomain id assemble the
280  // face.
281  if (own_cell && !own_neighbor
282  && (flags & assemble_ghost_faces_once)
283  && (neighbor_subdomain_id < current_subdomain_id))
284  continue;
285 
286  const unsigned int neighbor_face_no = periodic_neighbor?
287  cell->periodic_neighbor_face_no(face_no):
288  cell->neighbor_face_no(face_no);
289  Assert (periodic_neighbor || neighbor->face(neighbor_face_no) == face, ExcInternalError());
290 
291  face_worker(cell, face_no, numbers::invalid_unsigned_int,
292  neighbor, neighbor_face_no, numbers::invalid_unsigned_int,
293  scratch, copy);
294  }
295  }
296  } // faces
297 
298  // Execute the cell_worker if faces are handled before cells
299  if ((flags & cells_after_faces) &&
300  ( ((flags & assemble_own_cells) && own_cell) || ((flags & assemble_ghost_cells) && !own_cell)))
301  cell_worker(cell, scratch, copy);
302  };
303 
304  // Submit to workstream
305  WorkStream::run(begin, end,
306  cell_action, copier,
307  sample_scratch_data, sample_copy_data,
308  queue_length, chunk_size);
309  }
310 }
311 
312 DEAL_II_NAMESPACE_CLOSE
313 
314 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:173
void mesh_loop(const CellIteratorType &begin, const typename identity< CellIteratorType >::type &end, const typename identity< std::function< void(const CellIteratorType &, 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 CellIteratorType &, const unsigned int &, ScratchData &, CopyData &)>>::type &boundary_worker=std::function< void(const CellIteratorType &, const unsigned int &, ScratchData &, CopyData &)>(), const typename identity< std::function< void(const CellIteratorType &, const unsigned int &, const unsigned int &, const CellIteratorType &, const unsigned int &, const unsigned int &, ScratchData &, CopyData &)>>::type &face_worker=std::function< void(const CellIteratorType &, const unsigned int &, const unsigned int &, const CellIteratorType &, 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:116
const types::subdomain_id invalid_subdomain_id
Definition: types.h:248
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
bool is_active_iterator(const DI &)
Definition: loop.h:41
unsigned int subdomain_id
Definition: types.h:42
const types::subdomain_id artificial_subdomain_id
Definition: types.h:264
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:1106
static unsigned int n_threads()
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:178
static ::ExceptionBase & ExcInternalError()