PolyDEAL
 
Loading...
Searching...
No Matches
agglomeration_handler.h
Go to the documentation of this file.
1// -----------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
4// Copyright (C) XXXX - YYYY by the polyDEAL authors
5//
6// This file is part of the polyDEAL library.
7//
8// Detailed license information governing the source code
9// can be found in LICENSE.md at the top level directory.
10//
11// -----------------------------------------------------------------------------
12#ifndef agglomeration_handler_h
13#define agglomeration_handler_h
14
15#include <deal.II/base/mpi.h>
19
22
25
26#include <deal.II/fe/fe_dgp.h>
27#include <deal.II/fe/fe_dgq.h>
34
36#include <deal.II/grid/tria.h>
37
39
45#include <deal.II/lac/vector.h>
46
48
51
53#include <agglomerator.h>
54#include <mapping_box.h>
55#include <utils.h>
56
57#include <fstream>
58
59using namespace dealii;
60
61// Forward declarations
62template <int dim, int spacedim>
64
65namespace dealii
66{
67 namespace internal
68 {
72 template <int, int>
74 } // namespace internal
75} // namespace dealii
76
77
78
83namespace dealii
84{
85 namespace internal
86 {
87 template <int dim, int spacedim>
89 {
90 public:
94 PolytopeCache() = default;
95
99 ~PolytopeCache() = default;
100
101 void
103 {
104 // clear all the members
105 cell_face_at_boundary.clear();
106 interface.clear();
108 }
109
117 mutable std::set<std::pair<types::global_cell_index, unsigned int>>
119
120
121 mutable std::set<std::pair<CellId, unsigned int>>
123
124
125
134 mutable std::map<
135 std::pair<types::global_cell_index, unsigned int>,
136 std::pair<bool,
139
144 mutable std::map<std::pair<CellId, unsigned int>, CellId>
146
157 mutable std::map<
158 std::pair<CellId, CellId>,
159 std::vector<
160 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
161 unsigned int>>>
163 };
164 } // namespace internal
165} // namespace dealii
166
167
171template <int dim, int spacedim = dim>
173{
174public:
176
179
180
182 {
185 };
186
187
188
189 explicit AgglomerationHandler(
191
193
195 {
196 // disconnect the signal
197 tria_listener.disconnect();
198 }
199
204 begin() const;
205
211
216 end() const;
217
223
229
236
237 template <int, int>
239
240 template <int, int>
242
247 void
249
255 void
257 const Quadrature<dim> &cell_quadrature = QGauss<dim>(1),
258 const UpdateFlags &flags = UpdateFlags::update_default,
259 const Quadrature<dim - 1> &face_quadrature = QGauss<dim - 1>(1),
260 const UpdateFlags &face_flags = UpdateFlags::update_default);
261
268 template <typename SparsityPatternType, typename Number = double>
269 void
271 SparsityPatternType &sparsity_pattern,
273 const bool keep_constrained_dofs = true,
275
286
292 void
294
295 inline const Triangulation<dim, spacedim> &
297
298 inline const FiniteElement<dim, spacedim> &
299 get_fe() const;
300
301 inline const Mapping<dim> &
302 get_mapping() const;
303
304 inline const MappingBox<dim> &
306
307 inline const std::vector<BoundingBox<dim>> &
309
314 double
316
320 const;
321
322 inline decltype(auto)
324
328 template <typename CellIterator>
329 inline bool
330 is_master_cell(const CellIterator &cell) const;
331
336 inline const std::vector<
339
340
343
351 inline std::vector<
355 &master_cell) const;
356
361 template <class StreamType>
362 void
363 print_agglomeration(StreamType &out)
364 {
365 for (const auto &cell : tria->active_cell_iterators())
366 out << "Cell with index: " << cell->active_cell_index()
367 << " has associated value: "
368 << master_slave_relationships[cell->global_active_cell_index()]
369 << std::endl;
370 }
371
378 inline const DoFHandler<dim, spacedim> &
380
384 unsigned int
386
390 unsigned int
393 const;
394
399 reinit(const AgglomerationIterator<dim, spacedim> &polytope) const;
400
407 const unsigned int face_index) const;
408
414 std::pair<const FEValuesBase<dim, spacedim> &,
417 const AgglomerationIterator<dim, spacedim> &neigh_polytope,
418 const unsigned int local_in,
419 const unsigned int local_outside) const;
420
428 const AgglomerationContainer &cells,
430 &master_cell) const;
431
432
442 inline bool
445 const unsigned int f) const;
446
447 inline unsigned int
448 n_dofs_per_cell() const noexcept;
449
450 inline types::global_dof_index
451 n_dofs() const noexcept;
452
453
454
460 inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
462 const typename Triangulation<dim>::active_cell_iterator &cell);
463
464
468 DoFHandler<dim, spacedim> agglo_dh;
469
473 DoFHandler<dim, spacedim> output_dh;
474
475 std::unique_ptr<MappingBox<dim>> box_mapping;
476
485 void
487
488 void
490
491 // TODO: move it to private interface
492 mutable std::map<
493 types::subdomain_id,
494 std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>>
496
497 mutable std::map<
498 types::subdomain_id,
499 std::map<std::pair<CellId, unsigned int>, std::vector<double>>>
501
502 mutable std::map<
503 types::subdomain_id,
504 std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>>
506
507 mutable std::map<
508 types::subdomain_id,
509 std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>>
511
512 mutable std::map<types::subdomain_id,
513 std::map<std::pair<CellId, unsigned int>,
514 std::vector<std::vector<Tensor<1, spacedim>>>>>
516
521 inline const typename DoFHandler<dim, spacedim>::active_cell_iterator
522 polytope_to_dh_iterator(const types::global_cell_index polytope_index) const;
523
527 template <typename RtreeType>
528 void
529 connect_hierarchy(const CellsAgglomerator<dim, RtreeType> &agglomerator);
530
531private:
535 void
537 const std::unique_ptr<GridTools::Cache<dim, spacedim>> &cache_tria);
538
539 void
541 AgglomerationContainer &polytope,
542 const typename Triangulation<dim, spacedim>::active_cell_iterator
543 &master_cell);
544
548 void
550 {
551 // First disconnect existing connections
552 tria_listener.disconnect();
553 tria_listener = tria->signals.any_change.connect(
554 [&]() { this->initialize_agglomeration_data(this->cached_tria); });
555 }
556
562
566
571 void
573
574
578 const;
579
583 inline bool
587 &other_cell) const;
588
597 void
599
600
607 const unsigned int face_number,
609 &agglo_isv_ptr) const;
610
611
616 template <typename CellIterator>
617 inline bool
618 is_slave_cell(const CellIterator &cell) const;
619
620
626 void
628
629
633 unsigned int n_agglomerations;
634
635
643
651
653
654 mutable std::vector<types::global_cell_index> number_of_agglomerated_faces;
655
661 mutable std::map<
663 std::vector<typename Triangulation<dim>::active_face_iterator>>
665
666
673 std::vector<BoundingBox<spacedim>> bboxes;
674
676
677
678 // n_faces
679 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
681
682 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
684
685
686 // CellId (including slaves)
687 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
689
690 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
692
693
694 // send to neighborign rank the information that
695 // - current polytope id
696 // - face f
697 // has the following neighboring id.
698 mutable std::map<types::subdomain_id,
699 std::map<CellId, std::map<unsigned int, CellId>>>
701
702 mutable std::map<types::subdomain_id,
703 std::map<CellId, std::map<unsigned int, CellId>>>
705
706 // CellIds from neighboring rank
707 mutable std::map<types::subdomain_id,
708 std::map<CellId, std::map<unsigned int, bool>>>
710
711 mutable std::map<types::subdomain_id,
712 std::map<CellId, std::map<unsigned int, bool>>>
714
715 // Exchange neighboring bounding boxes
716 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
718
719 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
721
722 // Exchange DoF indices with ghosted polytopes
723 mutable std::map<types::subdomain_id,
724 std::map<CellId, std::vector<types::global_dof_index>>>
726
727 mutable std::map<types::subdomain_id,
728 std::map<CellId, std::vector<types::global_dof_index>>>
730
731 // Exchange qpoints
732 mutable std::map<
734 std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>>
736
737 // Exchange jxws
738 mutable std::map<
740 std::map<std::pair<CellId, unsigned int>, std::vector<double>>>
742
743 // Exchange normals
744 mutable std::map<
746 std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>>
748
749 // Exchange values
750 mutable std::map<
752 std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>>
754
755 mutable std::map<types::subdomain_id,
756 std::map<std::pair<CellId, unsigned int>,
757 std::vector<std::vector<Tensor<1, spacedim>>>>>
759
760
761
763
765
767
768 std::unique_ptr<GridTools::Cache<dim, spacedim>> cached_tria;
769
770 const MPI_Comm communicator;
771
772 // The FiniteElement space we have on each cell. Currently supported types are
773 // FE_DGQ and FE_DGP elements.
774 std::unique_ptr<FiniteElement<dim>> fe;
775
777
782
783
788 mutable std::unique_ptr<ScratchData> standard_scratch;
789
795 mutable std::unique_ptr<ScratchData> agglomerated_scratch;
796
797
798 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
800
801 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
803
804 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
806
807 boost::signals2::connection tria_listener;
808
810
814
816
820
822
824
825 // Associate the master cell to the slaves.
826 std::unordered_map<
828 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
830
831 // Map the master cell index with the polytope index
832 std::map<types::global_cell_index, types::global_cell_index> master2polygon;
833
834
835 std::vector<typename Triangulation<dim>::active_cell_iterator>
837
838 // Dummy FiniteElement objects needed only to generate quadratures
839
844
848 std::unique_ptr<FEValues<dim, spacedim>> no_values;
849
853 std::unique_ptr<FEFaceValues<dim, spacedim>> no_face_values;
854
858 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
860
861 friend class internal::AgglomerationHandlerImplementation<dim, spacedim>;
862
864
870
871 std::map<std::pair<types::global_cell_index, types::global_cell_index>,
872 std::vector<types::global_cell_index>>
874
876};
877
878
879
880// ------------------------------ inline functions -------------------------
881template <int dim, int spacedim>
882inline const FiniteElement<dim, spacedim> &
884{
885 return *fe;
886}
887
888
889
890template <int dim, int spacedim>
891inline const Mapping<dim> &
896
897
898
899template <int dim, int spacedim>
900inline const MappingBox<dim> &
905
906
907
908template <int dim, int spacedim>
909inline const Triangulation<dim, spacedim> &
914
915
916template <int dim, int spacedim>
917inline const std::vector<BoundingBox<dim>> &
922
923
924
925template <int dim, int spacedim>
929{
930 return master2polygon.at(cell->active_cell_index());
931}
932
933
934
935template <int dim, int spacedim>
936inline decltype(auto)
941
942
943
944template <int dim, int spacedim>
950
951
952
953template <int dim, int spacedim>
954inline std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
957 &master_cell) const
958{
959 Assert(is_master_cell(master_cell), ExcInternalError());
960 auto agglomeration = get_slaves_of_idx(master_cell->active_cell_index());
961 agglomeration.push_back(master_cell);
962 return agglomeration;
963}
964
965
966
967template <int dim, int spacedim>
968inline const DoFHandler<dim, spacedim> &
973
974
975
976template <int dim, int spacedim>
977inline const std::vector<
984
985
986
987template <int dim, int spacedim>
988template <typename CellIterator>
989inline bool
991 const CellIterator &cell) const
992{
993 return master_slave_relationships[cell->global_active_cell_index()] == -1;
994}
995
996
997
1002template <int dim, int spacedim>
1003template <typename CellIterator>
1004inline bool
1006 const CellIterator &cell) const
1007{
1008 return master_slave_relationships[cell->global_active_cell_index()] >= 0;
1009}
1010
1011
1012
1013template <int dim, int spacedim>
1014inline bool
1017 const unsigned int face_index) const
1018{
1019 Assert(!is_slave_cell(cell),
1020 ExcMessage("This function should not be called for a slave cell."));
1021
1022 return polytope_cache.cell_face_at_boundary
1023 .at({master2polygon.at(cell->active_cell_index()), face_index})
1024 .first;
1025}
1026
1027
1028template <int dim, int spacedim>
1029inline unsigned int
1031{
1032 return fe->n_dofs_per_cell();
1033}
1034
1035
1036
1037template <int dim, int spacedim>
1040{
1041 return agglo_dh.n_dofs();
1042}
1043
1044
1045
1046template <int dim, int spacedim>
1047inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
1053
1054
1055
1056template <int dim, int spacedim>
1063
1064
1065
1066template <int dim, int spacedim>
1069 const typename Triangulation<dim, spacedim>::active_cell_iterator &cell) const
1070{
1071 auto idx = master_slave_relationships[cell->global_active_cell_index()];
1072 if (idx == -1)
1073 return cell->global_active_cell_index();
1074 else
1075 return static_cast<types::global_cell_index>(idx);
1076}
1077
1078
1079
1080template <int dim, int spacedim>
1081inline bool
1084 const typename Triangulation<dim, spacedim>::active_cell_iterator &other_cell)
1085 const
1086{
1087 // if different subdomain, then **by construction** they will not be together
1088 // if (cell->subdomain_id() != other_cell->subdomain_id())
1089 // return false;
1090 // else
1091 return (get_master_idx_of_cell(cell) == get_master_idx_of_cell(other_cell));
1092}
1093
1094
1095
1096template <int dim, int spacedim>
1097inline unsigned int
1102
1103
1104
1105template <int dim, int spacedim>
1108 const types::global_cell_index polytope_index) const
1109{
1110 return master_cells_container[polytope_index]->as_dof_handler_iterator(
1111 agglo_dh);
1112}
1113
1114
1115
1116template <int dim, int spacedim>
1119{
1121 ExcMessage("No agglomeration has been performed."));
1122 return {*master_cells_container.begin(), this};
1123}
1124
1125
1126
1127template <int dim, int spacedim>
1130{
1132 ExcMessage("No agglomeration has been performed."));
1133 return {*master_cells_container.begin(), this};
1134}
1135
1136
1137
1138template <int dim, int spacedim>
1141{
1143 ExcMessage("No agglomeration has been performed."));
1144 return {*master_cells_container.end(), this};
1145}
1146
1147
1148
1149template <int dim, int spacedim>
1152{
1154 ExcMessage("No agglomeration has been performed."));
1155 return {*master_cells_container.end(), this};
1156}
1157
1158
1159
1160template <int dim, int spacedim>
1163{
1165 ExcMessage("No agglomeration has been performed."));
1166 return {master_cells_container.back(), this};
1167}
1168
1169
1170
1171template <int dim, int spacedim>
1180
1181template <int dim, int spacedim>
1182template <typename RtreeType>
1183void
1190
1191
1192#endif
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
MeshWorker::ScratchData< dim, spacedim > ScratchData
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
SmartPointer< const Triangulation< dim, spacedim > > tria
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
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
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
AgglomerationHandler(const GridTools::Cache< dim, spacedim > &cached_tria)
std::map< types::subdomain_id, std::map< CellId, std::vector< types::global_dof_index > > > recv_ghost_dofs
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)
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
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
hp::FECollection< dim, spacedim > fe_collection
std::map< types::subdomain_id, std::map< CellId, BoundingBox< dim > > > local_ghosted_bbox
DoFHandler< dim, spacedim > agglo_dh
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
void update_agglomerate(AgglomerationContainer &polytope, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell)
double get_mesh_size() const
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, CellId > > > recv_ghosted_master_id
SmartPointer< const Mapping< dim, spacedim > > mapping
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
EnableObserverPointer Subscriptor
Point< 2 > first
#define Assert(cond, exc)
UpdateFlags
update_values
update_normal_vectors
update_JxW_values
update_inverse_jacobians
update_gradients
update_quadrature_points
constexpr types::subdomain_id invalid_subdomain_id
unsigned int global_dof_index
unsigned int subdomain_id
unsigned int global_cell_index
ObserverPointer< T, P > SmartPointer