Reference documentation for deal.II version 8.4.2
loop.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2015 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/std_cxx11/function.h>
22 #include <deal.II/base/work_stream.h>
23 #include <deal.II/base/template_constraints.h>
24 #include <deal.II/grid/tria.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 
30 #define DEAL_II_MESHWORKER_PARALLEL 1
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 template <typename> class TriaActiveIterator;
35 template <typename> class FilteredIterator;
36 
37 namespace internal
38 {
42  template <class DI>
43  inline bool is_active_iterator(const DI &)
44  {
45  return false;
46  }
47 
48  template <class ACCESSOR>
50  {
51  return true;
52  }
53 
54  template <class ACCESSOR>
56  {
57  return true;
58  }
59 
60  template<int dim, class DOFINFO, class A>
61  void assemble(const MeshWorker::DoFInfoBox<dim, DOFINFO> &dinfo, A *assembler)
62  {
63  dinfo.assemble(*assembler);
64  }
65 }
66 
67 
68 
69 namespace MeshWorker
70 {
75  {
76  public:
77 
82  : own_cells(true), ghost_cells(false),
83  faces_to_ghost(LoopControl::one), own_faces(LoopControl::one),
84  cells_first(true)
85  {
86  }
87 
91  bool own_cells;
97 
98  enum FaceOption
99  {
100  never,
101  one,
102  both
103  };
104 
113  FaceOption faces_to_ghost;
114 
121  FaceOption own_faces;
122 
123 
129  };
130 
131 
132 
160  template<class INFOBOX, class DOFINFO, int dim, int spacedim, class ITERATOR>
162  ITERATOR cell,
163  DoFInfoBox<dim, DOFINFO> &dof_info,
164  INFOBOX &info,
165  const std_cxx11::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker,
166  const std_cxx11::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker,
167  const std_cxx11::function<void (DOFINFO &, DOFINFO &,
168  typename INFOBOX::CellInfo &,
169  typename INFOBOX::CellInfo &)> &face_worker,
170  const LoopControl &loop_control)
171  {
172  const bool ignore_subdomain = (cell->get_triangulation().locally_owned_subdomain()
174 
175  types::subdomain_id csid = (cell->is_level_cell())
176  ? cell->level_subdomain_id()
177  : cell->subdomain_id();
178 
179  const bool own_cell = ignore_subdomain || (csid == cell->get_triangulation().locally_owned_subdomain());
180 
181  dof_info.reset();
182 
183  if ((!ignore_subdomain) && (csid == numbers::artificial_subdomain_id))
184  return;
185 
186  dof_info.cell.reinit(cell);
187  dof_info.cell_valid = true;
188 
189  const bool integrate_cell = (cell_worker != 0);
190  const bool integrate_boundary = (boundary_worker != 0);
191  const bool integrate_interior_face = (face_worker != 0);
192 
193  if (integrate_cell)
194  info.cell.reinit(dof_info.cell);
195  // Execute this, if cells
196  // have to be dealt with
197  // before faces
198  if (integrate_cell && loop_control.cells_first &&
199  ((loop_control.own_cells && own_cell) || (loop_control.ghost_cells && !own_cell)))
200  cell_worker(dof_info.cell, info.cell);
201 
202  // Call the callback function in
203  // the info box to do
204  // computations between cell and
205  // face action.
206  info.post_cell(dof_info);
207 
208  if (integrate_interior_face || integrate_boundary)
209  for (unsigned int face_no=0; face_no < GeometryInfo<ITERATOR::AccessorType::Container::dimension>::faces_per_cell; ++face_no)
210  {
211  typename ITERATOR::AccessorType::Container::face_iterator face = cell->face(face_no);
212  if (cell->at_boundary(face_no))
213  {
214  // only integrate boundary faces of own cells
215  if (integrate_boundary && own_cell)
216  {
217  dof_info.interior_face_available[face_no] = true;
218  dof_info.interior[face_no].reinit(cell, face, face_no);
219  info.boundary.reinit(dof_info.interior[face_no]);
220  boundary_worker(dof_info.interior[face_no], info.boundary);
221  }
222  }
223  else if (integrate_interior_face)
224  {
225  // Interior face
226  TriaIterator<typename ITERATOR::AccessorType> neighbor = cell->neighbor(face_no);
227 
229  if (neighbor->is_level_cell())
230  neighbid = neighbor->level_subdomain_id();
231  //subdomain id is only valid for active cells
232  else if (neighbor->active())
233  neighbid = neighbor->subdomain_id();
234 
235  const bool own_neighbor = ignore_subdomain ||
236  (neighbid == cell->get_triangulation().locally_owned_subdomain());
237 
238  // skip all faces between two ghost cells
239  if (!own_cell && !own_neighbor)
240  continue;
241 
242  // skip if the user doesn't want faces between own cells
243  if (own_cell && own_neighbor && loop_control.own_faces==LoopControl::never)
244  continue;
245 
246  // skip face to ghost
247  if (own_cell != own_neighbor && loop_control.faces_to_ghost==LoopControl::never)
248  continue;
249 
250  // Deal with
251  // refinement edges
252  // from the refined
253  // side. Assuming
254  // one-irregular
255  // meshes, this
256  // situation should
257  // only occur if
258  // both cells are
259  // active.
260  if (cell->neighbor_is_coarser(face_no))
261  {
262  Assert(!cell->has_children(), ExcInternalError());
263  Assert(!neighbor->has_children(), ExcInternalError());
264 
265  // skip if only one processor needs to assemble the face
266  // to a ghost cell and the fine cell is not ours.
267  if (!own_cell
268  && loop_control.faces_to_ghost == LoopControl::one)
269  continue;
270 
271  const std::pair<unsigned int, unsigned int> neighbor_face_no
272  = cell->neighbor_of_coarser_neighbor(face_no);
273  const typename ITERATOR::AccessorType::Container::face_iterator nface
274  = neighbor->face(neighbor_face_no.first);
275 
276  dof_info.interior_face_available[face_no] = true;
277  dof_info.exterior_face_available[face_no] = true;
278  dof_info.interior[face_no].reinit(cell, face, face_no);
279  info.face.reinit(dof_info.interior[face_no]);
280  dof_info.exterior[face_no].reinit(
281  neighbor, nface, neighbor_face_no.first, neighbor_face_no.second);
282  info.subface.reinit(dof_info.exterior[face_no]);
283 
284  face_worker(dof_info.interior[face_no], dof_info.exterior[face_no],
285  info.face, info.subface);
286  }
287  else
288  {
289  // If iterator is active and neighbor is refined, skip
290  // internal face.
291  if (internal::is_active_iterator(cell) && neighbor->has_children())
292  {
293  Assert(loop_control.own_faces != LoopControl::both, ExcMessage(
294  "Assembling from both sides for own_faces is not "
295  "supported with hanging nodes!"));
296  continue;
297  }
298 
299  // Now neighbor is on same level, double-check this:
300  Assert(cell->level()==neighbor->level(), ExcInternalError());
301 
302  // If we own both cells only do faces from one side (unless
303  // LoopControl says otherwise). Here, we rely on cell comparison
304  // that will look at cell->index().
305  if (own_cell && own_neighbor
306  && loop_control.own_faces == LoopControl::one
307  && (neighbor < cell))
308  continue;
309 
310  // independent of loop_control.faces_to_ghost,
311  // we only look at faces to ghost on the same level once
312  // (only where own_cell=true and own_neighbor=false)
313  if (!own_cell)
314  continue;
315 
316  // now only one processor assembles faces_to_ghost. We let the
317  // processor with the smaller (level-)subdomain id assemble the
318  // face.
319  if (own_cell && !own_neighbor
320  && loop_control.faces_to_ghost == LoopControl::one
321  && (neighbid < csid))
322  continue;
323 
324  const unsigned int neighbor_face_no = cell->neighbor_face_no(face_no);
325  Assert (neighbor->face(neighbor_face_no) == face, ExcInternalError());
326  // Regular interior face
327  dof_info.interior_face_available[face_no] = true;
328  dof_info.exterior_face_available[face_no] = true;
329  dof_info.interior[face_no].reinit(cell, face, face_no);
330  info.face.reinit(dof_info.interior[face_no]);
331  dof_info.exterior[face_no].reinit(
332  neighbor, neighbor->face(neighbor_face_no), neighbor_face_no);
333  info.neighbor.reinit(dof_info.exterior[face_no]);
334 
335  face_worker(dof_info.interior[face_no], dof_info.exterior[face_no],
336  info.face, info.neighbor);
337  }
338  }
339  } // faces
340  // Call the callback function in
341  // the info box to do
342  // computations between face and
343  // cell action.
344  info.post_faces(dof_info);
345 
346  // Execute this, if faces
347  // have to be handled first
348  if (integrate_cell && !loop_control.cells_first &&
349  ((loop_control.own_cells && own_cell) || (loop_control.ghost_cells && !own_cell)))
350  cell_worker(dof_info.cell, info.cell);
351  }
352 
353 
369  template<int dim, int spacedim, class DOFINFO, class INFOBOX, class ASSEMBLER, class ITERATOR>
370  void loop(ITERATOR begin,
371  typename identity<ITERATOR>::type end,
372  DOFINFO &dinfo,
373  INFOBOX &info,
374  const std_cxx11::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker,
375  const std_cxx11::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker,
376  const std_cxx11::function<void (DOFINFO &, DOFINFO &,
377  typename INFOBOX::CellInfo &,
378  typename INFOBOX::CellInfo &)> &face_worker,
379  ASSEMBLER &assembler,
380  const LoopControl &lctrl = LoopControl())
381  {
382  DoFInfoBox<dim, DOFINFO> dof_info(dinfo);
383 
384  assembler.initialize_info(dof_info.cell, false);
385  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
386  {
387  assembler.initialize_info(dof_info.interior[i], true);
388  assembler.initialize_info(dof_info.exterior[i], true);
389  }
390 
391  // Loop over all cells
392 #ifdef DEAL_II_MESHWORKER_PARALLEL
393  WorkStream::run(begin, end,
394  std_cxx11::bind(&cell_action<INFOBOX, DOFINFO, dim, spacedim, ITERATOR>,
395  std_cxx11::_1, std_cxx11::_3, std_cxx11::_2,
396  cell_worker, boundary_worker, face_worker, lctrl),
397  std_cxx11::bind(&internal::assemble<dim,DOFINFO,ASSEMBLER>, std_cxx11::_1, &assembler),
398  info, dof_info);
399 #else
400  for (ITERATOR cell = begin; cell != end; ++cell)
401  {
402  cell_action<INFOBOX,DOFINFO,dim,spacedim>(cell, dof_info,
403  info, cell_worker,
404  boundary_worker, face_worker,
405  lctrl);
406  dof_info.assemble(assembler);
407  }
408 #endif
409  }
410 
411 
419  template<int dim, int spacedim, class ITERATOR, class ASSEMBLER>
420  void integration_loop(ITERATOR begin,
421  typename identity<ITERATOR>::type end,
422  DoFInfo<dim, spacedim> &dof_info,
424  const LocalIntegrator<dim, spacedim> &integrator,
425  ASSEMBLER &assembler,
426  const LoopControl &lctrl = LoopControl())
427  {
428  std_cxx11::function<void (DoFInfo<dim>&, IntegrationInfo<dim, spacedim>&)> cell_worker;
429  std_cxx11::function<void (DoFInfo<dim>&, IntegrationInfo<dim, spacedim>&)> boundary_worker;
430  std_cxx11::function<void (DoFInfo<dim> &, DoFInfo<dim> &,
432  IntegrationInfo<dim, spacedim> &)> face_worker;
433  if (integrator.use_cell)
434  cell_worker = std_cxx11::bind(&LocalIntegrator<dim, spacedim>::cell, &integrator, std_cxx11::_1, std_cxx11::_2);
435  if (integrator.use_boundary)
436  boundary_worker = std_cxx11::bind(&LocalIntegrator<dim, spacedim>::boundary, &integrator, std_cxx11::_1, std_cxx11::_2);
437  if (integrator.use_face)
438  face_worker = std_cxx11::bind(&LocalIntegrator<dim, spacedim>::face, &integrator, std_cxx11::_1, std_cxx11::_2, std_cxx11::_3, std_cxx11::_4);
439 
440  loop<dim, spacedim>
441  (begin, end,
442  dof_info,
443  box,
444  cell_worker,
445  boundary_worker,
446  face_worker,
447  assembler,
448  lctrl);
449  }
450 
451 }
452 
453 DEAL_II_NAMESPACE_CLOSE
454 
455 #endif
const types::subdomain_id invalid_subdomain_id
Definition: types.h:239
::ExceptionBase & ExcMessage(std::string arg1)
FaceOption own_faces
Definition: loop.h:121
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:264
void cell_action(ITERATOR cell, DoFInfoBox< dim, DOFINFO > &dof_info, INFOBOX &info, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std_cxx11::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, const LoopControl &loop_control)
Definition: loop.h:161
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:258
#define Assert(cond, exc)
Definition: exceptions.h:294
bool is_active_iterator(const DI &)
Definition: loop.h:43
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std_cxx11::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:370
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:453
unsigned int subdomain_id
Definition: types.h:42
const types::subdomain_id artificial_subdomain_id
Definition: types.h:255
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:248
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:1119
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:420
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:252
FaceOption faces_to_ghost
Definition: loop.h:113