DUNE-ACFEM (2.5.1)

subentitymapper.hh
Go to the documentation of this file.
1
5#ifndef __DUNE_ACFEM_SUBENTITYMAPPER_HH__
6#define __DUNE_ACFEM_SUBENTITYMAPPER_HH__
7
8#include <dune/geometry/referenceelements.hh>
9#include <dune/grid/common/intersection.hh>
10#include <array>
11#include <unordered_map>
12
13namespace Dune {
14
15 namespace ACFem {
16
21 namespace {
22
29 template<class DofMapper, int coDimension>
30 struct MapEntityHierarchyHelper
31 {
32 typedef DofMapper DofMapperType;
33 typedef typename DofMapperType::GridPartType GridPartType;
34 typedef typename GridPartType::ctype FieldType;
35 typedef typename DofMapperType::ElementType ElementType;
36 typedef typename DofMapperType::GlobalKeyType GlobalKeyType;
37 typedef ReferenceElements<FieldType, ElementType::dimension> ReferenceElementsType;
38 typedef ReferenceElement<FieldType, ElementType::dimension> ReferenceElementType;
39
40 // Recurse into the entity hierarchy.
41 template<int codim>
42 struct SubSubEntityFunctor
43 {
44 enum {
45 subCoDimension = codim,
46 gridDimension = GridPartType::dimension,
47 subSubEntityDimension = gridDimension-subCoDimension
48 };
49 template<class RefElem, class Functor>
50 static void apply(const DofMapperType& mapper, const ElementType& element, const RefElem& refElem, int subIndex, Functor f,
51 const std::unordered_map<GlobalKeyType, unsigned>& globalToLocal)
52 {
53 const int numSubSubEntities = refElem.size(subIndex, coDimension, subCoDimension);
54
55 for (int i = 0; i < numSubSubEntities; ++i) {
56 const unsigned subSubEntityIndex = refElem.subEntity(subIndex, coDimension, i, subCoDimension);
57 const auto subSubEntity = element.template subEntity<subCoDimension>(subSubEntityIndex);
58
59 mapper.mapEachEntityDof(subSubEntity, [f, globalToLocal] (unsigned localIndex, const GlobalKeyType& globalIndex) {
60 f(globalToLocal.at(globalIndex), globalIndex);
61 });
62 }
63 }
64 };
65
66 };
67
68 }
69
90 template<class DofMapper, int coDimension, class Functor>
91 void mapEntityHierarchy(const DofMapper& mapper,
92 const typename DofMapper::ElementType& element,
93 int subIndex,
94 const std::integral_constant<int, coDimension>& dummy,
95 Functor f)
96 {
97 typedef DofMapper DofMapperType;
98 typedef MapEntityHierarchyHelper<DofMapperType, coDimension> HelperType;
99 typedef typename HelperType::GridPartType GridPartType;
100 typedef typename HelperType::FieldType FieldType;
101 typedef typename HelperType::ElementType ElementType;
102 typedef typename HelperType::GlobalKeyType GlobalKeyType;
103 typedef typename HelperType::ReferenceElementType ReferenceElementType;
104 typedef typename HelperType::ReferenceElementsType ReferenceElementsType;
105
106 enum { dimension = GridPartType::dimension };
107
108 // first construct a reverse lookup for the global keys
109 std::unordered_map<GlobalKeyType, unsigned> globalToLocal;
110 mapper.mapEach(element, [&globalToLocal] (unsigned local, const GlobalKeyType& global) {
111 globalToLocal[global] = local;
112 });
113
114 // get codim-0 reference element
115 const auto& refElem = ReferenceElementsType::general(element.type());
116
117 // Use ForLoop in order to generate a constant co-dimension as
118 // we need element.template subEntity<cdim>
119#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 5)
120 Hybrid::forEach(Std::make_index_sequence<dimension-coDimension+1>{},
121 [&](auto i){
122 HelperType::template SubSubEntityFunctor<i+coDimension>::apply(mapper, element, refElem, subIndex, f, globalToLocal);
123 });
124#else
125 ForLoop<HelperType::template SubSubEntityFunctor, coDimension, dimension>::apply(mapper, element, refElem, subIndex, f, globalToLocal);
126#endif
127 }
128
148 template<class DofMapper, class Intersection, class Functor>
149 void mapEachFaceDof(const DofMapper& mapper,
150 const Intersection& intersection,
151 Functor f)
152 {
153 mapEntityHierarchy(mapper, intersection.inside(), intersection.indexInInside(),
154 std::integral_constant<int, 1>(), f);
155 }
156
178 template<class DofMapper, class Intersection, class LocalIndices, class GlobalIndices>
179 unsigned mapFaceDofs(const DofMapper& mapper, const Intersection& intersection,
180 LocalIndices& localIndices, GlobalIndices& globalIndices)
181 {
182 typedef typename DofMapper::GlobalKeyType GlobalKeyType;
183
184 unsigned count = 0;
185 mapEachFaceDof(mapper, intersection,
186 [&localIndices, &globalIndices, &count] (unsigned localIndex, GlobalKeyType globalIndex) {
187 localIndices[count] = localIndex;
188 globalIndices[count] = globalIndex;
189 ++count;
190 });
191
192 return count;
193 }
194
197 }
198
199}
200
201#endif // __DUNE_ACFEM_SUBENTITYMAPPER_HH__
void mapEachFaceDof(const DofMapper &mapper, const Intersection &intersection, Functor f)
Iterate over the codim-1 sub-entity hierarchy linked to the given intersection and call the dof-mappe...
Definition: subentitymapper.hh:149
void mapEntityHierarchy(const DofMapper &mapper, const typename DofMapper::ElementType &element, int subIndex, const std::integral_constant< int, coDimension > &dummy, Functor f)
Iterate over the entire sub-entity hierarchy below one given sub-entity (including the given one) and...
Definition: subentitymapper.hh:91
unsigned mapFaceDofs(const DofMapper &mapper, const Intersection &intersection, LocalIndices &localIndices, GlobalIndices &globalIndices)
Fetch all global DoFs of the codim-1 entity intersection belongs to.
Definition: subentitymapper.hh:179
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Mar 12, 23:28, 2025)