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\}}\)
loop.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 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_loop_h
18 #define dealii_mesh_worker_loop_h
19 
20 #include <deal.II/base/config.h>
21 
24 
26 #include <deal.II/grid/tria.h>
27 
31 
32 #include <functional>
33 
35 
36 // Forward declaration
37 #ifndef DOXYGEN
38 template <typename>
39 class TriaActiveIterator;
40 #endif
41 
42 namespace internal
43 {
47  template <class DI>
48  inline bool
49  is_active_iterator(const DI &)
50  {
51  return false;
52  }
53 
54  template <class ACCESSOR>
55  inline bool
57  {
58  return true;
59  }
60 
61  template <class ACCESSOR>
62  inline bool
64  const ::FilteredIterator<TriaActiveIterator<ACCESSOR>> &)
65  {
66  return true;
67  }
68 
69  template <int dim, class DOFINFO, class A>
70  void
72  {
73  dinfo.assemble(*assembler);
74  }
75 } // namespace internal
76 
77 
78 
79 namespace MeshWorker
80 {
85  {
86  public:
91  : own_cells(true)
92  , ghost_cells(false)
95  , cells_first(true)
96  {}
97 
101  bool own_cells;
102 
108 
115  {
128  };
129 
143 
154 
160  };
161 
162 
163 
191  template <class INFOBOX, class DOFINFO, int dim, int spacedim, class ITERATOR>
192  void
194  ITERATOR cell,
195  DoFInfoBox<dim, DOFINFO> &dof_info,
196  INFOBOX & info,
197  const std::function<void(DOFINFO &, typename INFOBOX::CellInfo &)>
198  &cell_worker,
199  const std::function<void(DOFINFO &, typename INFOBOX::CellInfo &)>
200  & boundary_worker,
201  const std::function<void(DOFINFO &,
202  DOFINFO &,
203  typename INFOBOX::CellInfo &,
204  typename INFOBOX::CellInfo &)> &face_worker,
205  const LoopControl & loop_control)
206  {
207  const bool ignore_subdomain =
208  (cell->get_triangulation().locally_owned_subdomain() ==
210 
211  types::subdomain_id csid = (cell->is_level_cell()) ?
212  cell->level_subdomain_id() :
213  cell->subdomain_id();
214 
215  const bool own_cell =
216  ignore_subdomain ||
217  (csid == cell->get_triangulation().locally_owned_subdomain());
218 
219  dof_info.reset();
220 
221  if ((!ignore_subdomain) && (csid == numbers::artificial_subdomain_id))
222  return;
223 
224  dof_info.cell.reinit(cell);
225  dof_info.cell_valid = true;
226 
227  const bool integrate_cell = (cell_worker != nullptr);
228  const bool integrate_boundary = (boundary_worker != nullptr);
229  const bool integrate_interior_face = (face_worker != nullptr);
230 
231  if (integrate_cell)
232  info.cell.reinit(dof_info.cell);
233  // Execute this, if cells
234  // have to be dealt with
235  // before faces
236  if (integrate_cell && loop_control.cells_first &&
237  ((loop_control.own_cells && own_cell) ||
238  (loop_control.ghost_cells && !own_cell)))
239  cell_worker(dof_info.cell, info.cell);
240 
241  // Call the callback function in
242  // the info box to do
243  // computations between cell and
244  // face action.
245  info.post_cell(dof_info);
246 
247  if (integrate_interior_face || integrate_boundary)
248  for (const unsigned int face_no : GeometryInfo<
249  ITERATOR::AccessorType::Container::dimension>::face_indices())
250  {
251  typename ITERATOR::AccessorType::Container::face_iterator face =
252  cell->face(face_no);
253  if (cell->at_boundary(face_no) &&
254  !cell->has_periodic_neighbor(face_no))
255  {
256  // only integrate boundary faces of own cells
257  if (integrate_boundary && own_cell)
258  {
259  dof_info.interior_face_available[face_no] = true;
260  dof_info.interior[face_no].reinit(cell, face, face_no);
261  info.boundary.reinit(dof_info.interior[face_no]);
262  boundary_worker(dof_info.interior[face_no], info.boundary);
263  }
264  }
265  else if (integrate_interior_face)
266  {
267  // Interior face
269  cell->neighbor_or_periodic_neighbor(face_no);
270 
272  if (neighbor->is_level_cell())
273  neighbid = neighbor->level_subdomain_id();
274  // subdomain id is only valid for active cells
275  else if (neighbor->is_active())
276  neighbid = neighbor->subdomain_id();
277 
278  const bool own_neighbor =
279  ignore_subdomain ||
280  (neighbid ==
281  cell->get_triangulation().locally_owned_subdomain());
282 
283  // skip all faces between two ghost cells
284  if (!own_cell && !own_neighbor)
285  continue;
286 
287  // skip if the user doesn't want faces between own cells
288  if (own_cell && own_neighbor &&
289  loop_control.own_faces == LoopControl::never)
290  continue;
291 
292  // skip face to ghost
293  if (own_cell != own_neighbor &&
294  loop_control.faces_to_ghost == LoopControl::never)
295  continue;
296 
297  // Deal with refinement edges from the refined side. Assuming
298  // one-irregular meshes, this situation should only occur if both
299  // cells are active.
300  const bool periodic_neighbor =
301  cell->has_periodic_neighbor(face_no);
302 
303  if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no)) ||
304  (periodic_neighbor &&
305  cell->periodic_neighbor_is_coarser(face_no)))
306  {
307  Assert(cell->is_active(), ExcInternalError());
308  Assert(neighbor->is_active(), ExcInternalError());
309 
310  // skip if only one processor needs to assemble the face
311  // to a ghost cell and the fine cell is not ours.
312  if (!own_cell &&
313  loop_control.faces_to_ghost == LoopControl::one)
314  continue;
315 
316  const std::pair<unsigned int, unsigned int> neighbor_face_no =
317  periodic_neighbor ?
318  cell->periodic_neighbor_of_coarser_periodic_neighbor(
319  face_no) :
320  cell->neighbor_of_coarser_neighbor(face_no);
321  const typename ITERATOR::AccessorType::Container::
322  face_iterator nface =
323  neighbor->face(neighbor_face_no.first);
324 
325  dof_info.interior_face_available[face_no] = true;
326  dof_info.exterior_face_available[face_no] = true;
327  dof_info.interior[face_no].reinit(cell, face, face_no);
328  info.face.reinit(dof_info.interior[face_no]);
329  dof_info.exterior[face_no].reinit(neighbor,
330  nface,
331  neighbor_face_no.first,
332  neighbor_face_no.second);
333  info.subface.reinit(dof_info.exterior[face_no]);
334 
335  face_worker(dof_info.interior[face_no],
336  dof_info.exterior[face_no],
337  info.face,
338  info.subface);
339  }
340  else
341  {
342  // If iterator is active and neighbor is refined, skip
343  // internal face.
344  if (::internal::is_active_iterator(cell) &&
345  neighbor->has_children())
346  {
347  Assert(
348  loop_control.own_faces != LoopControl::both,
349  ExcMessage(
350  "Assembling from both sides for own_faces is not "
351  "supported with hanging nodes!"));
352  continue;
353  }
354 
355  // Now neighbor is on same level, double-check this:
356  Assert(cell->level() == neighbor->level(),
357  ExcInternalError());
358 
359  // If we own both cells only do faces from one side (unless
360  // LoopControl says otherwise). Here, we rely on cell
361  // comparison that will look at cell->index().
362  if (own_cell && own_neighbor &&
363  loop_control.own_faces == LoopControl::one &&
364  (neighbor < cell))
365  continue;
366 
367  // independent of loop_control.faces_to_ghost,
368  // we only look at faces to ghost on the same level once
369  // (only where own_cell=true and own_neighbor=false)
370  if (!own_cell)
371  continue;
372 
373  // now only one processor assembles faces_to_ghost. We let the
374  // processor with the smaller (level-)subdomain id assemble
375  // the face.
376  if (own_cell && !own_neighbor &&
377  loop_control.faces_to_ghost == LoopControl::one &&
378  (neighbid < csid))
379  continue;
380 
381  const unsigned int neighbor_face_no =
382  periodic_neighbor ?
383  cell->periodic_neighbor_face_no(face_no) :
384  cell->neighbor_face_no(face_no);
385  Assert(periodic_neighbor ||
386  neighbor->face(neighbor_face_no) == face,
387  ExcInternalError());
388  // Regular interior face
389  dof_info.interior_face_available[face_no] = true;
390  dof_info.exterior_face_available[face_no] = true;
391  dof_info.interior[face_no].reinit(cell, face, face_no);
392  info.face.reinit(dof_info.interior[face_no]);
393  dof_info.exterior[face_no].reinit(neighbor,
394  neighbor->face(
395  neighbor_face_no),
396  neighbor_face_no);
397  info.neighbor.reinit(dof_info.exterior[face_no]);
398 
399  face_worker(dof_info.interior[face_no],
400  dof_info.exterior[face_no],
401  info.face,
402  info.neighbor);
403  }
404  }
405  } // faces
406  // Call the callback function in
407  // the info box to do
408  // computations between face and
409  // cell action.
410  info.post_faces(dof_info);
411 
412  // Execute this, if faces
413  // have to be handled first
414  if (integrate_cell && !loop_control.cells_first &&
415  ((loop_control.own_cells && own_cell) ||
416  (loop_control.ghost_cells && !own_cell)))
417  cell_worker(dof_info.cell, info.cell);
418  }
419 
420 
436  template <int dim,
437  int spacedim,
438  class DOFINFO,
439  class INFOBOX,
440  class ASSEMBLER,
441  class ITERATOR>
442  void
443  loop(ITERATOR begin,
444  typename identity<ITERATOR>::type end,
445  DOFINFO & dinfo,
446  INFOBOX & info,
447  const std::function<void(DOFINFO &, typename INFOBOX::CellInfo &)>
448  &cell_worker,
449  const std::function<void(DOFINFO &, typename INFOBOX::CellInfo &)>
450  & boundary_worker,
451  const std::function<void(DOFINFO &,
452  DOFINFO &,
453  typename INFOBOX::CellInfo &,
454  typename INFOBOX::CellInfo &)> &face_worker,
455  ASSEMBLER & assembler,
456  const LoopControl &lctrl = LoopControl())
457  {
458  DoFInfoBox<dim, DOFINFO> dof_info(dinfo);
459 
460  assembler.initialize_info(dof_info.cell, false);
461  for (unsigned int i : GeometryInfo<dim>::face_indices())
462  {
463  assembler.initialize_info(dof_info.interior[i], true);
464  assembler.initialize_info(dof_info.exterior[i], true);
465  }
466 
467  // Loop over all cells
469  begin,
470  end,
471  [&cell_worker, &boundary_worker, &face_worker, &lctrl](
472  ITERATOR cell, INFOBOX &info, DoFInfoBox<dim, DOFINFO> &dof_info) {
473  cell_action<INFOBOX, DOFINFO, dim, spacedim, ITERATOR>(cell,
474  dof_info,
475  info,
476  cell_worker,
477  boundary_worker,
478  face_worker,
479  lctrl);
480  },
481  [&assembler](const MeshWorker::DoFInfoBox<dim, DOFINFO> &dinfo) {
482  ::internal::assemble<dim, DOFINFO, ASSEMBLER>(dinfo, &assembler);
483  },
484  info,
485  dof_info);
486  }
487 
488 
496  template <int dim, int spacedim, class ITERATOR, class ASSEMBLER>
497  void
499  typename identity<ITERATOR>::type end,
500  DoFInfo<dim, spacedim> & dof_info,
502  const LocalIntegrator<dim, spacedim> &integrator,
503  ASSEMBLER & assembler,
504  const LoopControl & lctrl = LoopControl())
505  {
506  std::function<void(DoFInfo<dim, spacedim> &,
508  cell_worker;
509  std::function<void(DoFInfo<dim, spacedim> &,
511  boundary_worker;
512  std::function<void(DoFInfo<dim, spacedim> &,
516  face_worker;
517  if (integrator.use_cell)
518  cell_worker =
519  [&integrator](DoFInfo<dim, spacedim> & dof_info,
520  IntegrationInfo<dim, spacedim> &integration_info) {
521  integrator.cell(dof_info, integration_info);
522  };
523  if (integrator.use_boundary)
524  boundary_worker =
525  [&integrator](DoFInfo<dim, spacedim> & dof_info,
526  IntegrationInfo<dim, spacedim> &integration_info) {
527  integrator.boundary(dof_info, integration_info);
528  };
529  if (integrator.use_face)
530  face_worker =
531  [&integrator](DoFInfo<dim, spacedim> & dof_info_1,
532  DoFInfo<dim, spacedim> & dof_info_2,
533  IntegrationInfo<dim, spacedim> &integration_info_1,
534  IntegrationInfo<dim, spacedim> &integration_info_2) {
535  integrator.face(dof_info_1,
536  dof_info_2,
537  integration_info_1,
538  integration_info_2);
539  };
540 
541  loop<dim, spacedim>(begin,
542  end,
543  dof_info,
544  box,
545  cell_worker,
546  boundary_worker,
547  face_worker,
548  assembler,
549  lctrl);
550  }
551 
552 } // namespace MeshWorker
553 
555 
556 #endif
identity::type
T type
Definition: template_constraints.h:270
MeshWorker::LoopControl::never
@ never
Definition: loop.h:119
MeshWorker::DoFInfoBox::interior_face_available
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:276
MeshWorker::DoFInfoBox::cell
DOFINFO cell
Definition: dof_info.h:262
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::LocalIntegrator::boundary
virtual void boundary(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
Definition: mesh_worker.cc:91
tria.h
MeshWorker
Definition: assemble_flags.h:30
MeshWorker::DoFInfoBox::interior
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:266
MeshWorker::LoopControl::faces_to_ghost
FaceOption faces_to_ghost
Definition: loop.h:142
GeometryInfo
Definition: geometry_info.h:1224
MeshWorker::DoFInfoBox::exterior
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:270
integration_info.h
MeshWorker::DoFInfoBox::reset
void reset()
Definition: dof_info.h:481
MeshWorker::LocalIntegrator::use_cell
bool use_cell
Definition: local_integrator.h:105
MeshWorker::IntegrationInfo
Definition: integration_info.h:78
MeshWorker::LocalIntegrator::use_face
bool use_face
Definition: local_integrator.h:117
MeshWorker::DoFInfoBox::cell_valid
bool cell_valid
Definition: dof_info.h:287
MeshWorker::LoopControl
Definition: loop.h:84
internal::is_active_iterator
bool is_active_iterator(const DI &)
Definition: loop.h:49
MeshWorker::LocalIntegrator::cell
virtual void cell(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
Definition: mesh_worker.cc:81
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
MeshWorker::LoopControl::own_cells
bool own_cells
Definition: loop.h:101
work_stream.h
MeshWorker::LocalIntegrator::use_boundary
bool use_boundary
Definition: local_integrator.h:111
MeshWorker::LocalIntegrator
Definition: local_integrator.h:58
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MeshWorker::LoopControl::own_faces
FaceOption own_faces
Definition: loop.h:153
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
MeshWorker::LoopControl::both
@ both
Definition: loop.h:127
dof_info.h
TriaActiveIterator
Definition: tria_iterator.h:759
MeshWorker::LoopControl::cells_first
bool cells_first
Definition: loop.h:159
LAPACKSupport::A
static const char A
Definition: lapack_support.h:155
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
MeshWorker::LoopControl::ghost_cells
bool ghost_cells
Definition: loop.h:107
MeshWorker::LocalIntegrator::face
virtual void face(DoFInfo< dim, spacedim, number > &dinfo1, DoFInfo< dim, spacedim, number > &dinfo2, IntegrationInfo< dim, spacedim > &info1, IntegrationInfo< dim, spacedim > &info2) const
Definition: mesh_worker.cc:101
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
MeshWorker::DoFInfoBox::assemble
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:495
internal::assemble
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition: loop.h:71
MeshWorker::IntegrationInfoBox
Definition: integration_info.h:299
MeshWorker::integration_loop
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:498
local_integrator.h
MeshWorker::DoFInfoBox::exterior_face_available
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:282
template_constraints.h
config.h
MeshWorker::LoopControl::FaceOption
FaceOption
Definition: loop.h:114
MeshWorker::LoopControl::one
@ one
Definition: loop.h:123
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
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MeshWorker::loop
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:443
TriaIterator
Definition: tria_iterator.h:578
MeshWorker::LoopControl::LoopControl
LoopControl()
Definition: loop.h:90
MeshWorker::DoFInfo
Definition: dof_info.h:76
filtered_iterator.h
MeshWorker::DoFInfoBox
Definition: dof_info.h:223