12#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
13#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
21#include <dune/common/fvector.hh>
36template<
int dimworld,
class CoordTypeA,
class CoordTypeB = CoordTypeA>
40 using ctype =
typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
44 template<
class Corners>
48 , corners_(c.begin(), c.end())
60 const std::vector<GlobalPosition>&
corners()
const
67 bool cornersMatch(
const std::vector<GlobalPosition>& otherCorners)
const
69 if (otherCorners.size() != corners_.size())
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());
79 eps2 *= 1.5e-7*1.5e-7;
81 for (
int i = 0; i < corners_.size(); ++i)
83 if (std::none_of(otherCorners.begin(),
85 [&] (
const auto& other) { return (corners_[i] - other).two_norm2() < eps2; }))
93 std::vector<GlobalPosition> corners_;
100template<
class EntitySet,
class ctype,
int dimworld>
101inline std::vector<std::size_t>
104 bool isCartesianGrid =
false)
107 std::vector<std::size_t> entities;
116template<
class EntitySet,
class ctype,
int dimworld>
120 std::vector<std::size_t>& entities,
121 bool isCartesianGrid =
false)
131 else if (tree.
isLeaf(bBox, node))
133 const std::size_t entityIdx = bBox.child1;
136 entities.push_back(entityIdx);
139 const auto geometry = tree.
entitySet().entity(entityIdx).geometry();
142 entities.push_back(entityIdx);
158template<
class Geometry,
class EntitySet>
159inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
171template<
class Geometry,
class EntitySet,
class IntersectionPolicy>
172inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
175 IntersectionPolicy intersectionPolicy)
178 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
179 "Can only intersect geometry and bounding box tree of same world dimension");
182 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
184 static constexpr int dimworld = Geometry::coorddimension;
187 std::array<ctype, 2*Geometry::coorddimension> bBox;
188 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
191 auto corner = geometry.corner(0);
192 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
193 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
196 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
198 corner = geometry.corner(cornerIdx);
199 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
203 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
204 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
211 intersections, intersectionPolicy);
213 return intersections;
220template<
class Geometry,
class EntitySet>
223 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
226 typename Geometry::ctype,
227 typename EntitySet::ctype>>& intersections)
236template<
class Geometry,
class EntitySet,
class IntersectionPolicy>
239 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
242 typename Geometry::ctype,
243 typename EntitySet::ctype>>& intersections,
244 IntersectionPolicy intersectionPolicy)
247 static constexpr int dimworld = Geometry::coorddimension;
255 if (tree.
isLeaf(bBoxNode, nodeIdx))
258 const auto eIdxA = 0;
259 const auto eIdxB = bBoxNode.child1;
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;
267 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
269 static constexpr int dimIntersection = IntersectionPolicy::dimIntersection;
270 if constexpr (dimIntersection >= 2)
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]));
277 intersections.emplace_back(eIdxA, eIdxB, intersection);
293template<
class EntitySet0,
class EntitySet1>
294inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
306template<
class EntitySet0,
class EntitySet1,
class IntersectionPolicy>
307inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
310 IntersectionPolicy intersectionPolicy)
313 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
314 "Can only intersect bounding box trees of same world dimension");
317 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
323 intersections, intersectionPolicy);
325 return intersections;
332template<
class EntitySet0,
class EntitySet1>
335 std::size_t nodeA, std::size_t nodeB,
337 typename EntitySet0::ctype,
338 typename EntitySet1::ctype>>& intersections)
348template<
class EntitySet0,
class EntitySet1,
class IntersectionPolicy>
351 std::size_t nodeA, std::size_t nodeB,
353 typename EntitySet0::ctype,
354 typename EntitySet1::ctype>>& intersections,
355 IntersectionPolicy intersectionPolicy)
362 static constexpr int dimworld = EntitySet0::dimensionworld;
368 const bool isLeafA = treeA.
isLeaf(bBoxA, nodeA);
369 const bool isLeafB = treeB.
isLeaf(bBoxB, nodeB);
372 if (isLeafA && isLeafB)
374 const auto eIdxA = bBoxA.child1;
375 const auto eIdxB = bBoxB.child1;
377 const auto geometryA = treeA.
entitySet().entity(eIdxA).geometry();
378 const auto geometryB = treeB.
entitySet().entity(eIdxB).geometry();
380 using GeometryA = std::decay_t<
decltype(geometryA)>;
381 using GeometryB = std::decay_t<
decltype(geometryB)>;
383 using Intersection =
typename IntersectionAlgorithm::Intersection;
385 if (Intersection intersection; IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
387 static constexpr int dimIntersection = IntersectionPolicy::dimIntersection;
391 if constexpr (dimIntersection >= 2)
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]));
398 intersections.emplace_back(eIdxA, eIdxB, intersection);
418 else if (nodeA > nodeB)
438template<
class ctype,
int dimworld>
440 const Dune::FieldVector<ctype, dimworld>& min,
441 const Dune::FieldVector<ctype, dimworld>& max,
442 const std::array<
int, std::size_t(dimworld)>& cells)
444 std::size_t index = 0;
445 for (
int i = 0; i < dimworld; ++i)
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);
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
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.