14#ifndef agglomeration_accessor_h
15#define agglomeration_accessor_h
41template <
int dim,
int spacedim = dim>
49 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>;
101 const std::vector<typename Triangulation<dim>::active_face_iterator> &
178 inline const std::vector<types::global_cell_index> &
279template <
int dim,
int spacedim>
284 unsigned int n_neighbors = 0;
285 for (
const auto &f : cell->face_indices())
287 const auto &neighboring_cell = cell->neighbor(f);
288 if ((cell->face(f)->at_boundary()) ||
289 (neighboring_cell->is_active() &&
290 !
handler->are_cells_agglomerated(cell, neighboring_cell)))
300template <
int dim,
int spacedim>
305 ExcMessage(
"You cannot pass a slave cell."));
311template <
int dim,
int spacedim>
322 const unsigned int sender_rank =
master_cell->subdomain_id();
324 const CellId &master_id_ghosted_neighbor =
325 handler->recv_ghosted_master_id.at(sender_rank)
333 master_id_ghosted_neighbor,
354 handler->polytope_cache.cell_face_at_boundary.at({polytope_index, f})
358 if (neigh->is_locally_owned())
386 const CellId &master_id_neighbor =
391 return {neigh, master_id_neighbor,
handler};
402template <
int dim,
int spacedim>
405 const unsigned int f)
const
410 const auto &neigh_polytope =
416 unsigned int n_faces_agglomerated_neighbor;
419 if (neigh_polytope->is_locally_owned())
421 n_faces_agglomerated_neighbor = neigh_polytope->n_faces();
429 const CellId &master_id_neighbor = neigh_polytope->id();
432 const unsigned int sender_rank = neigh_polytope->subdomain_id();
436 n_faces_agglomerated_neighbor =
437 handler->recv_n_faces.at(sender_rank).at(master_id_neighbor);
442 for (
unsigned int f_out = 0; f_out < n_faces_agglomerated_neighbor;
447 if (neigh_polytope->neighbor(f_out)->id() ==
present_id)
461template <
int dim,
int spacedim>
467template <
int dim,
int spacedim>
473 if (&(*
handler->master_cells_container.end()) == std::addressof(cell))
491template <
int dim,
int spacedim>
494 const CellId &master_cell_id,
497 Assert(neigh_cell->is_ghost(), ExcInternalError());
510template <
int dim,
int spacedim>
513 std::vector<types::global_dof_index> &dof_indices)
const
515 Assert(dof_indices.size() > 0,
517 "The vector of DoFs indices must be already properly resized."));
523 master_cell_dh->get_dof_indices(dof_indices);
527 const std::vector<types::global_dof_index> &recv_dof_indices =
530 std::copy(recv_dof_indices.cbegin(),
531 recv_dof_indices.cend(),
532 dof_indices.begin());
538template <
int dim,
int spacedim>
544 return agglomeration;
549template <
int dim,
int spacedim>
550inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
558template <
int dim,
int spacedim>
563 ExcMessage(
"The present function cannot be called for slave cells."));
568 const auto &bdary_pts =
570 return (bdary_pts.second - bdary_pts.first).norm();
581template <
int dim,
int spacedim>
593template <
int dim,
int spacedim>
598 ExcMessage(
"The present function cannot be called for slave cells."));
612template <
int dim,
int spacedim>
620 if (present_index < handler->master_cells_container.size())
630template <
int dim,
int spacedim>
641template <
int dim,
int spacedim>
649template <
int dim,
int spacedim>
654 return !(*
this == other);
659template <
int dim,
int spacedim>
668template <
int dim,
int spacedim>
674 return master_cell->as_dof_handler_iterator(dof_handler);
679template <
int dim,
int spacedim>
689template <
int dim,
int spacedim>
699template <
int dim,
int spacedim>
704 unsigned int n_neighbors = 0;
705 for (
const auto &cell : agglomeration)
712template <
int dim,
int spacedim>
718 const unsigned int sender_rank =
master_cell->subdomain_id();
741 "This function should not be called for a slave cell."));
746 return handler->at_boundary(cell_dh, f);
752template <
int dim,
int spacedim>
761template <
int dim,
int spacedim>
770template <
int dim,
int spacedim>
777template <
int dim,
int spacedim>
778inline const std::vector<types::global_cell_index> &
781 Assert(!
handler->parent_child_info.empty(), ExcInternalError());
782 return handler->parent_child_info.at(
unsigned int n_agglomerated_faces() const
AgglomerationContainer get_agglomerate() const
types::global_cell_index index() const
~AgglomerationAccessor()=default
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > AgglomerationContainer
unsigned int n_agglomerated_faces_per_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
unsigned int n_background_cells() const
unsigned int neighbor_of_agglomerated_neighbor(const unsigned int f) const
AgglomerationAccessor(const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell, const AgglomerationHandler< dim, spacedim > *ah)
const std::vector< typename Triangulation< dim >::active_face_iterator > & polytope_boundary() const
bool is_locally_owned() const
const AgglomerationIterator< dim, spacedim > neighbor(const unsigned int f) const
bool at_boundary(const unsigned int f) const
const BoundingBox< dim > & get_bounding_box() const
friend class AgglomerationIterator
const std::vector< types::global_cell_index > & children() const
bool operator!=(const AgglomerationAccessor< dim, spacedim > &other) const
types::global_cell_index present_index
bool operator==(const AgglomerationAccessor< dim, spacedim > &other) const
AgglomerationHandler< dim, spacedim > * handler
DoFHandler< dim, spacedim >::active_cell_iterator as_dof_handler_iterator(const DoFHandler< dim, spacedim > &dof_handler) const
void get_dof_indices(std::vector< types::global_dof_index > &) const
types::subdomain_id subdomain_id() const
AgglomerationAccessor(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const CellId &cell_id, const AgglomerationHandler< dim, spacedim > *ah)
const AgglomerationContainer & get_slaves() const
types::subdomain_id present_subdomain_id
Triangulation< dim, spacedim >::active_cell_iterator master_cell
unsigned int n_faces() const
typename AgglomerationAccessor< dim, spacedim >::AgglomerationContainer AgglomerationContainer
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
constexpr unsigned int invalid_unsigned_int
constexpr types::subdomain_id invalid_subdomain_id
unsigned int subdomain_id
unsigned int global_cell_index