3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
intersectingentities.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 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
25#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
26
27#include <cmath>
28#include <type_traits>
29#include <vector>
30#include <algorithm>
31#include <limits>
32
33#include <dune/common/fvector.hh>
34
35#include <dumux/common/math.hh>
40
41namespace Dumux {
42
48template<int dimworld, class CoordTypeA, class CoordTypeB = CoordTypeA>
50{
51public:
52 using ctype = typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
53 static constexpr int dimensionworld = dimworld;
54 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
55
56 template<class Corners>
57 explicit IntersectionInfo(std::size_t a, std::size_t b, Corners&& c)
58 : a_(a)
59 , b_(b)
60 , corners_(c.begin(), c.end())
61 {}
62
64 std::size_t first() const
65 { return a_; }
66
68 std::size_t second() const
69 { return b_; }
70
72 const std::vector<GlobalPosition>& corners() const
73 { return corners_; }
74
79 bool cornersMatch(const std::vector<GlobalPosition>& otherCorners) const
80 {
81 if (otherCorners.size() != corners_.size())
82 return false;
83
84 using std::max;
85 ctype eps2 = std::numeric_limits<ctype>::min();
86 for (int i = 1; i < corners_.size(); ++i)
87 eps2 = max(eps2, (corners_[i] - corners_[0]).two_norm2());
88
89 // We use a base epsilon of 1.5e-7 for comparisons of lengths.
90 // Since here we compare squared lengths, we multiply by its square.
91 eps2 *= 1.5e-7*1.5e-7;
92
93 for (int i = 0; i < corners_.size(); ++i)
94 // early return if none of the other corners are equal to this corner
95 if (std::none_of(otherCorners.begin(),
96 otherCorners.end(),
97 [&] (const auto& other) { return (corners_[i] - other).two_norm2() < eps2; }))
98 return false;
99
100 return true;
101 }
102
103private:
104 std::size_t a_, b_;
105 std::vector<GlobalPosition> corners_;
106};
107
112template<class EntitySet, class ctype, int dimworld>
113inline std::vector<std::size_t>
114intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
115 const BoundingBoxTree<EntitySet>& tree,
116 bool isCartesianGrid = false)
117{
118 // Call the recursive find function to find candidates
119 std::vector<std::size_t> entities;
120 intersectingEntities(point, tree, tree.numBoundingBoxes() - 1, entities, isCartesianGrid);
121 return entities;
122}
123
128template<class EntitySet, class ctype, int dimworld>
129void intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
130 const BoundingBoxTree<EntitySet>& tree,
131 std::size_t node,
132 std::vector<std::size_t>& entities,
133 bool isCartesianGrid = false)
134{
135 // Get the bounding box for the current node
136 const auto& bBox = tree.getBoundingBoxNode(node);
137
138 // if the point is not in the bounding box we can stop
139 if (!intersectsPointBoundingBox(point, tree.getBoundingBoxCoordinates(node))) return;
140
141 // now we know it's inside
142 // if the box is a leaf do the primitive test.
143 else if (tree.isLeaf(bBox, node))
144 {
145 const std::size_t entityIdx = bBox.child1;
146 // for structured cube grids skip the primitive test
147 if (isCartesianGrid)
148 entities.push_back(entityIdx);
149 else
150 {
151 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
152 // if the primitive is positive it intersects the actual geometry, add the entity to the list
153 if (intersectsPointGeometry(point, geometry))
154 entities.push_back(entityIdx);
155 }
156 }
157
158 // No leaf. Check both children nodes.
159 else
160 {
161 intersectingEntities(point, tree, bBox.child0, entities, isCartesianGrid);
162 intersectingEntities(point, tree, bBox.child1, entities, isCartesianGrid);
163 }
164}
165
170template<class Geometry, class EntitySet>
171inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
172intersectingEntities(const Geometry& geometry,
173 const BoundingBoxTree<EntitySet>& tree)
174{
175 // check if the world dimensions match
176 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
177 "Can only intersect geometry and bounding box tree of same world dimension");
178
179 // Create data structure for return type
180 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
182 static constexpr int dimworld = Geometry::coorddimension;
183
184 // compute the bounding box of the given geometry
185 std::array<ctype, 2*Geometry::coorddimension> bBox;
186 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
187
188 // Get coordinates of first vertex
189 auto corner = geometry.corner(0);
190 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
191 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
192
193 // Compute the min and max over the remaining vertices
194 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
195 {
196 corner = geometry.corner(cornerIdx);
197 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
198 {
199 using std::max;
200 using std::min;
201 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
202 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
203 }
204 }
205
206 // Call the recursive find function to find candidates
207 intersectingEntities(geometry, tree,
208 bBox, tree.numBoundingBoxes() - 1,
209 intersections);
210
211 return intersections;
212}
213
218template<class Geometry, class EntitySet>
219void intersectingEntities(const Geometry& geometry,
220 const BoundingBoxTree<EntitySet>& tree,
221 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
222 std::size_t nodeIdx,
223 std::vector<IntersectionInfo<Geometry::coorddimension,
224 typename Geometry::ctype,
225 typename EntitySet::ctype>>& intersections)
226{
227 // if the two bounding boxes don't intersect we can stop searching
228 static constexpr int dimworld = Geometry::coorddimension;
229 if (!intersectsBoundingBoxBoundingBox<dimworld>(bBox.data(), tree.getBoundingBoxCoordinates(nodeIdx)))
230 return;
231
232 // get node info for current bounding box node
233 const auto& bBoxNode = tree.getBoundingBoxNode(nodeIdx);
234
235 // if the box is a leaf do the primitive test.
236 if (tree.isLeaf(bBoxNode, nodeIdx))
237 {
238 // eIdxA is always 0 since we intersect with exactly one geometry
239 const auto eIdxA = 0;
240 const auto eIdxB = bBoxNode.child1;
241
242 const auto geometryTree = tree.entitySet().entity(eIdxB).geometry();
243 using GeometryTree = std::decay_t<decltype(geometryTree)>;
245 using IntersectionAlgorithm = GeometryIntersection<Geometry, GeometryTree, Policy>;
246 using Intersection = typename IntersectionAlgorithm::Intersection;
247 Intersection intersection;
248
249 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
250 {
251 static constexpr int dimIntersection = Policy::dimIntersection;
252 if (dimIntersection >= 2)
253 {
254 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
255 for (unsigned int i = 0; i < triangulation.size(); ++i)
256 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
257 }
258 else
259 intersections.emplace_back(eIdxA, eIdxB, intersection);
260 }
261 }
262
263 // No leaf. Check both children nodes.
264 else
265 {
266 intersectingEntities(geometry, tree, bBox, bBoxNode.child0, intersections);
267 intersectingEntities(geometry, tree, bBox, bBoxNode.child1, intersections);
268 }
269}
270
275template<class EntitySet0, class EntitySet1>
276inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
278 const BoundingBoxTree<EntitySet1>& treeB)
279{
280 // check if the world dimensions match
281 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
282 "Can only intersect bounding box trees of same world dimension");
283
284 // Create data structure for return type
285 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
286
287 // Call the recursive find function to find candidates
288 intersectingEntities(treeA, treeB,
289 treeA.numBoundingBoxes() - 1,
290 treeB.numBoundingBoxes() - 1,
291 intersections);
292
293 return intersections;
294}
295
300template<class EntitySet0, class EntitySet1>
302 const BoundingBoxTree<EntitySet1>& treeB,
303 std::size_t nodeA, std::size_t nodeB,
304 std::vector<IntersectionInfo<EntitySet0::dimensionworld,
305 typename EntitySet0::ctype,
306 typename EntitySet1::ctype>>& intersections)
307{
308 // Get the bounding box for the current node
309 const auto& bBoxA = treeA.getBoundingBoxNode(nodeA);
310 const auto& bBoxB = treeB.getBoundingBoxNode(nodeB);
311
312 // if the two bounding boxes of the current nodes don't intersect we can stop searching
313 static constexpr int dimworld = EntitySet0::dimensionworld;
314 if (!intersectsBoundingBoxBoundingBox<dimworld>(treeA.getBoundingBoxCoordinates(nodeA),
315 treeB.getBoundingBoxCoordinates(nodeB)))
316 return;
317
318 // Check if we have a leaf in treeA or treeB
319 const bool isLeafA = treeA.isLeaf(bBoxA, nodeA);
320 const bool isLeafB = treeB.isLeaf(bBoxB, nodeB);
321
322 // If both boxes are leaves do the primitive test
323 if (isLeafA && isLeafB)
324 {
325 const auto eIdxA = bBoxA.child1;
326 const auto eIdxB = bBoxB.child1;
327
328 const auto geometryA = treeA.entitySet().entity(eIdxA).geometry();
329 const auto geometryB = treeB.entitySet().entity(eIdxB).geometry();
330
331 using GeometryA = std::decay_t<decltype(geometryA)>;
332 using GeometryB = std::decay_t<decltype(geometryB)>;
334 using IntersectionAlgorithm = GeometryIntersection<GeometryA, GeometryB, Policy>;
335 using Intersection = typename IntersectionAlgorithm::Intersection;
336
337 if (Intersection intersection; IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
338 {
339 static constexpr int dimIntersection = Policy::dimIntersection;
340
341 // intersection is returned as a point cloud for dim >= 2
342 // so we have to triangulate first
343 if (dimIntersection >= 2)
344 {
345 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
346 for (unsigned int i = 0; i < triangulation.size(); ++i)
347 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
348 }
349 else
350 intersections.emplace_back(eIdxA, eIdxB, intersection);
351 }
352 }
353
354 // if we reached the leaf in treeA, just continue in treeB
355 else if (isLeafA)
356 {
357 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
358 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
359 }
360
361 // if we reached the leaf in treeB, just continue in treeA
362 else if (isLeafB)
363 {
364 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
365 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
366 }
367
368 // we know now that both trees didn't reach the leaf yet so
369 // we continue with the larger tree first (bigger node number)
370 else if (nodeA > nodeB)
371 {
372 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
373 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
374 }
375 else
376 {
377 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
378 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
379 }
380}
381
390template<class ctype, int dimworld>
391inline std::size_t intersectingEntityCartesianGrid(const Dune::FieldVector<ctype, dimworld>& point,
392 const Dune::FieldVector<ctype, dimworld>& min,
393 const Dune::FieldVector<ctype, dimworld>& max,
394 const std::array<int, std::size_t(dimworld)>& cells)
395{
396 std::size_t index = 0;
397 for (int i = 0; i < dimworld; ++i)
398 {
399 using std::clamp; using std::floor;
400 ctype dimOffset = clamp<ctype>(floor((point[i]-min[i])*cells[i]/(max[i]-min[i])), 0.0, cells[i]-1);
401 for (int j = 0; j < i; ++j)
402 dimOffset *= cells[j];
403 index += static_cast<std::size_t>(dimOffset);
404 }
405 return index;
406}
407
408} // end namespace Dumux
409
410#endif
Define some often used mathematical functions.
A class for collision detection of two geometries and computation of intersection corners.
An axis-aligned bounding box volume hierarchy for dune grids.
Detect if a point intersects a geometry.
Functionality to triangulate point clouds.
std::size_t intersectingEntityCartesianGrid(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &min, const Dune::FieldVector< ctype, dimworld > &max, const std::array< int, std::size_t(dimworld)> &cells)
Compute the index of the intersecting element of a Cartesian grid with a point The grid is given by t...
Definition: intersectingentities.hh:391
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:40
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: intersectingentities.hh:114
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: boundingboxtree.hh:306
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:116
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:68
const ctype * getBoundingBoxCoordinates(std::size_t nodeIdx) const
Get an existing bounding box for a given node.
Definition: boundingboxtree.hh:147
bool isLeaf(const BoundingBoxNode &node, std::size_t nodeIdx) const
Definition: boundingboxtree.hh:156
const EntitySet & entitySet() const
the entity set this tree was built with
Definition: boundingboxtree.hh:135
std::size_t numBoundingBoxes() const
Get the number of bounding boxes currently in the tree.
Definition: boundingboxtree.hh:151
const BoundingBoxNode & getBoundingBoxNode(std::size_t nodeIdx) const
Interface to be used by other bounding box trees.
Definition: boundingboxtree.hh:143
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:219
An intersection object resulting from the intersection of two primitives in an entity set.
Definition: intersectingentities.hh:50
static constexpr int dimensionworld
Definition: intersectingentities.hh:53
bool cornersMatch(const std::vector< GlobalPosition > &otherCorners) const
Check if the corners of this intersection match with the given corners.
Definition: intersectingentities.hh:79
const std::vector< GlobalPosition > & corners() const
Get the corners of the intersection geometry.
Definition: intersectingentities.hh:72
IntersectionInfo(std::size_t a, std::size_t b, Corners &&c)
Definition: intersectingentities.hh:57
Dune::FieldVector< ctype, dimworld > GlobalPosition
Definition: intersectingentities.hh:54
std::size_t second() const
Get the index of the intersecting entity belonging to the other grid.
Definition: intersectingentities.hh:68
std::size_t first() const
Get the index of the intersecting entity belonging to this grid.
Definition: intersectingentities.hh:64
typename Dune::PromotionTraits< CoordTypeA, CoordTypeB >::PromotedType ctype
Definition: intersectingentities.hh:52