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>
18
21
24
25#include <deal.II/fe/fe_dgp.h>
26#include <deal.II/fe/fe_dgq.h>
33
35#include <deal.II/grid/tria.h>
36
38
44#include <deal.II/lac/vector.h>
45
47
50
52#include <agglomerator.h>
53#include <mapping_box.h>
54#include <utils.h>
55
56#include <fstream>
57
58using namespace dealii;
59
60// Forward declarations
61template <int dim, int spacedim>
63
64namespace dealii
65{
66 namespace internal
67 {
71 template <int, int>
73 } // namespace internal
74} // namespace dealii
75
76
77
82namespace dealii
83{
84 namespace internal
85 {
86 template <int dim, int spacedim>
88 {
89 public:
93 PolytopeCache() = default;
94
98 ~PolytopeCache() = default;
99
100 void
102 {
103 // clear all the members
104 cell_face_at_boundary.clear();
105 interface.clear();
107 }
108
116 mutable std::set<std::pair<types::global_cell_index, unsigned int>>
118
119
120 mutable std::set<std::pair<CellId, unsigned int>>
122
123
124
133 mutable std::map<
134 std::pair<types::global_cell_index, unsigned int>,
135 std::pair<bool,
138
143 mutable std::map<std::pair<CellId, unsigned int>, CellId>
145
156 mutable std::map<
157 std::pair<CellId, CellId>,
158 std::vector<
159 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
160 unsigned int>>>
162 };
163 } // namespace internal
164} // namespace dealii
165
166
170template <int dim, int spacedim = dim>
172{
173public:
175
178
179
181 {
183 slave = 1
184 };
185
186
187
188 explicit AgglomerationHandler(
190
192
194 {
195 // disconnect the signal
196 tria_listener.disconnect();
197 }
198
203 begin() const;
204
210
215 end() const;
216
222
228
235
236 template <int, int>
238
239 template <int, int>
241
246 void
248
254 void
256 const Quadrature<dim> &cell_quadrature = QGauss<dim>(1),
257 const UpdateFlags &flags = UpdateFlags::update_default,
258 const Quadrature<dim - 1> &face_quadrature = QGauss<dim - 1>(1),
259 const UpdateFlags &face_flags = UpdateFlags::update_default);
260
267 template <typename SparsityPatternType, typename Number = double>
268 void
270 SparsityPatternType &sparsity_pattern,
272 const bool keep_constrained_dofs = true,
274
285
291 void
293
294 inline const Triangulation<dim, spacedim> &
296
297 inline const FiniteElement<dim, spacedim> &
298 get_fe() const;
299
300 inline const Mapping<dim> &
301 get_mapping() const;
302
303 inline const MappingBox<dim> &
305
306 inline const std::vector<BoundingBox<dim>> &
308
313 double
315
319 const;
320
321 inline decltype(auto)
323
327 template <typename CellIterator>
328 inline bool
329 is_master_cell(const CellIterator &cell) const;
330
335 inline const std::vector<
338
339
342
350 inline std::vector<
354 &master_cell) const;
355
360 template <class StreamType>
361 void
362 print_agglomeration(StreamType &out)
363 {
364 for (const auto &cell : tria->active_cell_iterators())
365 out << "Cell with index: " << cell->active_cell_index()
366 << " has associated value: "
367 << master_slave_relationships[cell->global_active_cell_index()]
368 << std::endl;
369 }
370
377 inline const DoFHandler<dim, spacedim> &
379
383 unsigned int
385
389 unsigned int
392 const;
393
398 reinit(const AgglomerationIterator<dim, spacedim> &polytope) const;
399
406 const unsigned int face_index) const;
407
413 std::pair<const FEValuesBase<dim, spacedim> &,
416 const AgglomerationIterator<dim, spacedim> &neigh_polytope,
417 const unsigned int local_in,
418 const unsigned int local_outside) const;
419
427 const AgglomerationContainer &cells,
429 &master_cell) const;
430
431
441 inline bool
444 const unsigned int f) const;
445
446 inline unsigned int
447 n_dofs_per_cell() const noexcept;
448
449 inline types::global_dof_index
450 n_dofs() const noexcept;
451
452
453
459 inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
461 const typename Triangulation<dim>::active_cell_iterator &cell);
462
463
467 DoFHandler<dim, spacedim> agglo_dh;
468
472 DoFHandler<dim, spacedim> output_dh;
473
474 std::unique_ptr<MappingBox<dim>> box_mapping;
475
484 void
486
487 void
489
490 // TODO: move it to private interface
491 mutable std::map<
492 types::subdomain_id,
493 std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>>
495
496 mutable std::map<
497 types::subdomain_id,
498 std::map<std::pair<CellId, unsigned int>, std::vector<double>>>
500
501 mutable std::map<
502 types::subdomain_id,
503 std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>>
505
506 mutable std::map<
507 types::subdomain_id,
508 std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>>
510
511 mutable std::map<types::subdomain_id,
512 std::map<std::pair<CellId, unsigned int>,
513 std::vector<std::vector<Tensor<1, spacedim>>>>>
515
520 inline const typename DoFHandler<dim, spacedim>::active_cell_iterator
521 polytope_to_dh_iterator(const types::global_cell_index polytope_index) const;
522
526 template <typename RtreeType>
527 void
528 connect_hierarchy(const CellsAgglomerator<dim, RtreeType> &agglomerator);
529
530private:
534 void
536 const std::unique_ptr<GridTools::Cache<dim, spacedim>> &cache_tria);
537
538 void
540 AgglomerationContainer &polytope,
541 const typename Triangulation<dim, spacedim>::active_cell_iterator
542 &master_cell);
543
547 void
549 {
550 // First disconnect existing connections
551 tria_listener.disconnect();
552 tria_listener = tria->signals.any_change.connect(
553 [&]() { this->initialize_agglomeration_data(this->cached_tria); });
554 }
555
565
570 void
572
573
577 const;
578
582 inline bool
586 &other_cell) const;
587
596 void
598
599
606 const unsigned int face_number,
608 &agglo_isv_ptr) const;
609
610
615 template <typename CellIterator>
616 inline bool
617 is_slave_cell(const CellIterator &cell) const;
618
619
625 void
627
628
632 unsigned int n_agglomerations;
633
634
642
650
652
653 mutable std::vector<types::global_cell_index> number_of_agglomerated_faces;
654
660 mutable std::map<
662 std::vector<typename Triangulation<dim>::active_face_iterator>>
664
665
672 std::vector<BoundingBox<spacedim>> bboxes;
673
675
676
677 // n_faces
678 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
680
681 mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>>
683
684
685 // CellId (including slaves)
686 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
688
689 mutable std::map<types::subdomain_id, std::map<CellId, CellId>>
691
692
693 // send to neighborign rank the information that
694 // - current polytope id
695 // - face f
696 // has the following neighboring id.
697 mutable std::map<types::subdomain_id,
698 std::map<CellId, std::map<unsigned int, CellId>>>
700
701 mutable std::map<types::subdomain_id,
702 std::map<CellId, std::map<unsigned int, CellId>>>
704
705 // CellIds from neighboring rank
706 mutable std::map<types::subdomain_id,
707 std::map<CellId, std::map<unsigned int, bool>>>
709
710 mutable std::map<types::subdomain_id,
711 std::map<CellId, std::map<unsigned int, bool>>>
713
714 // Exchange neighboring bounding boxes
715 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
717
718 mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>>
720
721 // Exchange DoF indices with ghosted polytopes
722 mutable std::map<types::subdomain_id,
723 std::map<CellId, std::vector<types::global_dof_index>>>
725
726 mutable std::map<types::subdomain_id,
727 std::map<CellId, std::vector<types::global_dof_index>>>
729
730 // Exchange qpoints
731 mutable std::map<
733 std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>>
735
736 // Exchange jxws
737 mutable std::map<
739 std::map<std::pair<CellId, unsigned int>, std::vector<double>>>
741
742 // Exchange normals
743 mutable std::map<
745 std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>>
747
748 // Exchange values
749 mutable std::map<
751 std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>>
753
754 mutable std::map<types::subdomain_id,
755 std::map<std::pair<CellId, unsigned int>,
756 std::vector<std::vector<Tensor<1, spacedim>>>>>
758
759
760
762
764
766
767 std::unique_ptr<GridTools::Cache<dim, spacedim>> cached_tria;
768
769 const MPI_Comm communicator;
770
771 // The FiniteElement space we have on each cell. Currently supported types are
772 // FE_DGQ and FE_DGP elements.
773 std::unique_ptr<FiniteElement<dim>> fe;
774
776
781
782
787 mutable std::unique_ptr<ScratchData> standard_scratch;
788
794 mutable std::unique_ptr<ScratchData> agglomerated_scratch;
795
796
797 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
799
800 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
802
803 mutable std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim>>
805
806 boost::signals2::connection tria_listener;
807
809
813
815
819
821
823
824 // Associate the master cell to the slaves.
825 std::unordered_map<
827 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
829
830 // Map the master cell index with the polytope index
831 std::map<types::global_cell_index, types::global_cell_index> master2polygon;
832
833
834 std::vector<typename Triangulation<dim>::active_cell_iterator>
836
837 // Dummy FiniteElement objects needed only to generate quadratures
838
843
847 std::unique_ptr<FEValues<dim, spacedim>> no_values;
848
852 std::unique_ptr<FEFaceValues<dim, spacedim>> no_face_values;
853
857 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
859
860 friend class internal::AgglomerationHandlerImplementation<dim, spacedim>;
861
863
869
870 std::map<std::pair<types::global_cell_index, types::global_cell_index>,
871 std::vector<types::global_cell_index>>
873
875};
876
877
878
879// ------------------------------ inline functions -------------------------
880template <int dim, int spacedim>
881inline const FiniteElement<dim, spacedim> &
883{
884 return *fe;
885}
886
887
888
889template <int dim, int spacedim>
890inline const Mapping<dim> &
892{
893 return *mapping;
894}
895
896
897
898template <int dim, int spacedim>
899inline const MappingBox<dim> &
901{
902 return *box_mapping;
903}
904
905
906
907template <int dim, int spacedim>
908inline const Triangulation<dim, spacedim> &
913
914
915template <int dim, int spacedim>
916inline const std::vector<BoundingBox<dim>> &
918{
919 return bboxes;
920}
921
922
923
924template <int dim, int spacedim>
928{
929 return master2polygon.at(cell->active_cell_index());
930}
931
932
933
934template <int dim, int spacedim>
935inline decltype(auto)
937{
938 return polytope_cache.interface;
939}
940
941
942
943template <int dim, int spacedim>
946{
947 return master_slave_relationships;
948}
949
950
951
952template <int dim, int spacedim>
953inline std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
956 &master_cell) const
957{
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;
962}
963
964
965
966template <int dim, int spacedim>
967inline const DoFHandler<dim, spacedim> &
969{
970 return agglo_dh;
971}
972
973
974
975template <int dim, int spacedim>
976inline const std::vector<
979 types::global_cell_index idx) const
980{
981 return master2slaves.at(idx);
982}
983
984
985
986template <int dim, int spacedim>
987template <typename CellIterator>
988inline bool
990 const CellIterator &cell) const
991{
992 return master_slave_relationships[cell->global_active_cell_index()] == -1;
993}
994
995
996
1001template <int dim, int spacedim>
1002template <typename CellIterator>
1003inline bool
1005 const CellIterator &cell) const
1006{
1007 return master_slave_relationships[cell->global_active_cell_index()] >= 0;
1008}
1009
1010
1011
1012template <int dim, int spacedim>
1013inline bool
1016 const unsigned int face_index) const
1017{
1018 Assert(!is_slave_cell(cell),
1019 ExcMessage("This function should not be called for a slave cell."));
1020
1021 return polytope_cache.cell_face_at_boundary
1022 .at({master2polygon.at(cell->active_cell_index()), face_index})
1023 .first;
1024}
1025
1026
1027template <int dim, int spacedim>
1028inline unsigned int
1030{
1031 return fe->n_dofs_per_cell();
1032}
1033
1034
1035
1036template <int dim, int spacedim>
1039{
1040 return agglo_dh.n_dofs();
1041}
1042
1043
1044
1045template <int dim, int spacedim>
1046inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
1048 const typename Triangulation<dim>::active_cell_iterator &cell)
1049{
1050 return polygon_boundary[cell];
1051}
1052
1053
1054
1055template <int dim, int spacedim>
1059{
1060 return master_slave_relationships_iterators.at(cell->active_cell_index());
1061}
1062
1063
1064
1065template <int dim, int spacedim>
1068 const typename Triangulation<dim, spacedim>::active_cell_iterator &cell) const
1069{
1070 auto idx = master_slave_relationships[cell->global_active_cell_index()];
1071 if (idx == -1)
1072 return cell->global_active_cell_index();
1073 else
1074 return static_cast<types::global_cell_index>(idx);
1075}
1076
1077
1078
1079template <int dim, int spacedim>
1080inline bool
1083 const typename Triangulation<dim, spacedim>::active_cell_iterator &other_cell)
1084 const
1085{
1086 // if different subdomain, then **by construction** they will not be together
1087 // if (cell->subdomain_id() != other_cell->subdomain_id())
1088 // return false;
1089 // else
1090 return (get_master_idx_of_cell(cell) == get_master_idx_of_cell(other_cell));
1091}
1092
1093
1094
1095template <int dim, int spacedim>
1096inline unsigned int
1098{
1099 return n_agglomerations;
1100}
1101
1102
1103
1104template <int dim, int spacedim>
1107 const types::global_cell_index polytope_index) const
1108{
1109 return master_cells_container[polytope_index]->as_dof_handler_iterator(
1110 agglo_dh);
1111}
1112
1113
1114
1115template <int dim, int spacedim>
1118{
1119 Assert(n_agglomerations > 0,
1120 ExcMessage("No agglomeration has been performed."));
1121 return {*master_cells_container.begin(), this};
1122}
1123
1124
1125
1126template <int dim, int spacedim>
1129{
1130 Assert(n_agglomerations > 0,
1131 ExcMessage("No agglomeration has been performed."));
1132 return {*master_cells_container.begin(), this};
1133}
1134
1135
1136
1137template <int dim, int spacedim>
1140{
1141 Assert(n_agglomerations > 0,
1142 ExcMessage("No agglomeration has been performed."));
1143 return {*master_cells_container.end(), this};
1144}
1145
1146
1147
1148template <int dim, int spacedim>
1151{
1152 Assert(n_agglomerations > 0,
1153 ExcMessage("No agglomeration has been performed."));
1154 return {*master_cells_container.end(), this};
1155}
1156
1157
1158
1159template <int dim, int spacedim>
1162{
1163 Assert(n_agglomerations > 0,
1164 ExcMessage("No agglomeration has been performed."));
1165 return {master_cells_container.back(), this};
1166}
1167
1168
1169
1170template <int dim, int spacedim>
1179
1180template <int dim, int spacedim>
1181template <typename RtreeType>
1182void
1184 const CellsAgglomerator<dim, RtreeType> &agglomerator)
1185{
1186 parent_child_info = agglomerator.parent_node_to_children_nodes;
1187 present_extraction_level = agglomerator.extraction_level;
1188}
1189
1190
1191#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
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
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
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
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
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
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
Point< 2 > first
#define Assert(cond, exc)
UpdateFlags
update_values
update_normal_vectors
update_JxW_values
update_inverse_jacobians
update_gradients
update_quadrature_points
const types::subdomain_id invalid_subdomain_id
unsigned int global_dof_index
unsigned int subdomain_id
unsigned int global_cell_index