version 3.8
enrichmenthelper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
7
13#ifndef DUMUX_VERTEX_ENRICHMENT_HELPER_HH
14#define DUMUX_VERTEX_ENRICHMENT_HELPER_HH
15
16#include <cmath>
17#include <vector>
18#include <unordered_map>
19#include <unordered_set>
20
21#include <dune/common/exceptions.hh>
22#include <dune/common/reservedvector.hh>
23
24#include <dune/grid/common/mcmgmapper.hh>
25
26#include <dumux/common/math.hh>
28
29namespace Dumux {
30
42template< class GridView, class CodimOneGridView >
44{
45 static constexpr int dim = GridView::dimension;
46 static constexpr int dimWorld = GridView::dimensionworld;
47 static_assert(dim == 2 || dim == 3, "Grid dimension must be two or three");
48 static_assert(dimWorld == int(CodimOneGridView::dimensionworld), "world dimension mismatch");
49
50 using Intersection = typename GridView::Intersection;
51 using MCMGMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
52 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
53 using Element = typename GridView::template Codim<0>::Entity;
54 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
55
56 using ElementPath = std::vector< GridIndexType >;
57 using NodalElementPaths = std::vector< std::vector<ElementPath> >;
58
59public:
60
76 template< class IndexMap, class CodimOneGridAdapter >
77 static std::size_t enrich(IndexMap& indexMap,
78 const std::vector<bool>& vertexMarkers,
79 const GridView& gridView,
80 const MCMGMapper& vertexMapper,
81 const MCMGMapper& elementMapper,
82 const CodimOneGridView& codimOneGridView,
83 const CodimOneGridAdapter& codimOneGridAdapter)
84 {
85 // determine the bulk element paths around vertices marked for enrichment
86 NodalElementPaths nodalPaths(gridView.size(dim));
87 for (const auto& e : elements(gridView))
88 {
89 const auto eIdx = elementMapper.index(e);
90
91 std::vector<unsigned int> handledFacets;
92 const auto refElement = referenceElement(e);
93 for (const auto& is : intersections(gridView, e))
94 {
95 // do not start path search on boundary intersections
96 if (!is.neighbor())
97 continue;
98
99 // if facet has been handled already, skip rest (necessary for e.g. Dune::FoamGrid)
100 if (std::find(handledFacets.begin(), handledFacets.end(), is.indexInInside()) != handledFacets.end())
101 continue;
102
103 // first, get indices of all vertices of this intersection
104 const auto numCorners = is.geometry().corners();
105 std::vector<GridIndexType> faceVertexIndices(numCorners);
106 for (int i = 0; i < numCorners; ++i)
107 faceVertexIndices[i] = vertexMapper.subIndex( e,
108 refElement.subEntity(is.indexInInside(), 1, i, dim),
109 dim );
110
111 // start path search on intersections that lie on a facet grid element
112 if (codimOneGridAdapter.composeFacetElement(faceVertexIndices))
113 {
114 handledFacets.push_back(is.indexInInside());
115 for (int i = 0; i < numCorners; ++i)
116 {
117 // set up path only around those vertices that are to be enriched
118 const auto vIdxGlobal = faceVertexIndices[i];
119 if (!vertexMarkers[vIdxGlobal])
120 continue;
121
122 // construct the path only if element is not yet contained in another path
123 bool found = false;
124 for (const auto& path : nodalPaths[vIdxGlobal])
125 if (std::find(path.begin(), path.end(), eIdx) != path.end())
126 { found = true; break; }
127
128 if (!found)
129 {
130 ElementPath path;
131
132 // Reserve enough memory so that reallocation during the recursive search
133 // does not happen and references are invalidated. Memory is not an issue
134 // here as after enrichment the container is thrown away again.
135 // TODO improve this!?
136 path.reserve(150);
137 path.push_back(eIdx);
138 continuePathSearch_(path, gridView, elementMapper, vertexMapper, codimOneGridAdapter, e, refElement, is, vIdxGlobal);
139 nodalPaths[vIdxGlobal].emplace_back(std::move(path));
140 }
141 }
142 }
143 }
144 }
145
146 // determine the offsets for each bulk vertex index on the basis of the paths found per vertex
147 std::vector<std::size_t> bulkVertexIndexOffsets(gridView.size(dim), 0);
148 for (const auto& v : vertices(gridView))
149 {
150 const auto vIdx = vertexMapper.index(v);
151 if (vertexMarkers[vIdx])
152 bulkVertexIndexOffsets[vIdx] = nodalPaths[vIdx].size()-1;
153 }
154
155 // ... and accumulate the offsets
156 std::size_t sumOffset = 0;
157 std::size_t size = 0;
158 for (auto& nodalOffset : bulkVertexIndexOffsets)
159 {
160 const auto os = nodalOffset;
161 nodalOffset = sumOffset;
162 sumOffset += os;
163 size += (os == 0) ? 1 : os + 1;
164 }
165
166 // Now, finally set up the new index map
167 for (const auto& e : elements(gridView))
168 {
169 const auto& eg = e.geometry();
170 const auto eIdx = elementMapper.index(e);
171 for (int i = 0; i < eg.corners(); ++i)
172 {
173 const auto origVIdx = vertexMapper.subIndex(e, i, dim);
174
175 // it the node itself is not enriched, simply add offset
176 if (!vertexMarkers[origVIdx])
177 indexMap[eIdx][i] += bulkVertexIndexOffsets[origVIdx];
178
179 // find the local index of the path the element belongs to
180 else
181 {
182 bool found = false;
183 const auto& paths = nodalPaths[ origVIdx ];
184 for (int pathIdx = 0; pathIdx < paths.size(); ++pathIdx)
185 {
186 const auto& curPath = paths[pathIdx];
187 if ( std::find(curPath.begin(), curPath.end(), eIdx) != curPath.end() )
188 {
189 indexMap[eIdx][i] += bulkVertexIndexOffsets[origVIdx] + pathIdx;
190 found = true;
191 break;
192 }
193 }
194
195 if (!found)
196 DUNE_THROW(Dune::InvalidStateException, "Element not found in any path");
197 }
198 }
199 }
200
201 // return the number of dofs after enrichment
202 return size;
203 }
204
205private:
206 template< class ReferenceElement, class CodimOneGridAdapter >
207 static void continuePathSearch_(ElementPath& path,
208 const GridView& gridView,
209 const MCMGMapper& elementMapper,
210 const MCMGMapper& vertexMapper,
211 const CodimOneGridAdapter& codimOneGridAdapter,
212 const Element& element,
213 const ReferenceElement& refElement,
214 const Intersection& prevIntersection,
215 GridIndexType vIdxGlobal)
216 {
217 // on 2d/3d grids, we need to find 1/2 more intersections at vertex
218 static constexpr int numIsToFind = dim == 3 ? 2 : 1;
219
220 // keep track of facets handled already while searching
221 unsigned int foundCounter = 0;
222 std::vector<unsigned int> handledFacets;
223 for (const auto& is : intersections(gridView, element))
224 {
225 // skip if on previous intersection again
226 if (is.indexInInside() == prevIntersection.indexInInside())
227 continue;
228
229 // skip intersection if handled already (necessary for e.g. Dune::FoamGrid)
230 if (std::count(handledFacets.begin(), handledFacets.end(), is.indexInInside()))
231 continue;
232
233 // determine all vertex indices of this face
234 const auto numCorners = is.geometry().corners();
235 std::vector<GridIndexType> faceVertexIndices(numCorners);
236 for (int i = 0; i < numCorners; ++i)
237 faceVertexIndices[i] = vertexMapper.subIndex( element,
238 refElement.subEntity(is.indexInInside(), 1, i, dim),
239 dim );
240
241 // we found another intersection at the vertex if it contains the vertex around which we rotate
242 if (std::find(faceVertexIndices.begin(), faceVertexIndices.end(), vIdxGlobal) != faceVertexIndices.end())
243 {
244 foundCounter++;
245 handledFacets.push_back(is.indexInInside());
246
247 // keep searching in the outside element only if ...
248 // ... this is a not (processor) boundary
249 if (!is.neighbor())
250 continue;
251
252 // ... this face does not lie on the facet grid
253 if (codimOneGridAdapter.composeFacetElement(faceVertexIndices))
254 continue;
255
256 // ... outside element is not contained yet in path
257 const auto outsideElement = is.outside();
258 const auto outsideElemIdx = elementMapper.index(outsideElement);
259 if (std::find(path.begin(), path.end(), outsideElemIdx) == path.end())
260 {
261 const auto idxInOutside = is.indexInOutside();
262 const auto outsideRefElement = referenceElement(outsideElement);
263 path.push_back(outsideElemIdx);
264
265 // on 2d grids, there is only going to be one more
266 for (const auto& outsideIs : intersections(gridView, outsideElement))
267 {
268 if (outsideIs.indexInInside() == idxInOutside)
269 {
270 // if this is the last intersection to find, return
271 if (foundCounter == numIsToFind)
272 return continuePathSearch_(path,
273 gridView,
274 elementMapper,
275 vertexMapper,
276 codimOneGridAdapter,
277 outsideElement,
278 outsideRefElement,
279 outsideIs,
280 vIdxGlobal);
281
282 // otherwise, put one more search on stack
283 else
284 continuePathSearch_(path,
285 gridView,
286 elementMapper,
287 vertexMapper,
288 codimOneGridAdapter,
289 outsideElement,
290 outsideRefElement,
291 outsideIs,
292 vIdxGlobal);
293 }
294 }
295 }
296 }
297 }
298
299 if (foundCounter != numIsToFind)
300 DUNE_THROW(Dune::InvalidStateException, "Found " << foundCounter << " instead of " << numIsToFind << " intersections around vertex");
301 }
302};
303
304} // end namespace Dumux
305
306#endif
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:40
bool composeFacetElement(const IndexStorage &bulkVertexIndices) const
Returns true if a given set of bulk vertex indices make up a facet grid element.
Definition: codimonegridadapter.hh:180
Specialization of the enrichment helper class for 2d grids. In this case, we look for two-dimensional...
Definition: enrichmenthelper.hh:44
static std::size_t enrich(IndexMap &indexMap, const std::vector< bool > &vertexMarkers, const GridView &gridView, const MCMGMapper &vertexMapper, const MCMGMapper &elementMapper, const CodimOneGridView &codimOneGridView, const CodimOneGridAdapter &codimOneGridAdapter)
Enriches the dof map subject to a (dim-1)-dimensional grid.
Definition: enrichmenthelper.hh:77
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27