17 #ifndef dealii_mesh_worker_loop_h 18 #define dealii_mesh_worker_loop_h 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> 31 DEAL_II_NAMESPACE_OPEN
46 template <
class ACCESSOR>
52 template <
class ACCESSOR>
58 template <
int dim,
class DOFINFO,
class A>
177 template <
class INFOBOX,
class DOFINFO,
int dim,
int spacedim,
class ITERATOR>
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,
189 const bool ignore_subdomain = (cell->get_triangulation().locally_owned_subdomain()
193 ? cell->level_subdomain_id()
194 : cell->subdomain_id();
196 const bool own_cell = ignore_subdomain || (csid == cell->get_triangulation().locally_owned_subdomain());
203 dof_info.
cell.reinit(cell);
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);
211 info.cell.reinit(dof_info.
cell);
217 cell_worker(dof_info.
cell, info.cell);
223 info.post_cell(dof_info);
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)
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))
232 if (integrate_boundary && own_cell)
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);
240 else if (integrate_interior_face)
246 if (neighbor->is_level_cell())
247 neighbid = neighbor->level_subdomain_id();
249 else if (neighbor->active())
250 neighbid = neighbor->subdomain_id();
252 const bool own_neighbor = ignore_subdomain ||
253 (neighbid == cell->get_triangulation().locally_owned_subdomain());
256 if (!own_cell && !own_neighbor)
269 const bool periodic_neighbor = cell->has_periodic_neighbor(face_no);
271 if ((!periodic_neighbor && cell->neighbor_is_coarser(face_no))
272 || (periodic_neighbor && cell->periodic_neighbor_is_coarser(face_no)))
283 const std::pair<unsigned int, unsigned int> neighbor_face_no
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);
292 dof_info.
interior[face_no].reinit(cell, face, face_no);
293 info.face.reinit(dof_info.
interior[face_no]);
295 neighbor, nface, neighbor_face_no.first, neighbor_face_no.second);
296 info.subface.reinit(dof_info.
exterior[face_no]);
299 info.face, info.subface);
308 "Assembling from both sides for own_faces is not " 309 "supported with hanging nodes!"));
319 if (own_cell && own_neighbor
321 && (neighbor < cell))
333 if (own_cell && !own_neighbor
335 && (neighbid < csid))
338 const unsigned int neighbor_face_no = periodic_neighbor?
339 cell->periodic_neighbor_face_no(face_no):
340 cell->neighbor_face_no(face_no);
345 dof_info.
interior[face_no].reinit(cell, face, face_no);
346 info.face.reinit(dof_info.
interior[face_no]);
348 neighbor, neighbor->face(neighbor_face_no), neighbor_face_no);
349 info.neighbor.reinit(dof_info.
exterior[face_no]);
352 info.face, info.neighbor);
360 info.post_faces(dof_info);
366 cell_worker(dof_info.
cell, info.cell);
385 template <
int dim,
int spacedim,
class DOFINFO,
class INFOBOX,
class ASSEMBLER,
class ITERATOR>
387 typename identity<ITERATOR>::type end,
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,
400 assembler.initialize_info(dof_info.
cell,
false);
401 for (
unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
403 assembler.initialize_info(dof_info.
interior[i],
true);
404 assembler.initialize_info(dof_info.
exterior[i],
true);
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),
424 template <
int dim,
int spacedim,
class ITERATOR,
class ASSEMBLER>
426 typename identity<ITERATOR>::type end,
430 ASSEMBLER &assembler,
458 DEAL_II_NAMESPACE_CLOSE
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
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 &)
void assemble(ASSEMBLER &ass) const
unsigned int subdomain_id
const types::subdomain_id artificial_subdomain_id
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]
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)
static ::ExceptionBase & ExcInternalError()