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>
191 template <
class INFOBOX,
class DOFINFO,
int dim,
int spacedim,
class ITERATOR>
197 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
199 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
201 const std::function<
void(DOFINFO &,
203 typename INFOBOX::CellInfo &,
204 typename INFOBOX::CellInfo &)> &face_worker,
207 const bool ignore_subdomain =
208 (cell->get_triangulation().locally_owned_subdomain() ==
212 cell->level_subdomain_id() :
213 cell->subdomain_id();
215 const bool own_cell =
217 (csid == cell->get_triangulation().locally_owned_subdomain());
224 dof_info.
cell.reinit(cell);
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);
232 info.cell.reinit(dof_info.
cell);
239 cell_worker(dof_info.
cell, info.cell);
245 info.post_cell(dof_info);
247 if (integrate_interior_face || integrate_boundary)
249 ITERATOR::AccessorType::Container::dimension>::face_indices())
251 typename ITERATOR::AccessorType::Container::face_iterator face =
253 if (cell->at_boundary(face_no) &&
254 !cell->has_periodic_neighbor(face_no))
257 if (integrate_boundary && own_cell)
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);
265 else if (integrate_interior_face)
269 cell->neighbor_or_periodic_neighbor(face_no);
272 if (neighbor->is_level_cell())
273 neighbid = neighbor->level_subdomain_id();
275 else if (neighbor->is_active())
276 neighbid = neighbor->subdomain_id();
278 const bool own_neighbor =
281 cell->get_triangulation().locally_owned_subdomain());
284 if (!own_cell && !own_neighbor)
288 if (own_cell && own_neighbor &&
289 loop_control.
own_faces == LoopControl::never)
293 if (own_cell != own_neighbor &&
300 const bool periodic_neighbor =
301 cell->has_periodic_neighbor(face_no);
303 if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no)) ||
304 (periodic_neighbor &&
305 cell->periodic_neighbor_is_coarser(face_no)))
316 const std::pair<unsigned int, unsigned int> neighbor_face_no =
318 cell->periodic_neighbor_of_coarser_periodic_neighbor(
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);
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,
331 neighbor_face_no.first,
332 neighbor_face_no.second);
333 info.subface.reinit(dof_info.
exterior[face_no]);
335 face_worker(dof_info.
interior[face_no],
345 neighbor->has_children())
348 loop_control.
own_faces != LoopControl::both,
350 "Assembling from both sides for own_faces is not " 351 "supported with hanging nodes!"));
356 Assert(cell->level() == neighbor->level(),
362 if (own_cell && own_neighbor &&
376 if (own_cell && !own_neighbor &&
381 const unsigned int neighbor_face_no =
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,
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,
397 info.neighbor.reinit(dof_info.
exterior[face_no]);
399 face_worker(dof_info.
interior[face_no],
410 info.post_faces(dof_info);
417 cell_worker(dof_info.
cell, info.cell);
447 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
449 const std::function<
void(DOFINFO &,
typename INFOBOX::CellInfo &)>
451 const std::function<
void(DOFINFO &,
453 typename INFOBOX::CellInfo &,
454 typename INFOBOX::CellInfo &)> &face_worker,
455 ASSEMBLER & assembler,
460 assembler.initialize_info(dof_info.
cell,
false);
463 assembler.initialize_info(dof_info.
interior[i],
true);
464 assembler.initialize_info(dof_info.
exterior[i],
true);
471 [&cell_worker, &boundary_worker, &face_worker, &lctrl](
473 cell_action<INFOBOX, DOFINFO, dim, spacedim, ITERATOR>(cell,
482 ::internal::assemble<dim, DOFINFO, ASSEMBLER>(dinfo, &assembler);
496 template <
int dim,
int spacedim,
class ITERATOR,
class ASSEMBLER>
503 ASSEMBLER & assembler,
506 std::function<void(DoFInfo<dim, spacedim> &,
509 std::function<void(DoFInfo<dim, spacedim> &,
512 std::function<void(DoFInfo<dim, spacedim> &,
515 IntegrationInfo<dim, spacedim> &)>
519 [&integrator](DoFInfo<dim, spacedim> & dof_info,
520 IntegrationInfo<dim, spacedim> &integration_info) {
521 integrator.
cell(dof_info, integration_info);
525 [&integrator](DoFInfo<dim, spacedim> & dof_info,
526 IntegrationInfo<dim, spacedim> &integration_info) {
527 integrator.
boundary(dof_info, integration_info);
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,
541 loop<dim, spacedim>(
begin,
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())
const types::subdomain_id invalid_subdomain_id
virtual void cell(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
static const types::blas_int one
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
static ::ExceptionBase & ExcMessage(std::string arg1)
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
#define Assert(cond, exc)
bool is_active_iterator(const DI &)
#define DEAL_II_NAMESPACE_CLOSE
VectorType::value_type * end(VectorType &V)
void assemble(ASSEMBLER &ass) const
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
const types::subdomain_id artificial_subdomain_id
#define DEAL_II_NAMESPACE_OPEN
VectorType::value_type * begin(VectorType &V)
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
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 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())
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
virtual void face(DoFInfo< dim, spacedim, number > &dinfo1, DoFInfo< dim, spacedim, number > &dinfo2, IntegrationInfo< dim, spacedim > &info1, IntegrationInfo< dim, spacedim > &info2) const
FaceOption faces_to_ghost
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)
virtual void boundary(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
static ::ExceptionBase & ExcInternalError()