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