3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
19
25#ifndef DUMUX_VERTEX_ENRICHMENT_HELPER_HH
26#define DUMUX_VERTEX_ENRICHMENT_HELPER_HH
27
28#include <cmath>
29#include <vector>
30#include <unordered_map>
31#include <unordered_set>
32
33#include <dune/common/exceptions.hh>
34#include <dune/common/reservedvector.hh>
35
36#include <dune/grid/common/mcmgmapper.hh>
37
38#include <dumux/common/math.hh>
40
41namespace Dumux {
42
54template< class GridView, class CodimOneGridView >
56{
57 static constexpr int dim = GridView::dimension;
58 static constexpr int dimWorld = GridView::dimensionworld;
59 static_assert(dim == 2 || dim == 3, "Grid dimension must be two or three");
60 static_assert(dimWorld == int(CodimOneGridView::dimensionworld), "world dimension mismatch");
61
62 using Intersection = typename GridView::Intersection;
63 using MCMGMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
64 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
65 using Element = typename GridView::template Codim<0>::Entity;
66 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
67
68 using ElementPath = std::vector< GridIndexType >;
69 using NodalElementPaths = std::vector< std::vector<ElementPath> >;
70
71public:
72
88 template< class IndexMap, class CodimOneGridAdapter >
89 static std::size_t enrich(IndexMap& indexMap,
90 const std::vector<bool>& vertexMarkers,
91 const GridView& gridView,
92 const MCMGMapper& vertexMapper,
93 const MCMGMapper& elementMapper,
94 const CodimOneGridView& codimOneGridView,
95 const CodimOneGridAdapter& codimOneGridAdapter)
96 {
97 // determine the bulk element paths around vertices marked for enrichment
98 NodalElementPaths nodalPaths(gridView.size(dim));
99 for (const auto& e : elements(gridView))
100 {
101 const auto eIdx = elementMapper.index(e);
102
103 std::vector<unsigned int> handledFacets;
104 const auto refElement = referenceElement(e);
105 for (const auto& is : intersections(gridView, e))
106 {
107 // do not start path search on boundary intersections
108 if (!is.neighbor())
109 continue;
110
111 // if facet has been handled already, skip rest (necesary for e.g. Dune::FoamGrid)
112 if (std::find(handledFacets.begin(), handledFacets.end(), is.indexInInside()) != handledFacets.end())
113 continue;
114
115 // first, get indices of all vertices of this intersection
116 const auto numCorners = is.geometry().corners();
117 std::vector<GridIndexType> faceVertexIndices(numCorners);
118 for (int i = 0; i < numCorners; ++i)
119 faceVertexIndices[i] = vertexMapper.subIndex( e,
120 refElement.subEntity(is.indexInInside(), 1, i, dim),
121 dim );
122
123 // start path search on intersections that lie on a facet grid element
124 if (codimOneGridAdapter.composeFacetElement(faceVertexIndices))
125 {
126 handledFacets.push_back(is.indexInInside());
127 for (int i = 0; i < numCorners; ++i)
128 {
129 // set up path only around those vertices that are to be enriched
130 const auto vIdxGlobal = faceVertexIndices[i];
131 if (!vertexMarkers[vIdxGlobal])
132 continue;
133
134 // construct the path only if element is not yet contained in another path
135 bool found = false;
136 for (const auto& path : nodalPaths[vIdxGlobal])
137 if (std::find(path.begin(), path.end(), eIdx) != path.end())
138 { found = true; break; }
139
140 if (!found)
141 {
142 ElementPath path;
143
144 // Reserve enough memory so that reallocation during the recursive search
145 // does not happen and references are invalidated. Memory is not an issue
146 // here as after enrichment the container is thrown away again.
147 // TODO improve this!?
148 path.reserve(150);
149 path.push_back(eIdx);
150 continuePathSearch_(path, gridView, elementMapper, vertexMapper, codimOneGridAdapter, e, refElement, is, vIdxGlobal);
151 nodalPaths[vIdxGlobal].emplace_back(std::move(path));
152 }
153 }
154 }
155 }
156 }
157
158 // determine the offsets for each bulk vertex index on the basis of the paths found per vertex
159 std::vector<std::size_t> bulkVertexIndexOffsets(gridView.size(dim), 0);
160 for (const auto& v : vertices(gridView))
161 {
162 const auto vIdx = vertexMapper.index(v);
163 if (vertexMarkers[vIdx])
164 bulkVertexIndexOffsets[vIdx] = nodalPaths[vIdx].size()-1;
165 }
166
167 // ... and accumulate the offsets
168 std::size_t sumOffset = 0;
169 std::size_t size = 0;
170 for (auto& nodalOffset : bulkVertexIndexOffsets)
171 {
172 const auto os = nodalOffset;
173 nodalOffset = sumOffset;
174 sumOffset += os;
175 size += (os == 0) ? 1 : os + 1;
176 }
177
178 // Now, finally set up the new index map
179 for (const auto& e : elements(gridView))
180 {
181 const auto& eg = e.geometry();
182 const auto eIdx = elementMapper.index(e);
183 for (int i = 0; i < eg.corners(); ++i)
184 {
185 const auto origVIdx = vertexMapper.subIndex(e, i, dim);
186
187 // it the node itself is not enriched, simply add offset
188 if (!vertexMarkers[origVIdx])
189 indexMap[eIdx][i] += bulkVertexIndexOffsets[origVIdx];
190
191 // find the local index of the path the element belongs to
192 else
193 {
194 bool found = false;
195 const auto& paths = nodalPaths[ origVIdx ];
196 for (int pathIdx = 0; pathIdx < paths.size(); ++pathIdx)
197 {
198 const auto& curPath = paths[pathIdx];
199 if ( std::find(curPath.begin(), curPath.end(), eIdx) != curPath.end() )
200 {
201 indexMap[eIdx][i] += bulkVertexIndexOffsets[origVIdx] + pathIdx;
202 found = true;
203 break;
204 }
205 }
206
207 if (!found)
208 DUNE_THROW(Dune::InvalidStateException, "Element not found in any path");
209 }
210 }
211 }
212
213 // return the number of dofs after enrichment
214 return size;
215 }
216
217private:
218 template< class ReferenceElement, class CodimOneGridAdapter >
219 static void continuePathSearch_(ElementPath& path,
220 const GridView& gridView,
221 const MCMGMapper& elementMapper,
222 const MCMGMapper& vertexMapper,
223 const CodimOneGridAdapter& codimOneGridAdapter,
224 const Element& element,
225 const ReferenceElement& refElement,
226 const Intersection& prevIntersection,
227 GridIndexType vIdxGlobal)
228 {
229 // on 2d/3d grids, we need to find 1/2 more intersections at vertex
230 static constexpr int numIsToFind = dim == 3 ? 2 : 1;
231
232 // keep track of facets handled already while searching
233 unsigned int foundCounter = 0;
234 std::vector<unsigned int> handledFacets;
235 for (const auto& is : intersections(gridView, element))
236 {
237 // skip if on previous intersection again
238 if (is.indexInInside() == prevIntersection.indexInInside())
239 continue;
240
241 // skip intersection if handled already (necessary for e.g. Dune::FoamGrid)
242 if (std::count(handledFacets.begin(), handledFacets.end(), is.indexInInside()))
243 continue;
244
245 // determine all vertex indices of this face
246 const auto numCorners = is.geometry().corners();
247 std::vector<GridIndexType> faceVertexIndices(numCorners);
248 for (int i = 0; i < numCorners; ++i)
249 faceVertexIndices[i] = vertexMapper.subIndex( element,
250 refElement.subEntity(is.indexInInside(), 1, i, dim),
251 dim );
252
253 // we found another intersection at the vertex if it contains the vertex around which we rotate
254 if (std::find(faceVertexIndices.begin(), faceVertexIndices.end(), vIdxGlobal) != faceVertexIndices.end())
255 {
256 foundCounter++;
257 handledFacets.push_back(is.indexInInside());
258
259 // keep searching in the outside element only if ...
260 // ... this is a not (processor) boundary
261 if (!is.neighbor())
262 continue;
263
264 // ... this face does not lie on the facet grid
265 if (codimOneGridAdapter.composeFacetElement(faceVertexIndices))
266 continue;
267
268 // ... outside element is not contained yet in path
269 const auto outsideElement = is.outside();
270 const auto outsideElemIdx = elementMapper.index(outsideElement);
271 if (std::find(path.begin(), path.end(), outsideElemIdx) == path.end())
272 {
273 const auto idxInOutside = is.indexInOutside();
274 const auto outsideRefElement = referenceElement(outsideElement);
275 path.push_back(outsideElemIdx);
276
277 // on 2d grids, there is only going to be one more
278 for (const auto& outsideIs : intersections(gridView, outsideElement))
279 {
280 if (outsideIs.indexInInside() == idxInOutside)
281 {
282 // if this is the last intersection to find, return
283 if (foundCounter == numIsToFind)
284 return continuePathSearch_(path,
285 gridView,
286 elementMapper,
287 vertexMapper,
288 codimOneGridAdapter,
289 outsideElement,
290 outsideRefElement,
291 outsideIs,
292 vIdxGlobal);
293
294 // otherwise, put one more search on stack
295 else
296 continuePathSearch_(path,
297 gridView,
298 elementMapper,
299 vertexMapper,
300 codimOneGridAdapter,
301 outsideElement,
302 outsideRefElement,
303 outsideIs,
304 vIdxGlobal);
305 }
306 }
307 }
308 }
309 }
310
311 if (foundCounter != numIsToFind)
312 DUNE_THROW(Dune::InvalidStateException, "Found " << foundCounter << " instead of " << numIsToFind << " intersections around vertex");
313 }
314};
315
316} // end namespace Dumux
317
318#endif
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
Definition: adapt.hh:29
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:52
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:192
Specialization of the enrichment helper class for 2d grids. In this case, we look for two-dimensional...
Definition: enrichmenthelper.hh:56
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:89