Reference documentation for deal.II version 9.0.0
loop.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 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_loop_h
18 #define dealii_mesh_worker_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/local_integrator.h>
26 #include <deal.II/meshworker/dof_info.h>
27 #include <deal.II/meshworker/integration_info.h>
28 
29 #include <functional>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 template <typename> class TriaActiveIterator;
34 
35 namespace internal
36 {
40  template <class DI>
41  inline bool is_active_iterator(const DI &)
42  {
43  return false;
44  }
45 
46  template <class ACCESSOR>
48  {
49  return true;
50  }
51 
52  template <class ACCESSOR>
53  inline bool is_active_iterator(const ::FilteredIterator<TriaActiveIterator<ACCESSOR> > &)
54  {
55  return true;
56  }
57 
58  template <int dim, class DOFINFO, class A>
59  void assemble(const MeshWorker::DoFInfoBox<dim, DOFINFO> &dinfo, A *assembler)
60  {
61  dinfo.assemble(*assembler);
62  }
63 }
64 
65 
66 
67 namespace MeshWorker
68 {
73  {
74  public:
79  : own_cells(true), ghost_cells(false),
81  cells_first(true)
82  {
83  }
84 
88  bool own_cells;
89 
95 
102  {
115  };
116 
129 
140 
146  };
147 
148 
149 
177  template <class INFOBOX, class DOFINFO, int dim, int spacedim, class ITERATOR>
179  ITERATOR cell,
180  DoFInfoBox<dim, DOFINFO> &dof_info,
181  INFOBOX &info,
182  const std::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker,
183  const std::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker,
184  const std::function<void (DOFINFO &, DOFINFO &,
185  typename INFOBOX::CellInfo &,
186  typename INFOBOX::CellInfo &)> &face_worker,
187  const LoopControl &loop_control)
188  {
189  const bool ignore_subdomain = (cell->get_triangulation().locally_owned_subdomain()
191 
192  types::subdomain_id csid = (cell->is_level_cell())
193  ? cell->level_subdomain_id()
194  : cell->subdomain_id();
195 
196  const bool own_cell = ignore_subdomain || (csid == cell->get_triangulation().locally_owned_subdomain());
197 
198  dof_info.reset();
199 
200  if ((!ignore_subdomain) && (csid == numbers::artificial_subdomain_id))
201  return;
202 
203  dof_info.cell.reinit(cell);
204  dof_info.cell_valid = true;
205 
206  const bool integrate_cell = (cell_worker != nullptr);
207  const bool integrate_boundary = (boundary_worker != nullptr);
208  const bool integrate_interior_face = (face_worker != nullptr);
209 
210  if (integrate_cell)
211  info.cell.reinit(dof_info.cell);
212  // Execute this, if cells
213  // have to be dealt with
214  // before faces
215  if (integrate_cell && loop_control.cells_first &&
216  ((loop_control.own_cells && own_cell) || (loop_control.ghost_cells && !own_cell)))
217  cell_worker(dof_info.cell, info.cell);
218 
219  // Call the callback function in
220  // the info box to do
221  // computations between cell and
222  // face action.
223  info.post_cell(dof_info);
224 
225  if (integrate_interior_face || integrate_boundary)
226  for (unsigned int face_no=0; face_no < GeometryInfo<ITERATOR::AccessorType::Container::dimension>::faces_per_cell; ++face_no)
227  {
228  typename ITERATOR::AccessorType::Container::face_iterator face = cell->face(face_no);
229  if (cell->at_boundary(face_no) && !cell->has_periodic_neighbor(face_no))
230  {
231  // only integrate boundary faces of own cells
232  if (integrate_boundary && own_cell)
233  {
234  dof_info.interior_face_available[face_no] = true;
235  dof_info.interior[face_no].reinit(cell, face, face_no);
236  info.boundary.reinit(dof_info.interior[face_no]);
237  boundary_worker(dof_info.interior[face_no], info.boundary);
238  }
239  }
240  else if (integrate_interior_face)
241  {
242  // Interior face
243  TriaIterator<typename ITERATOR::AccessorType> neighbor = cell->neighbor_or_periodic_neighbor(face_no);
244 
246  if (neighbor->is_level_cell())
247  neighbid = neighbor->level_subdomain_id();
248  //subdomain id is only valid for active cells
249  else if (neighbor->active())
250  neighbid = neighbor->subdomain_id();
251 
252  const bool own_neighbor = ignore_subdomain ||
253  (neighbid == cell->get_triangulation().locally_owned_subdomain());
254 
255  // skip all faces between two ghost cells
256  if (!own_cell && !own_neighbor)
257  continue;
258 
259  // skip if the user doesn't want faces between own cells
260  if (own_cell && own_neighbor && loop_control.own_faces==LoopControl::never)
261  continue;
262 
263  // skip face to ghost
264  if (own_cell != own_neighbor && loop_control.faces_to_ghost==LoopControl::never)
265  continue;
266 
267  // Deal with refinement edges from the refined side. Assuming one-irregular
268  // meshes, this situation should only occur if both cells are active.
269  const bool periodic_neighbor = cell->has_periodic_neighbor(face_no);
270 
271  if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no))
272  || (periodic_neighbor && cell->periodic_neighbor_is_coarser(face_no)))
273  {
274  Assert(!cell->has_children(), ExcInternalError());
275  Assert(!neighbor->has_children(), ExcInternalError());
276 
277  // skip if only one processor needs to assemble the face
278  // to a ghost cell and the fine cell is not ours.
279  if (!own_cell
280  && loop_control.faces_to_ghost == LoopControl::one)
281  continue;
282 
283  const std::pair<unsigned int, unsigned int> neighbor_face_no
284  = periodic_neighbor?
285  cell->periodic_neighbor_of_coarser_periodic_neighbor(face_no):
286  cell->neighbor_of_coarser_neighbor(face_no);
287  const typename ITERATOR::AccessorType::Container::face_iterator nface
288  = neighbor->face(neighbor_face_no.first);
289 
290  dof_info.interior_face_available[face_no] = true;
291  dof_info.exterior_face_available[face_no] = true;
292  dof_info.interior[face_no].reinit(cell, face, face_no);
293  info.face.reinit(dof_info.interior[face_no]);
294  dof_info.exterior[face_no].reinit(
295  neighbor, nface, neighbor_face_no.first, neighbor_face_no.second);
296  info.subface.reinit(dof_info.exterior[face_no]);
297 
298  face_worker(dof_info.interior[face_no], dof_info.exterior[face_no],
299  info.face, info.subface);
300  }
301  else
302  {
303  // If iterator is active and neighbor is refined, skip
304  // internal face.
305  if (internal::is_active_iterator(cell) && neighbor->has_children())
306  {
307  Assert(loop_control.own_faces != LoopControl::both, ExcMessage(
308  "Assembling from both sides for own_faces is not "
309  "supported with hanging nodes!"));
310  continue;
311  }
312 
313  // Now neighbor is on same level, double-check this:
314  Assert(cell->level()==neighbor->level(), ExcInternalError());
315 
316  // If we own both cells only do faces from one side (unless
317  // LoopControl says otherwise). Here, we rely on cell comparison
318  // that will look at cell->index().
319  if (own_cell && own_neighbor
320  && loop_control.own_faces == LoopControl::one
321  && (neighbor < cell))
322  continue;
323 
324  // independent of loop_control.faces_to_ghost,
325  // we only look at faces to ghost on the same level once
326  // (only where own_cell=true and own_neighbor=false)
327  if (!own_cell)
328  continue;
329 
330  // now only one processor assembles faces_to_ghost. We let the
331  // processor with the smaller (level-)subdomain id assemble the
332  // face.
333  if (own_cell && !own_neighbor
334  && loop_control.faces_to_ghost == LoopControl::one
335  && (neighbid < csid))
336  continue;
337 
338  const unsigned int neighbor_face_no = periodic_neighbor?
339  cell->periodic_neighbor_face_no(face_no):
340  cell->neighbor_face_no(face_no);
341  Assert (periodic_neighbor || neighbor->face(neighbor_face_no) == face, ExcInternalError());
342  // Regular interior face
343  dof_info.interior_face_available[face_no] = true;
344  dof_info.exterior_face_available[face_no] = true;
345  dof_info.interior[face_no].reinit(cell, face, face_no);
346  info.face.reinit(dof_info.interior[face_no]);
347  dof_info.exterior[face_no].reinit(
348  neighbor, neighbor->face(neighbor_face_no), neighbor_face_no);
349  info.neighbor.reinit(dof_info.exterior[face_no]);
350 
351  face_worker(dof_info.interior[face_no], dof_info.exterior[face_no],
352  info.face, info.neighbor);
353  }
354  }
355  } // faces
356  // Call the callback function in
357  // the info box to do
358  // computations between face and
359  // cell action.
360  info.post_faces(dof_info);
361 
362  // Execute this, if faces
363  // have to be handled first
364  if (integrate_cell && !loop_control.cells_first &&
365  ((loop_control.own_cells && own_cell) || (loop_control.ghost_cells && !own_cell)))
366  cell_worker(dof_info.cell, info.cell);
367  }
368 
369 
385  template <int dim, int spacedim, class DOFINFO, class INFOBOX, class ASSEMBLER, class ITERATOR>
386  void loop(ITERATOR begin,
387  typename identity<ITERATOR>::type end,
388  DOFINFO &dinfo,
389  INFOBOX &info,
390  const std::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker,
391  const std::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker,
392  const std::function<void (DOFINFO &, DOFINFO &,
393  typename INFOBOX::CellInfo &,
394  typename INFOBOX::CellInfo &)> &face_worker,
395  ASSEMBLER &assembler,
396  const LoopControl &lctrl = LoopControl())
397  {
398  DoFInfoBox<dim, DOFINFO> dof_info(dinfo);
399 
400  assembler.initialize_info(dof_info.cell, false);
401  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
402  {
403  assembler.initialize_info(dof_info.interior[i], true);
404  assembler.initialize_info(dof_info.exterior[i], true);
405  }
406 
407  // Loop over all cells
408  WorkStream::run(begin, end,
409  std::bind(&cell_action<INFOBOX, DOFINFO, dim, spacedim, ITERATOR>,
410  std::placeholders::_1, std::placeholders::_3, std::placeholders::_2,
411  cell_worker, boundary_worker, face_worker, lctrl),
412  std::bind(&internal::assemble<dim,DOFINFO,ASSEMBLER>, std::placeholders::_1, &assembler),
413  info, dof_info);
414  }
415 
416 
424  template <int dim, int spacedim, class ITERATOR, class ASSEMBLER>
425  void integration_loop(ITERATOR begin,
426  typename identity<ITERATOR>::type end,
427  DoFInfo<dim, spacedim> &dof_info,
429  const LocalIntegrator<dim, spacedim> &integrator,
430  ASSEMBLER &assembler,
431  const LoopControl &lctrl = LoopControl())
432  {
433  std::function<void (DoFInfo<dim, spacedim>&, IntegrationInfo<dim, spacedim>&)> cell_worker;
434  std::function<void (DoFInfo<dim, spacedim>&, IntegrationInfo<dim, spacedim>&)> boundary_worker;
435  std::function<void (DoFInfo<dim, spacedim>&, DoFInfo<dim, spacedim> &,
437  IntegrationInfo<dim, spacedim> &)> face_worker;
438  if (integrator.use_cell)
439  cell_worker = std::bind(&LocalIntegrator<dim, spacedim>::cell, &integrator, std::placeholders::_1, std::placeholders::_2);
440  if (integrator.use_boundary)
441  boundary_worker = std::bind(&LocalIntegrator<dim, spacedim>::boundary, &integrator, std::placeholders::_1, std::placeholders::_2);
442  if (integrator.use_face)
443  face_worker = std::bind(&LocalIntegrator<dim, spacedim>::face, &integrator, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4);
444 
445  loop<dim, spacedim>
446  (begin, end,
447  dof_info,
448  box,
449  cell_worker,
450  boundary_worker,
451  face_worker,
452  assembler,
453  lctrl);
454  }
455 
456 }
457 
458 DEAL_II_NAMESPACE_CLOSE
459 
460 #endif
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, 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, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:386
const types::subdomain_id invalid_subdomain_id
Definition: types.h:248
FaceOption own_faces
Definition: loop.h:139
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:260
static ::ExceptionBase & ExcMessage(std::string arg1)
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:254
#define Assert(cond, exc)
Definition: exceptions.h:1142
bool is_active_iterator(const DI &)
Definition: loop.h:41
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:461
unsigned int subdomain_id
Definition: types.h:42
const types::subdomain_id artificial_subdomain_id
Definition: types.h:264
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:244
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
void integration_loop(ITERATOR begin, typename identity< ITERATOR >::type end, DoFInfo< dim, spacedim > &dof_info, IntegrationInfoBox< dim, spacedim > &box, const LocalIntegrator< dim, spacedim > &integrator, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:425
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:248
FaceOption faces_to_ghost
Definition: loop.h:128
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()