3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
geometry/intersectingentities.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
22#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
23#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
24
25#include <cmath>
26#include <type_traits>
27#include <vector>
28#include <algorithm>
29
30#include <dune/common/fvector.hh>
31
32#include <dumux/common/math.hh>
37
38namespace Dumux {
39
45template<int dimworld, class CoordTypeA, class CoordTypeB = CoordTypeA>
47{
48public:
49 using ctype = typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
50 static constexpr int dimensionworld = dimworld;
51 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
52
53 template<class Corners>
54 explicit IntersectionInfo(std::size_t a, std::size_t b, Corners&& c)
55 : a_(a)
56 , b_(b)
57 , corners_(c.begin(), c.end())
58 {}
59
61 std::size_t first() const
62 { return a_; }
63
65 std::size_t second() const
66 { return b_; }
67
69 std::vector<GlobalPosition> corners() const
70 { return corners_; }
71
76 bool cornersMatch(const std::vector<GlobalPosition>& otherCorners) const
77 {
78 if (otherCorners.size() != corners_.size())
79 return false;
80
81 const auto eps = 1.5e-7*(corners_[1] - corners_[0]).two_norm();
82 for (int i = 0; i < corners_.size(); ++i)
83 if ((corners_[i] - otherCorners[i]).two_norm() > eps)
84 return false;
85
86 return true;
87 }
88
89private:
90 std::size_t a_, b_;
91 std::vector<GlobalPosition> corners_;
92};
93
98template<class EntitySet, class ctype, int dimworld>
99inline std::vector<std::size_t>
100intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
101 const BoundingBoxTree<EntitySet>& tree,
102 bool isCartesianGrid = false)
103{
104 // Call the recursive find function to find candidates
105 std::vector<std::size_t> entities;
106 intersectingEntities(point, tree, tree.numBoundingBoxes() - 1, entities, isCartesianGrid);
107 return entities;
108}
109
114template<class EntitySet, class ctype, int dimworld>
115void intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
116 const BoundingBoxTree<EntitySet>& tree,
117 std::size_t node,
118 std::vector<std::size_t>& entities,
119 bool isCartesianGrid = false)
120{
121 // Get the bounding box for the current node
122 const auto& bBox = tree.getBoundingBoxNode(node);
123
124 // if the point is not in the bounding box we can stop
125 if (!intersectsPointBoundingBox(point, tree.getBoundingBoxCoordinates(node))) return;
126
127 // now we know it's inside
128 // if the box is a leaf do the primitive test.
129 else if (tree.isLeaf(bBox, node))
130 {
131 const std::size_t entityIdx = bBox.child1;
132 // for structured cube grids skip the primitive test
133 if (isCartesianGrid)
134 entities.push_back(entityIdx);
135 else
136 {
137 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
138 // if the primitive is positive it intersects the actual geometry, add the entity to the list
139 if (intersectsPointGeometry(point, geometry))
140 entities.push_back(entityIdx);
141 }
142 }
143
144 // No leaf. Check both children nodes.
145 else
146 {
147 intersectingEntities(point, tree, bBox.child0, entities, isCartesianGrid);
148 intersectingEntities(point, tree, bBox.child1, entities, isCartesianGrid);
149 }
150}
151
156template<class Geometry, class EntitySet>
157inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
158intersectingEntities(const Geometry& geometry,
159 const BoundingBoxTree<EntitySet>& tree)
160{
161 // check if the world dimensions match
162 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
163 "Can only intersect geometry and bounding box tree of same world dimension");
164
165 // Create data structure for return type
166 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
168 static constexpr int dimworld = Geometry::coorddimension;
169
170 // compute the bounding box of the given geometry
171 std::array<ctype, 2*Geometry::coorddimension> bBox;
172 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
173
174 // Get coordinates of first vertex
175 auto corner = geometry.corner(0);
176 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
177 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
178
179 // Compute the min and max over the remaining vertices
180 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
181 {
182 corner = geometry.corner(cornerIdx);
183 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
184 {
185 using std::max;
186 using std::min;
187 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
188 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
189 }
190 }
191
192 // Call the recursive find function to find candidates
193 intersectingEntities(geometry, tree,
194 bBox, tree.numBoundingBoxes() - 1,
195 intersections);
196
197 return intersections;
198}
199
204template<class Geometry, class EntitySet>
205void intersectingEntities(const Geometry& geometry,
206 const BoundingBoxTree<EntitySet>& tree,
207 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
208 std::size_t nodeIdx,
209 std::vector<IntersectionInfo<Geometry::coorddimension,
210 typename Geometry::ctype,
211 typename EntitySet::ctype>>& intersections)
212{
213 // if the two bounding boxes don't intersect we can stop searching
214 static constexpr int dimworld = Geometry::coorddimension;
215 if (!intersectsBoundingBoxBoundingBox<dimworld>(bBox.data(), tree.getBoundingBoxCoordinates(nodeIdx)))
216 return;
217
218 // get node info for current bounding box node
219 const auto& bBoxNode = tree.getBoundingBoxNode(nodeIdx);
220
221 // if the box is a leaf do the primitive test.
222 if (tree.isLeaf(bBoxNode, nodeIdx))
223 {
224 // eIdxA is always 0 since we intersect with exactly one geometry
225 const auto eIdxA = 0;
226 const auto eIdxB = bBoxNode.child1;
227
228 const auto geometryTree = tree.entitySet().entity(eIdxB).geometry();
229 using GeometryTree = std::decay_t<decltype(geometryTree)>;
231 using IntersectionAlgorithm = GeometryIntersection<Geometry, GeometryTree, Policy>;
232 using Intersection = typename IntersectionAlgorithm::Intersection;
233 Intersection intersection;
234
235 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
236 {
237 static constexpr int dimIntersection = Policy::dimIntersection;
238 if (dimIntersection >= 2)
239 {
240 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
241 for (unsigned int i = 0; i < triangulation.size(); ++i)
242 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
243 }
244 else
245 intersections.emplace_back(eIdxA, eIdxB, intersection);
246 }
247 }
248
249 // No leaf. Check both children nodes.
250 else
251 {
252 intersectingEntities(geometry, tree, bBox, bBoxNode.child0, intersections);
253 intersectingEntities(geometry, tree, bBox, bBoxNode.child1, intersections);
254 }
255}
256
261template<class EntitySet0, class EntitySet1>
262inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
264 const BoundingBoxTree<EntitySet1>& treeB)
265{
266 // check if the world dimensions match
267 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
268 "Can only intersect bounding box trees of same world dimension");
269
270 // Create data structure for return type
271 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
272
273 // Call the recursive find function to find candidates
274 intersectingEntities(treeA, treeB,
275 treeA.numBoundingBoxes() - 1,
276 treeB.numBoundingBoxes() - 1,
277 intersections);
278
279 return intersections;
280}
281
286template<class EntitySet0, class EntitySet1>
288 const BoundingBoxTree<EntitySet1>& treeB,
289 std::size_t nodeA, std::size_t nodeB,
290 std::vector<IntersectionInfo<EntitySet0::dimensionworld,
291 typename EntitySet0::ctype,
292 typename EntitySet1::ctype>>& intersections)
293{
294 // Get the bounding box for the current node
295 const auto& bBoxA = treeA.getBoundingBoxNode(nodeA);
296 const auto& bBoxB = treeB.getBoundingBoxNode(nodeB);
297
298 // if the two bounding boxes of the current nodes don't intersect we can stop searching
299 static constexpr int dimworld = EntitySet0::dimensionworld;
300 if (!intersectsBoundingBoxBoundingBox<dimworld>(treeA.getBoundingBoxCoordinates(nodeA),
301 treeB.getBoundingBoxCoordinates(nodeB)))
302 return;
303
304 // Check if we have a leaf in treeA or treeB
305 const bool isLeafA = treeA.isLeaf(bBoxA, nodeA);
306 const bool isLeafB = treeB.isLeaf(bBoxB, nodeB);
307
308 // If both boxes are leaves do the primitive test
309 if (isLeafA && isLeafB)
310 {
311 const auto eIdxA = bBoxA.child1;
312 const auto eIdxB = bBoxB.child1;
313
314 const auto geometryA = treeA.entitySet().entity(eIdxA).geometry();
315 const auto geometryB = treeB.entitySet().entity(eIdxB).geometry();
316
317 using GeometryA = std::decay_t<decltype(geometryA)>;
318 using GeometryB = std::decay_t<decltype(geometryB)>;
320 using IntersectionAlgorithm = GeometryIntersection<GeometryA, GeometryB, Policy>;
321 using Intersection = typename IntersectionAlgorithm::Intersection;
322 Intersection intersection;
323
324 if (IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
325 {
326 static constexpr int dimIntersection = Policy::dimIntersection;
327
328 if (dimIntersection >= 2)
329 {
330 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
331 for (unsigned int i = 0; i < triangulation.size(); ++i)
332 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
333 }
334 else
335 intersections.emplace_back(eIdxA, eIdxB, intersection);
336 }
337 }
338
339 // if we reached the leaf in treeA, just continue in treeB
340 else if (isLeafA)
341 {
342 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
343 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
344 }
345
346 // if we reached the leaf in treeB, just continue in treeA
347 else if (isLeafB)
348 {
349 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
350 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
351 }
352
353 // we know now that both trees didn't reach the leaf yet so
354 // we continue with the larger tree first (bigger node number)
355 else if (nodeA > nodeB)
356 {
357 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
358 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
359 }
360 else
361 {
362 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
363 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
364 }
365}
366
375template<class ctype, int dimworld>
376inline std::size_t intersectingEntityCartesianGrid(const Dune::FieldVector<ctype, dimworld>& point,
377 const Dune::FieldVector<ctype, dimworld>& min,
378 const Dune::FieldVector<ctype, dimworld>& max,
379 const std::array<int, std::size_t(dimworld)>& cells)
380{
381 std::size_t index = 0;
382 for (int i = 0; i < dimworld; ++i)
383 {
384 using std::clamp; using std::floor;
385 ctype dimOffset = clamp<ctype>(floor((point[i]-min[i])*cells[i]/(max[i]-min[i])), 0.0, cells[i]-1);
386 for (int j = 0; j < i; ++j)
387 dimOffset *= cells[j];
388 index += static_cast<std::size_t>(dimOffset);
389 }
390 return index;
391}
392
393} // end namespace Dumux
394
395#endif
Define some often used mathematical functions.
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: geometry/intersectingentities.hh:376
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: geometry/intersectspointgeometry.hh:38
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: geometry/intersectingentities.hh:100
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: geometry/boundingboxtree.hh:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometry/geometryintersection.hh:113
An axis-aligned bounding box volume tree implementation.
Definition: geometry/boundingboxtree.hh:66
const ctype * getBoundingBoxCoordinates(std::size_t nodeIdx) const
Get an existing bounding box for a given node.
Definition: geometry/boundingboxtree.hh:145
bool isLeaf(const BoundingBoxNode &node, std::size_t nodeIdx) const
Definition: geometry/boundingboxtree.hh:154
const EntitySet & entitySet() const
the entity set this tree was built with
Definition: geometry/boundingboxtree.hh:133
std::size_t numBoundingBoxes() const
Get the number of bounding boxes currently in the tree.
Definition: geometry/boundingboxtree.hh:149
const BoundingBoxNode & getBoundingBoxNode(std::size_t nodeIdx) const
Interface to be used by other bounding box trees.
Definition: geometry/boundingboxtree.hh:141
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometry/geometryintersection.hh:216
An intersection object resulting from the intersection of two primitives in an entity set.
Definition: geometry/intersectingentities.hh:47
static constexpr int dimensionworld
Definition: geometry/intersectingentities.hh:50
bool cornersMatch(const std::vector< GlobalPosition > &otherCorners) const
Check if the corners of this intersection match with the given corners.
Definition: geometry/intersectingentities.hh:76
IntersectionInfo(std::size_t a, std::size_t b, Corners &&c)
Definition: geometry/intersectingentities.hh:54
Dune::FieldVector< ctype, dimworld > GlobalPosition
Definition: geometry/intersectingentities.hh:51
std::vector< GlobalPosition > corners() const
Get the corners of the intersection geometry.
Definition: geometry/intersectingentities.hh:69
std::size_t second() const
Get the index of the intersecting entity belonging to the other grid.
Definition: geometry/intersectingentities.hh:65
std::size_t first() const
Get the index of the intersecting entity belonging to this grid.
Definition: geometry/intersectingentities.hh:61
typename Dune::PromotionTraits< CoordTypeA, CoordTypeB >::PromotedType ctype
Definition: geometry/intersectingentities.hh:49
An axis-aligned bounding box volume hierarchy for dune grids.
A class for collision detection of two geometries and computation of intersection corners.
Detect if a point intersects a geometry.
Functionality to triangulate point clouds.