12#ifndef agglomeration_handler_h
13#define agglomeration_handler_h
61template <
int dim,
int spacedim>
86 template <
int dim,
int spacedim>
116 mutable std::set<std::pair<types::global_cell_index, unsigned int>>
120 mutable std::set<std::pair<CellId, unsigned int>>
134 std::pair<types::global_cell_index, unsigned int>,
143 mutable std::map<std::pair<CellId, unsigned int>,
CellId>
157 std::pair<CellId, CellId>,
159 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
170template <
int dim,
int spacedim = dim>
257 const UpdateFlags &flags = UpdateFlags::update_default,
259 const UpdateFlags &face_flags = UpdateFlags::update_default);
267 template <
typename SparsityPatternType,
typename Number =
double>
270 SparsityPatternType &sparsity_pattern,
272 const bool keep_constrained_dofs =
true,
306 inline const std::vector<BoundingBox<dim>> &
321 inline decltype(
auto)
327 template <
typename CellIterator>
335 inline const std::vector<
360 template <
class StreamType>
364 for (
const auto &cell :
tria->active_cell_iterators())
365 out <<
"Cell with index: " << cell->active_cell_index()
366 <<
" has associated value: "
406 const unsigned int face_index)
const;
413 std::pair<const FEValuesBase<dim, spacedim> &,
417 const unsigned int local_in,
418 const unsigned int local_outside)
const;
444 const unsigned int f)
const;
449 inline
types::global_dof_index
459 inline const
std::vector<typename
Triangulation<dim>::active_face_iterator> &
461 const typename
Triangulation<dim>::active_cell_iterator &cell);
511 mutable
std::map<
types::subdomain_id,
520 inline const typename
DoFHandler<dim, spacedim>::active_cell_iterator
526 template <typename RtreeType>
536 const
std::unique_ptr<
GridTools::Cache<dim, spacedim>> &cache_tria);
541 const typename
Triangulation<dim, spacedim>::active_cell_iterator
606 const unsigned int face_number,
608 &agglo_isv_ptr)
const;
615 template <
typename CellIterator>
662 std::vector<typename Triangulation<dim>::active_face_iterator>>
672 std::vector<BoundingBox<spacedim>>
bboxes;
678 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
681 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
686 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
689 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
698 std::map<CellId, std::map<unsigned int, CellId>>>
702 std::map<CellId, std::map<unsigned int, CellId>>>
707 std::map<CellId, std::map<unsigned int, bool>>>
711 std::map<CellId, std::map<unsigned int, bool>>>
715 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
718 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
723 std::map<CellId, std::vector<types::global_dof_index>>>
727 std::map<CellId, std::vector<types::global_dof_index>>>
733 std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>>
739 std::map<std::pair<CellId, unsigned int>, std::vector<double>>>
745 std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>>
751 std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>>
755 std::map<std::pair<CellId, unsigned int>,
756 std::vector<std::vector<Tensor<1, spacedim>>>>>
773 std::unique_ptr<FiniteElement<dim>>
fe;
797 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
800 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
803 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
827 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
834 std::vector<typename Triangulation<dim>::active_cell_iterator>
857 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
870 std::map<std::pair<types::global_cell_index, types::global_cell_index>,
871 std::vector<types::global_cell_index>>
880template <
int dim,
int spacedim>
889template <
int dim,
int spacedim>
898template <
int dim,
int spacedim>
907template <
int dim,
int spacedim>
915template <
int dim,
int spacedim>
916inline const std::vector<BoundingBox<dim>> &
924template <
int dim,
int spacedim>
929 return master2polygon.at(cell->active_cell_index());
934template <
int dim,
int spacedim>
938 return polytope_cache.interface;
943template <
int dim,
int spacedim>
947 return master_slave_relationships;
952template <
int dim,
int spacedim>
953inline std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
958 Assert(is_master_cell(master_cell), ExcInternalError());
959 auto agglomeration = get_slaves_of_idx(master_cell->active_cell_index());
960 agglomeration.push_back(master_cell);
961 return agglomeration;
966template <
int dim,
int spacedim>
975template <
int dim,
int spacedim>
976inline const std::vector<
981 return master2slaves.at(idx);
986template <
int dim,
int spacedim>
987template <
typename CellIterator>
990 const CellIterator &cell)
const
992 return master_slave_relationships[cell->global_active_cell_index()] == -1;
1001template <
int dim,
int spacedim>
1002template <
typename CellIterator>
1005 const CellIterator &cell)
const
1007 return master_slave_relationships[cell->global_active_cell_index()] >= 0;
1012template <
int dim,
int spacedim>
1016 const unsigned int face_index)
const
1018 Assert(!is_slave_cell(cell),
1019 ExcMessage(
"This function should not be called for a slave cell."));
1021 return polytope_cache.cell_face_at_boundary
1022 .at({master2polygon.at(cell->active_cell_index()), face_index})
1027template <
int dim,
int spacedim>
1031 return fe->n_dofs_per_cell();
1036template <
int dim,
int spacedim>
1040 return agglo_dh.n_dofs();
1045template <
int dim,
int spacedim>
1046inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
1050 return polygon_boundary[cell];
1055template <
int dim,
int spacedim>
1060 return master_slave_relationships_iterators.at(cell->active_cell_index());
1065template <
int dim,
int spacedim>
1070 auto idx = master_slave_relationships[cell->global_active_cell_index()];
1072 return cell->global_active_cell_index();
1079template <
int dim,
int spacedim>
1090 return (get_master_idx_of_cell(cell) == get_master_idx_of_cell(other_cell));
1095template <
int dim,
int spacedim>
1099 return n_agglomerations;
1104template <
int dim,
int spacedim>
1109 return master_cells_container[polytope_index]->as_dof_handler_iterator(
1115template <
int dim,
int spacedim>
1119 Assert(n_agglomerations > 0,
1120 ExcMessage(
"No agglomeration has been performed."));
1121 return {*master_cells_container.begin(),
this};
1126template <
int dim,
int spacedim>
1130 Assert(n_agglomerations > 0,
1131 ExcMessage(
"No agglomeration has been performed."));
1132 return {*master_cells_container.begin(),
this};
1137template <
int dim,
int spacedim>
1141 Assert(n_agglomerations > 0,
1142 ExcMessage(
"No agglomeration has been performed."));
1143 return {*master_cells_container.end(),
this};
1148template <
int dim,
int spacedim>
1152 Assert(n_agglomerations > 0,
1153 ExcMessage(
"No agglomeration has been performed."));
1154 return {*master_cells_container.end(),
this};
1159template <
int dim,
int spacedim>
1163 Assert(n_agglomerations > 0,
1164 ExcMessage(
"No agglomeration has been performed."));
1165 return {master_cells_container.back(),
this};
1170template <
int dim,
int spacedim>
1180template <
int dim,
int spacedim>
1181template <
typename RtreeType>
unsigned int n_dofs_per_cell() const noexcept
AgglomerationHandler()=default
types::global_cell_index cell_to_polytope_index(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
std::map< types::subdomain_id, std::map< CellId, unsigned int > > local_n_faces
agglomeration_iterator end() const
const DoFHandler< dim, spacedim > & get_dof_handler() const
std::unique_ptr< GridTools::Cache< dim, spacedim > > cached_tria
boost::signals2::connection tria_listener
const std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > & get_slaves_of_idx(types::global_cell_index idx) const
agglomeration_iterator last()
LinearAlgebra::distributed::Vector< float > master_slave_relationships
types::global_cell_index get_master_idx_of_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Tensor< 1, spacedim > > > > local_normals
const FEValues< dim, spacedim > & reinit(const AgglomerationIterator< dim, spacedim > &polytope) const
std::map< types::subdomain_id, std::map< CellId, CellId > > recv_cell_ids_neigh_cell
FE_Nothing< dim, spacedim > dummy_fe
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_bdary
bool is_slave_cell(const CellIterator &cell) const
LinearAlgebra::distributed::Vector< double > euler_vector
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, bool > > > local_bdary_info
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< double > > > > local_values
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv
std::unique_ptr< ScratchData > standard_scratch
Quadrature< dim > agglomeration_quad
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< double > > > recv_jxws
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Tensor< 1, spacedim > > > > recv_normals
Quadrature< dim - 1 > agglomeration_face_quad
std::unique_ptr< FEValues< dim, spacedim > > no_values
const DoFHandler< dim, spacedim >::active_cell_iterator polytope_to_dh_iterator(const types::global_cell_index polytope_index) const
std::map< types::subdomain_id, std::map< CellId, unsigned int > > recv_n_faces
std::vector< BoundingBox< spacedim > > bboxes
void exchange_interface_values()
unsigned int present_extraction_level
types::global_dof_index n_dofs() const noexcept
std::map< types::subdomain_id, std::map< CellId, BoundingBox< dim > > > recv_ghosted_bbox
bool are_cells_agglomerated(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const typename Triangulation< dim, spacedim >::active_cell_iterator &other_cell) const
unsigned int n_agglomerations
Quadrature< dim > agglomerated_quadrature(const AgglomerationContainer &cells, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
AgglomerationIterator< dim, spacedim > agglomeration_iterator
std::map< const typename Triangulation< dim, spacedim >::active_cell_iterator, std::vector< typename Triangulation< dim >::active_face_iterator > > polygon_boundary
ObserverPointer< const Triangulation< dim, spacedim > > tria
agglomeration_iterator begin() const
std::unordered_map< types::global_cell_index, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > > master2slaves
std::vector< typename Triangulation< dim >::active_cell_iterator > master_disconnected
void create_bounding_box(const AgglomerationContainer &polytope)
std::map< types::global_cell_index, types::global_cell_index > master2polygon
DoFHandler< dim, spacedim > output_dh
typename AgglomerationIterator< dim, spacedim >::AgglomerationContainer AgglomerationContainer
std::map< std::pair< types::global_cell_index, types::global_cell_index >, std::vector< types::global_cell_index > > parent_child_info
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Point< spacedim > > > > recv_qpoints
std::pair< const FEValuesBase< dim, spacedim > &, const FEValuesBase< dim, spacedim > & > reinit_interface(const AgglomerationIterator< dim, spacedim > &polytope_in, const AgglomerationIterator< dim, spacedim > &neigh_polytope, const unsigned int local_in, const unsigned int local_outside) const
const MPI_Comm communicator
std::map< types::subdomain_id, std::map< CellId, std::vector< types::global_dof_index > > > recv_ghost_dofs
void setup_connectivity_of_agglomeration()
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< Tensor< 1, spacedim > > > > > recv_gradients
const Triangulation< dim, spacedim > & get_triangulation() const
void create_agglomeration_sparsity_pattern(SparsityPatternType &sparsity_pattern, const AffineConstraints< Number > &constraints=AffineConstraints< Number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
const Mapping< dim > & get_mapping() const
unsigned int n_agglomerates() const
const std::vector< typename Triangulation< dim >::active_face_iterator > & polytope_boundary(const typename Triangulation< dim >::active_cell_iterator &cell)
void setup_ghost_polytopes()
std::map< types::global_cell_index, typename Triangulation< dim, spacedim >::active_cell_iterator > master_slave_relationships_iterators
const std::vector< BoundingBox< dim > > & get_local_bboxes() const
const LinearAlgebra::distributed::Vector< float > & get_relationships() const
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > get_agglomerate(const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
unsigned int n_agglomerated_faces_per_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
const UpdateFlags internal_agglomeration_face_flags
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< Tensor< 1, spacedim > > > > > local_gradients
void distribute_agglomerated_dofs(const FiniteElement< dim > &fe_space)
std::vector< types::global_cell_index > number_of_agglomerated_faces
const UpdateFlags internal_agglomeration_flags
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, bool > > > recv_bdary_info
void connect_to_tria_signals()
const FiniteElement< dim, spacedim > & get_fe() const
std::unique_ptr< FEFaceValues< dim, spacedim > > no_face_values
const FEValuesBase< dim, spacedim > & reinit_master(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_number, std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > &agglo_isv_ptr) const
std::map< types::subdomain_id, std::map< CellId, std::vector< types::global_dof_index > > > local_ghost_dofs
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > master_cells_container
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< double > > > > recv_values
void connect_hierarchy(const CellsAgglomerator< dim, RtreeType > &agglomerator)
bool at_boundary(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int f) const
bool is_master_cell(const CellIterator &cell) const
IteratorRange< agglomeration_iterator > polytope_iterators() const
decltype(auto) get_interface() const
void initialize_agglomeration_data(const std::unique_ptr< GridTools::Cache< dim, spacedim > > &cache_tria)
void print_agglomeration(StreamType &out)
agglomeration_iterator define_agglomerate(const AgglomerationContainer &cells)
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< double > > > local_jxws
void initialize_fe_values(const Quadrature< dim > &cell_quadrature=QGauss< dim >(1), const UpdateFlags &flags=UpdateFlags::update_default, const Quadrature< dim - 1 > &face_quadrature=QGauss< dim - 1 >(1), const UpdateFlags &face_flags=UpdateFlags::update_default)
agglomeration_iterator end()
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Point< spacedim > > > > local_qpoints
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_neigh
void define_agglomerate_with_check(const AgglomerationContainer &cells)
const MappingBox< dim > & get_agglomeration_mapping() const
std::unique_ptr< ScratchData > agglomerated_scratch
UpdateFlags agglomeration_face_flags
hp::FECollection< dim, spacedim > fe_collection
std::map< types::subdomain_id, std::map< CellId, BoundingBox< dim > > > local_ghosted_bbox
DoFHandler< dim, spacedim > agglo_dh
void initialize_hp_structure()
Triangulation< dim, spacedim >::active_cell_iterator & is_slave_cell_of(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< MappingBox< dim > > box_mapping
agglomeration_iterator begin()
std::unique_ptr< FiniteElement< dim > > fe
ObserverPointer< const Mapping< dim, spacedim > > mapping
void update_agglomerate(AgglomerationContainer &polytope, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell)
double get_mesh_size() const
UpdateFlags agglomeration_flags
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, CellId > > > recv_ghosted_master_id
std::map< types::subdomain_id, std::map< CellId, CellId > > local_cell_ids_neigh_cell
internal::PolytopeCache< dim, spacedim > polytope_cache
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, CellId > > > local_ghosted_master_id
typename AgglomerationAccessor< dim, spacedim >::AgglomerationContainer AgglomerationContainer
const unsigned int extraction_level
std::map< std::pair< types::global_cell_index, types::global_cell_index >, std::vector< types::global_cell_index > > parent_node_to_children_nodes
std::map< std::pair< types::global_cell_index, unsigned int >, std::pair< bool, typename Triangulation< dim, spacedim >::active_cell_iterator > > cell_face_at_boundary
std::map< std::pair< CellId, unsigned int >, CellId > ghosted_master_id
std::set< std::pair< types::global_cell_index, unsigned int > > visited_cell_and_faces
std::set< std::pair< CellId, unsigned int > > visited_cell_and_faces_id
std::map< std::pair< CellId, CellId >, std::vector< std::pair< typename Triangulation< dim, spacedim >::active_cell_iterator, unsigned int > > > interface
#define Assert(cond, exc)
const types::subdomain_id invalid_subdomain_id
unsigned int global_dof_index
unsigned int subdomain_id
unsigned int global_cell_index