version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
13#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
14
15#include <cmath>
16#include <type_traits>
17#include <vector>
18#include <algorithm>
19#include <limits>
20
21#include <dune/common/fvector.hh>
22
23#include <dumux/common/math.hh>
28
29namespace Dumux {
30
36template<int dimworld, class CoordTypeA, class CoordTypeB = CoordTypeA>
38{
39public:
40 using ctype = typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
41 static constexpr int dimensionworld = dimworld;
42 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
43
44 template<class Corners>
45 explicit IntersectionInfo(std::size_t a, std::size_t b, Corners&& c)
46 : a_(a)
47 , b_(b)
48 , corners_(c.begin(), c.end())
49 {}
50
52 std::size_t first() const
53 { return a_; }
54
56 std::size_t second() const
57 { return b_; }
58
60 const std::vector<GlobalPosition>& corners() const
61 { return corners_; }
62
67 bool cornersMatch(const std::vector<GlobalPosition>& otherCorners) const
68 {
69 if (otherCorners.size() != corners_.size())
70 return false;
71
72 using std::max;
73 ctype eps2 = std::numeric_limits<ctype>::min();
74 for (int i = 1; i < corners_.size(); ++i)
75 eps2 = max(eps2, (corners_[i] - corners_[0]).two_norm2());
76
77 // We use a base epsilon of 1.5e-7 for comparisons of lengths.
78 // Since here we compare squared lengths, we multiply by its square.
79 eps2 *= 1.5e-7*1.5e-7;
80
81 for (int i = 0; i < corners_.size(); ++i)
82 // early return if none of the other corners are equal to this corner
83 if (std::none_of(otherCorners.begin(),
84 otherCorners.end(),
85 [&] (const auto& other) { return (corners_[i] - other).two_norm2() < eps2; }))
86 return false;
87
88 return true;
89 }
90
91private:
92 std::size_t a_, b_;
93 std::vector<GlobalPosition> corners_;
94};
95
100template<class EntitySet, class ctype, int dimworld>
101inline std::vector<std::size_t>
102intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
103 const BoundingBoxTree<EntitySet>& tree,
104 bool isCartesianGrid = false)
105{
106 // Call the recursive find function to find candidates
107 std::vector<std::size_t> entities;
108 intersectingEntities(point, tree, tree.numBoundingBoxes() - 1, entities, isCartesianGrid);
109 return entities;
110}
111
116template<class EntitySet, class ctype, int dimworld>
117void intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
118 const BoundingBoxTree<EntitySet>& tree,
119 std::size_t node,
120 std::vector<std::size_t>& entities,
121 bool isCartesianGrid = false)
122{
123 // Get the bounding box for the current node
124 const auto& bBox = tree.getBoundingBoxNode(node);
125
126 // if the point is not in the bounding box we can stop
127 if (!intersectsPointBoundingBox(point, tree.getBoundingBoxCoordinates(node))) return;
128
129 // now we know it's inside
130 // if the box is a leaf do the primitive test.
131 else if (tree.isLeaf(bBox, node))
132 {
133 const std::size_t entityIdx = bBox.child1;
134 // for structured cube grids skip the primitive test
135 if (isCartesianGrid)
136 entities.push_back(entityIdx);
137 else
138 {
139 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
140 // if the primitive is positive it intersects the actual geometry, add the entity to the list
141 if (intersectsPointGeometry(point, geometry))
142 entities.push_back(entityIdx);
143 }
144 }
145
146 // No leaf. Check both children nodes.
147 else
148 {
149 intersectingEntities(point, tree, bBox.child0, entities, isCartesianGrid);
150 intersectingEntities(point, tree, bBox.child1, entities, isCartesianGrid);
151 }
152}
153
158template<class Geometry, class EntitySet>
159inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
160intersectingEntities(const Geometry& geometry,
161 const BoundingBoxTree<EntitySet>& tree)
162{
164 return intersectingEntities(geometry, tree, IP{});
165}
166
171template<class Geometry, class EntitySet, class IntersectionPolicy>
172inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
173intersectingEntities(const Geometry& geometry,
174 const BoundingBoxTree<EntitySet>& tree,
175 IntersectionPolicy intersectionPolicy)
176{
177 // check if the world dimensions match
178 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
179 "Can only intersect geometry and bounding box tree of same world dimension");
180
181 // Create data structure for return type
182 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
184 static constexpr int dimworld = Geometry::coorddimension;
185
186 // compute the bounding box of the given geometry
187 std::array<ctype, 2*Geometry::coorddimension> bBox;
188 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
189
190 // Get coordinates of first vertex
191 auto corner = geometry.corner(0);
192 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
193 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
194
195 // Compute the min and max over the remaining vertices
196 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
197 {
198 corner = geometry.corner(cornerIdx);
199 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
200 {
201 using std::max;
202 using std::min;
203 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
204 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
205 }
206 }
207
208 // Call the recursive find function to find candidates
209 intersectingEntities(geometry, tree,
210 bBox, tree.numBoundingBoxes() - 1,
211 intersections, intersectionPolicy);
212
213 return intersections;
214}
215
220template<class Geometry, class EntitySet>
221void intersectingEntities(const Geometry& geometry,
222 const BoundingBoxTree<EntitySet>& tree,
223 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
224 std::size_t nodeIdx,
225 std::vector<IntersectionInfo<Geometry::coorddimension,
226 typename Geometry::ctype,
227 typename EntitySet::ctype>>& intersections)
228{
230 intersectingEntities(geometry, tree, bBox, nodeIdx, intersections, IP{});
231}
236template<class Geometry, class EntitySet, class IntersectionPolicy>
237void intersectingEntities(const Geometry& geometry,
238 const BoundingBoxTree<EntitySet>& tree,
239 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
240 std::size_t nodeIdx,
241 std::vector<IntersectionInfo<Geometry::coorddimension,
242 typename Geometry::ctype,
243 typename EntitySet::ctype>>& intersections,
244 IntersectionPolicy intersectionPolicy)
245{
246 // if the two bounding boxes don't intersect we can stop searching
247 static constexpr int dimworld = Geometry::coorddimension;
248 if (!intersectsBoundingBoxBoundingBox<dimworld>(bBox.data(), tree.getBoundingBoxCoordinates(nodeIdx)))
249 return;
250
251 // get node info for current bounding box node
252 const auto& bBoxNode = tree.getBoundingBoxNode(nodeIdx);
253
254 // if the box is a leaf do the primitive test.
255 if (tree.isLeaf(bBoxNode, nodeIdx))
256 {
257 // eIdxA is always 0 since we intersect with exactly one geometry
258 const auto eIdxA = 0;
259 const auto eIdxB = bBoxNode.child1;
260
261 const auto geometryTree = tree.entitySet().entity(eIdxB).geometry();
262 using GeometryTree = std::decay_t<decltype(geometryTree)>;
264 using Intersection = typename IntersectionAlgorithm::Intersection;
265 Intersection intersection;
266
267 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
268 {
269 static constexpr int dimIntersection = IntersectionPolicy::dimIntersection;
270 if constexpr (dimIntersection >= 2)
271 {
272 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
273 for (unsigned int i = 0; i < triangulation.size(); ++i)
274 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
275 }
276 else
277 intersections.emplace_back(eIdxA, eIdxB, intersection);
278 }
279 }
280
281 // No leaf. Check both children nodes.
282 else
283 {
284 intersectingEntities(geometry, tree, bBox, bBoxNode.child0, intersections);
285 intersectingEntities(geometry, tree, bBox, bBoxNode.child1, intersections);
286 }
287}
288
293template<class EntitySet0, class EntitySet1>
294inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
296 const BoundingBoxTree<EntitySet1>& treeB)
297{
299 return intersectingEntities(treeA, treeB, IP{});
300}
301
306template<class EntitySet0, class EntitySet1, class IntersectionPolicy>
307inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
309 const BoundingBoxTree<EntitySet1>& treeB,
310 IntersectionPolicy intersectionPolicy)
311{
312 // check if the world dimensions match
313 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
314 "Can only intersect bounding box trees of same world dimension");
315
316 // Create data structure for return type
317 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
318
319 // Call the recursive find function to find candidates
320 intersectingEntities(treeA, treeB,
321 treeA.numBoundingBoxes() - 1,
322 treeB.numBoundingBoxes() - 1,
323 intersections, intersectionPolicy);
324
325 return intersections;
326}
327
332template<class EntitySet0, class EntitySet1>
334 const BoundingBoxTree<EntitySet1>& treeB,
335 std::size_t nodeA, std::size_t nodeB,
336 std::vector<IntersectionInfo<EntitySet0::dimensionworld,
337 typename EntitySet0::ctype,
338 typename EntitySet1::ctype>>& intersections)
339{
341 intersectingEntities(treeA, treeB, nodeA, nodeB, intersections, IP{});
342}
343
348template<class EntitySet0, class EntitySet1, class IntersectionPolicy>
350 const BoundingBoxTree<EntitySet1>& treeB,
351 std::size_t nodeA, std::size_t nodeB,
352 std::vector<IntersectionInfo<EntitySet0::dimensionworld,
353 typename EntitySet0::ctype,
354 typename EntitySet1::ctype>>& intersections,
355 IntersectionPolicy intersectionPolicy)
356{
357 // Get the bounding box for the current node
358 const auto& bBoxA = treeA.getBoundingBoxNode(nodeA);
359 const auto& bBoxB = treeB.getBoundingBoxNode(nodeB);
360
361 // if the two bounding boxes of the current nodes don't intersect we can stop searching
362 static constexpr int dimworld = EntitySet0::dimensionworld;
363 if (!intersectsBoundingBoxBoundingBox<dimworld>(treeA.getBoundingBoxCoordinates(nodeA),
364 treeB.getBoundingBoxCoordinates(nodeB)))
365 return;
366
367 // Check if we have a leaf in treeA or treeB
368 const bool isLeafA = treeA.isLeaf(bBoxA, nodeA);
369 const bool isLeafB = treeB.isLeaf(bBoxB, nodeB);
370
371 // If both boxes are leaves do the primitive test
372 if (isLeafA && isLeafB)
373 {
374 const auto eIdxA = bBoxA.child1;
375 const auto eIdxB = bBoxB.child1;
376
377 const auto geometryA = treeA.entitySet().entity(eIdxA).geometry();
378 const auto geometryB = treeB.entitySet().entity(eIdxB).geometry();
379
380 using GeometryA = std::decay_t<decltype(geometryA)>;
381 using GeometryB = std::decay_t<decltype(geometryB)>;
383 using Intersection = typename IntersectionAlgorithm::Intersection;
384
385 if (Intersection intersection; IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
386 {
387 static constexpr int dimIntersection = IntersectionPolicy::dimIntersection;
388
389 // intersection is returned as a point cloud for dim >= 2
390 // so we have to triangulate first
391 if constexpr (dimIntersection >= 2)
392 {
393 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
394 for (unsigned int i = 0; i < triangulation.size(); ++i)
395 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
396 }
397 else
398 intersections.emplace_back(eIdxA, eIdxB, intersection);
399 }
400 }
401
402 // if we reached the leaf in treeA, just continue in treeB
403 else if (isLeafA)
404 {
405 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
406 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
407 }
408
409 // if we reached the leaf in treeB, just continue in treeA
410 else if (isLeafB)
411 {
412 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
413 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
414 }
415
416 // we know now that both trees didn't reach the leaf yet so
417 // we continue with the larger tree first (bigger node number)
418 else if (nodeA > nodeB)
419 {
420 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
421 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
422 }
423 else
424 {
425 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
426 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
427 }
428}
429
438template<class ctype, int dimworld>
439inline std::size_t intersectingEntityCartesianGrid(const Dune::FieldVector<ctype, dimworld>& point,
440 const Dune::FieldVector<ctype, dimworld>& min,
441 const Dune::FieldVector<ctype, dimworld>& max,
442 const std::array<int, std::size_t(dimworld)>& cells)
443{
444 std::size_t index = 0;
445 for (int i = 0; i < dimworld; ++i)
446 {
447 using std::clamp; using std::floor;
448 ctype dimOffset = clamp<ctype>(floor((point[i]-min[i])*cells[i]/(max[i]-min[i])), 0.0, cells[i]-1);
449 for (int j = 0; j < i; ++j)
450 dimOffset *= cells[j];
451 index += static_cast<std::size_t>(dimOffset);
452 }
453 return index;
454}
455
456} // end namespace Dumux
457
458#endif
An axis-aligned bounding box volume hierarchy for dune grids.
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:56
const ctype * getBoundingBoxCoordinates(std::size_t nodeIdx) const
Get an existing bounding box for a given node.
Definition: boundingboxtree.hh:135
bool isLeaf(const BoundingBoxNode &node, std::size_t nodeIdx) const
Definition: boundingboxtree.hh:144
const EntitySet & entitySet() const
the entity set this tree was built with
Definition: boundingboxtree.hh:123
std::size_t numBoundingBoxes() const
Get the number of bounding boxes currently in the tree.
Definition: boundingboxtree.hh:139
const BoundingBoxNode & getBoundingBoxNode(std::size_t nodeIdx) const
Interface to be used by other bounding box trees.
Definition: boundingboxtree.hh:131
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:207
An intersection object resulting from the intersection of two primitives in an entity set.
Definition: intersectingentities.hh:38
static constexpr int dimensionworld
Definition: intersectingentities.hh:41
bool cornersMatch(const std::vector< GlobalPosition > &otherCorners) const
Check if the corners of this intersection match with the given corners.
Definition: intersectingentities.hh:67
const std::vector< GlobalPosition > & corners() const
Get the corners of the intersection geometry.
Definition: intersectingentities.hh:60
IntersectionInfo(std::size_t a, std::size_t b, Corners &&c)
Definition: intersectingentities.hh:45
Dune::FieldVector< ctype, dimworld > GlobalPosition
Definition: intersectingentities.hh:42
std::size_t second() const
Get the index of the intersecting entity belonging to the other grid.
Definition: intersectingentities.hh:56
std::size_t first() const
Get the index of the intersecting entity belonging to this grid.
Definition: intersectingentities.hh:52
typename Dune::PromotionTraits< CoordTypeA, CoordTypeB >::PromotedType ctype
Definition: intersectingentities.hh:40
A class for collision detection of two geometries and computation of intersection corners.
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:439
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:28
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:102
Detect if a point intersects a geometry.
Define some often used mathematical functions.
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:104
Definition: adapt.hh:17
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:294
Functionality to triangulate point clouds.