3.2-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_INTERSECTING_ENTITIES_HH
23#define DUMUX_INTERSECTING_ENTITIES_HH
24
25#include <cmath>
26#include <type_traits>
27#include <vector>
28
29#include <dune/common/fvector.hh>
30
31#include <dumux/common/math.hh>
36
37namespace Dumux {
38
44template<int dimworld, class CoordTypeA, class CoordTypeB = CoordTypeA>
46{
47public:
48 using ctype = typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
49 static constexpr int dimensionworld = dimworld;
50 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
51
52 template<class Corners>
53 explicit IntersectionInfo(std::size_t a, std::size_t b, Corners&& c)
54 : a_(a)
55 , b_(b)
56 , corners_(c.begin(), c.end())
57 {}
58
60 std::size_t first() const
61 { return a_; }
62
64 std::size_t second() const
65 { return b_; }
66
68 std::vector<GlobalPosition> corners() const
69 { return corners_; }
70
75 bool cornersMatch(const std::vector<GlobalPosition>& otherCorners) const
76 {
77 if (otherCorners.size() != corners_.size())
78 return false;
79
80 const auto eps = 1.5e-7*(corners_[1] - corners_[0]).two_norm();
81 for (int i = 0; i < corners_.size(); ++i)
82 if ((corners_[i] - otherCorners[i]).two_norm() > eps)
83 return false;
84
85 return true;
86 }
87
88private:
89 std::size_t a_, b_;
90 std::vector<GlobalPosition> corners_;
91};
92
97template<class EntitySet, class ctype, int dimworld>
98inline std::vector<std::size_t>
99intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
100 const BoundingBoxTree<EntitySet>& tree,
101 bool isCartesianGrid = false)
102{
103 // Call the recursive find function to find candidates
104 std::vector<std::size_t> entities;
105 intersectingEntities(point, tree, tree.numBoundingBoxes() - 1, entities, isCartesianGrid);
106 return entities;
107}
108
113template<class EntitySet, class ctype, int dimworld>
114void intersectingEntities(const Dune::FieldVector<ctype, dimworld>& point,
115 const BoundingBoxTree<EntitySet>& tree,
116 std::size_t node,
117 std::vector<std::size_t>& entities,
118 bool isCartesianGrid = false)
119{
120 // Get the bounding box for the current node
121 const auto& bBox = tree.getBoundingBoxNode(node);
122
123 // if the point is not in the bounding box we can stop
124 if (!intersectsPointBoundingBox(point, tree.getBoundingBoxCoordinates(node))) return;
125
126 // now we know it's inside
127 // if the box is a leaf do the primitive test.
128 else if (tree.isLeaf(bBox, node))
129 {
130 const std::size_t entityIdx = bBox.child1;
131 // for structured cube grids skip the primitive test
132 if (isCartesianGrid)
133 entities.push_back(entityIdx);
134 else
135 {
136 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
137 // if the primitive is positive it intersects the actual geometry, add the entity to the list
138 if (intersectsPointGeometry(point, geometry))
139 entities.push_back(entityIdx);
140 }
141 }
142
143 // No leaf. Check both children nodes.
144 else
145 {
146 intersectingEntities(point, tree, bBox.child0, entities, isCartesianGrid);
147 intersectingEntities(point, tree, bBox.child1, entities, isCartesianGrid);
148 }
149}
150
155template<class Geometry, class EntitySet>
156inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
157intersectingEntities(const Geometry& geometry,
158 const BoundingBoxTree<EntitySet>& tree)
159{
160 // check if the world dimensions match
161 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
162 "Can only intersect geometry and bounding box tree of same world dimension");
163
164 // Create data structure for return type
165 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
167 static constexpr int dimworld = Geometry::coorddimension;
168
169 // compute the bounding box of the given geometry
170 std::array<ctype, 2*Geometry::coorddimension> bBox;
171 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
172
173 // Get coordinates of first vertex
174 auto corner = geometry.corner(0);
175 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
176 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
177
178 // Compute the min and max over the remaining vertices
179 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
180 {
181 corner = geometry.corner(cornerIdx);
182 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
183 {
184 using std::max;
185 using std::min;
186 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
187 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
188 }
189 }
190
191 // Call the recursive find function to find candidates
192 intersectingEntities(geometry, tree,
193 bBox, tree.numBoundingBoxes() - 1,
194 intersections);
195
196 return intersections;
197}
198
203template<class Geometry, class EntitySet>
204void intersectingEntities(const Geometry& geometry,
205 const BoundingBoxTree<EntitySet>& tree,
206 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
207 std::size_t nodeIdx,
208 std::vector<IntersectionInfo<Geometry::coorddimension,
209 typename Geometry::ctype,
210 typename EntitySet::ctype>>& intersections)
211{
212 // if the two bounding boxes don't intersect we can stop searching
213 static constexpr int dimworld = Geometry::coorddimension;
214 if (!intersectsBoundingBoxBoundingBox<dimworld>(bBox.data(), tree.getBoundingBoxCoordinates(nodeIdx)))
215 return;
216
217 // get node info for current bounding box node
218 const auto& bBoxNode = tree.getBoundingBoxNode(nodeIdx);
219
220 // if the box is a leaf do the primitive test.
221 if (tree.isLeaf(bBoxNode, nodeIdx))
222 {
223 // eIdxA is always 0 since we intersect with exactly one geometry
224 const auto eIdxA = 0;
225 const auto eIdxB = bBoxNode.child1;
226
227 const auto geometryTree = tree.entitySet().entity(eIdxB).geometry();
228 using GeometryTree = std::decay_t<decltype(geometryTree)>;
230 using IntersectionAlgorithm = GeometryIntersection<Geometry, GeometryTree, Policy>;
231 using Intersection = typename IntersectionAlgorithm::Intersection;
232 Intersection intersection;
233
234 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
235 {
236 static constexpr int dimIntersection = Policy::dimIntersection;
237 if (dimIntersection >= 2)
238 {
239 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
240 for (unsigned int i = 0; i < triangulation.size(); ++i)
241 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
242 }
243 else
244 intersections.emplace_back(eIdxA, eIdxB, intersection);
245 }
246 }
247
248 // No leaf. Check both children nodes.
249 else
250 {
251 intersectingEntities(geometry, tree, bBox, bBoxNode.child0, intersections);
252 intersectingEntities(geometry, tree, bBox, bBoxNode.child1, intersections);
253 }
254}
255
260template<class EntitySet0, class EntitySet1>
261inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
263 const BoundingBoxTree<EntitySet1>& treeB)
264{
265 // check if the world dimensions match
266 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
267 "Can only intersect bounding box trees of same world dimension");
268
269 // Create data structure for return type
270 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
271
272 // Call the recursive find function to find candidates
273 intersectingEntities(treeA, treeB,
274 treeA.numBoundingBoxes() - 1,
275 treeB.numBoundingBoxes() - 1,
276 intersections);
277
278 return intersections;
279}
280
285template<class EntitySet0, class EntitySet1>
287 const BoundingBoxTree<EntitySet1>& treeB,
288 std::size_t nodeA, std::size_t nodeB,
289 std::vector<IntersectionInfo<EntitySet0::dimensionworld,
290 typename EntitySet0::ctype,
291 typename EntitySet1::ctype>>& intersections)
292{
293 // Get the bounding box for the current node
294 const auto& bBoxA = treeA.getBoundingBoxNode(nodeA);
295 const auto& bBoxB = treeB.getBoundingBoxNode(nodeB);
296
297 // if the two bounding boxes of the current nodes don't intersect we can stop searching
298 static constexpr int dimworld = EntitySet0::dimensionworld;
299 if (!intersectsBoundingBoxBoundingBox<dimworld>(treeA.getBoundingBoxCoordinates(nodeA),
300 treeB.getBoundingBoxCoordinates(nodeB)))
301 return;
302
303 // Check if we have a leaf in treeA or treeB
304 const bool isLeafA = treeA.isLeaf(bBoxA, nodeA);
305 const bool isLeafB = treeB.isLeaf(bBoxB, nodeB);
306
307 // If both boxes are leaves do the primitive test
308 if (isLeafA && isLeafB)
309 {
310 const auto eIdxA = bBoxA.child1;
311 const auto eIdxB = bBoxB.child1;
312
313 const auto geometryA = treeA.entitySet().entity(eIdxA).geometry();
314 const auto geometryB = treeB.entitySet().entity(eIdxB).geometry();
315
316 using GeometryA = std::decay_t<decltype(geometryA)>;
317 using GeometryB = std::decay_t<decltype(geometryB)>;
319 using IntersectionAlgorithm = GeometryIntersection<GeometryA, GeometryB, Policy>;
320 using Intersection = typename IntersectionAlgorithm::Intersection;
321 Intersection intersection;
322
323 if (IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
324 {
325 static constexpr int dimIntersection = Policy::dimIntersection;
326
327 if (dimIntersection >= 2)
328 {
329 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
330 for (unsigned int i = 0; i < triangulation.size(); ++i)
331 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
332 }
333 else
334 intersections.emplace_back(eIdxA, eIdxB, intersection);
335 }
336 }
337
338 // if we reached the leaf in treeA, just continue in treeB
339 else if (isLeafA)
340 {
341 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
342 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
343 }
344
345 // if we reached the leaf in treeB, just continue in treeA
346 else if (isLeafB)
347 {
348 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
349 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
350 }
351
352 // we know now that both trees didn't reach the leaf yet so
353 // we continue with the larger tree first (bigger node number)
354 else if (nodeA > nodeB)
355 {
356 intersectingEntities(treeA, treeB, bBoxA.child0, nodeB, intersections);
357 intersectingEntities(treeA, treeB, bBoxA.child1, nodeB, intersections);
358 }
359 else
360 {
361 intersectingEntities(treeA, treeB, nodeA, bBoxB.child0, intersections);
362 intersectingEntities(treeA, treeB, nodeA, bBoxB.child1, intersections);
363 }
364}
365
366} // end namespace Dumux
367
368#endif
An axis-aligned bounding box volume hierarchy for dune grids.
Detect if a point intersects a geometry.
Functionality to triangulate point clouds.
A class for collision detection of two geometries and computation of intersection corners.
Define some often used mathematical functions.
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:99
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:114
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:217
An intersection object resulting from the intersection of two primitives in an entity set.
Definition: intersectingentities.hh:46
static constexpr int dimensionworld
Definition: intersectingentities.hh:49
bool cornersMatch(const std::vector< GlobalPosition > &otherCorners) const
Check if the corners of this intersection match with the given corners.
Definition: intersectingentities.hh:75
IntersectionInfo(std::size_t a, std::size_t b, Corners &&c)
Definition: intersectingentities.hh:53
Dune::FieldVector< ctype, dimworld > GlobalPosition
Definition: intersectingentities.hh:50
std::vector< GlobalPosition > corners() const
Get the corners of the intersection geometry.
Definition: intersectingentities.hh:68
std::size_t second() const
Get the index of the intersecting entity belonging to the other grid.
Definition: intersectingentities.hh:64
std::size_t first() const
Get the index of the intersecting entity belonging to this grid.
Definition: intersectingentities.hh:60
typename Dune::PromotionTraits< CoordTypeA, CoordTypeB >::PromotedType ctype
Definition: intersectingentities.hh:48