17#ifndef dealii_mesh_worker_loop_h
18#define dealii_mesh_worker_loop_h
54 template <
class ACCESSOR>
61 template <
class ACCESSOR>
69 template <
int dim,
class DOFINFO,
class A>
189 template <
class INFOBOX,
class DOFINFO,
int dim,
int spacedim,
class ITERATOR>
195 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
197 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
199 const std::function<
void(DOFINFO &,
201 typename INFOBOX::CellInfo &,
202 typename INFOBOX::CellInfo &)> &face_worker,
205 const bool ignore_subdomain =
206 (cell->get_triangulation().locally_owned_subdomain() ==
210 cell->level_subdomain_id() :
211 cell->subdomain_id();
213 const bool own_cell =
215 (csid == cell->get_triangulation().locally_owned_subdomain());
222 dof_info.
cell.reinit(cell);
225 const bool integrate_cell = (cell_worker !=
nullptr);
226 const bool integrate_boundary = (boundary_worker !=
nullptr);
227 const bool integrate_interior_face = (face_worker !=
nullptr);
230 info.cell.reinit(dof_info.
cell);
237 cell_worker(dof_info.
cell, info.cell);
243 info.post_cell(dof_info);
245 if (integrate_interior_face || integrate_boundary)
246 for (
const unsigned int face_no : cell->face_indices())
248 typename ITERATOR::AccessorType::Container::face_iterator face =
250 if (cell->at_boundary(face_no) &&
251 !cell->has_periodic_neighbor(face_no))
254 if (integrate_boundary && own_cell)
257 dof_info.
interior[face_no].reinit(cell, face, face_no);
258 info.boundary.reinit(dof_info.
interior[face_no]);
259 boundary_worker(dof_info.
interior[face_no], info.boundary);
262 else if (integrate_interior_face)
266 cell->neighbor_or_periodic_neighbor(face_no);
269 if (neighbor->is_level_cell())
270 neighbid = neighbor->level_subdomain_id();
272 else if (neighbor->is_active())
273 neighbid = neighbor->subdomain_id();
275 const bool own_neighbor =
278 cell->get_triangulation().locally_owned_subdomain());
281 if (!own_cell && !own_neighbor)
285 if (own_cell && own_neighbor &&
290 if (own_cell != own_neighbor &&
297 const bool periodic_neighbor =
298 cell->has_periodic_neighbor(face_no);
300 if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no)) ||
301 (periodic_neighbor &&
302 cell->periodic_neighbor_is_coarser(face_no)))
313 const std::pair<unsigned int, unsigned int> neighbor_face_no =
315 cell->periodic_neighbor_of_coarser_periodic_neighbor(
317 cell->neighbor_of_coarser_neighbor(face_no);
318 const typename ITERATOR::AccessorType::Container::
319 face_iterator nface =
320 neighbor->face(neighbor_face_no.first);
324 dof_info.
interior[face_no].reinit(cell, face, face_no);
325 info.face.reinit(dof_info.
interior[face_no]);
326 dof_info.
exterior[face_no].reinit(neighbor,
328 neighbor_face_no.first,
329 neighbor_face_no.second);
330 info.subface.reinit(dof_info.
exterior[face_no]);
332 face_worker(dof_info.
interior[face_no],
342 neighbor->has_children())
347 "Assembling from both sides for own_faces is not "
348 "supported with hanging nodes!"));
353 Assert(cell->level() == neighbor->level(),
359 if (own_cell && own_neighbor &&
373 if (own_cell && !own_neighbor &&
378 const unsigned int neighbor_face_no =
380 cell->periodic_neighbor_face_no(face_no) :
381 cell->neighbor_face_no(face_no);
382 Assert(periodic_neighbor ||
383 neighbor->face(neighbor_face_no) == face,
388 dof_info.
interior[face_no].reinit(cell, face, face_no);
389 info.face.reinit(dof_info.
interior[face_no]);
390 dof_info.
exterior[face_no].reinit(neighbor,
394 info.neighbor.reinit(dof_info.
exterior[face_no]);
396 face_worker(dof_info.
interior[face_no],
407 info.post_faces(dof_info);
414 cell_worker(dof_info.
cell, info.cell);
443 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
445 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
447 const std::function<
void(DOFINFO &,
449 typename INFOBOX::CellInfo &,
450 typename INFOBOX::CellInfo &)> &face_worker,
451 ASSEMBLER & assembler,
456 assembler.initialize_info(dof_info.
cell,
false);
459 assembler.initialize_info(dof_info.
interior[i],
true);
460 assembler.initialize_info(dof_info.
exterior[i],
true);
467 [&cell_worker, &boundary_worker, &face_worker, &lctrl](
469 cell_action<INFOBOX, DOFINFO, dim, spacedim, ITERATOR>(cell,
478 ::internal::assemble<dim, DOFINFO, ASSEMBLER>(dinfo, &assembler);
491 template <
int dim,
int spacedim,
class ITERATOR,
class ASSEMBLER>
498 ASSEMBLER & assembler,
516 integrator.
cell(dof_info, integration_info);
522 integrator.
boundary(dof_info, integration_info);
530 integrator.
face(dof_info_1,
536 loop<dim, spacedim>(begin,
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
void assemble(ASSEMBLER &ass) const
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
virtual void cell(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
virtual void face(DoFInfo< dim, spacedim, number > &dinfo1, DoFInfo< dim, spacedim, number > &dinfo2, IntegrationInfo< dim, spacedim > &info1, IntegrationInfo< dim, spacedim > &info2) const
virtual void boundary(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
FaceOption faces_to_ghost
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void integration_loop(ITERATOR begin, std_cxx20::type_identity_t< ITERATOR > end, DoFInfo< dim, spacedim > &dof_info, IntegrationInfoBox< dim, spacedim > &box, const LocalIntegrator< dim, spacedim > &integrator, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
void loop(ITERATOR begin, std_cxx20::type_identity_t< ITERATOR > 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())
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)
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)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
bool is_active_iterator(const DI &)
const types::subdomain_id artificial_subdomain_id
const types::subdomain_id invalid_subdomain_id
typename type_identity< T >::type type_identity_t
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()