22#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
23#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
31#include <dune/common/fvector.hh>
46template<
int dimworld,
class CoordTypeA,
class CoordTypeB = CoordTypeA>
50 using ctype =
typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
54 template<
class Corners>
58 , corners_(c.begin(), c.end())
70 const std::vector<GlobalPosition>&
corners()
const
77 bool cornersMatch(
const std::vector<GlobalPosition>& otherCorners)
const
79 if (otherCorners.size() != corners_.size())
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());
89 eps2 *= 1.5e-7*1.5e-7;
91 for (
int i = 0; i < corners_.size(); ++i)
93 if (std::none_of(otherCorners.begin(),
95 [&] (
const auto& other) { return (corners_[i] - other).two_norm2() < eps2; }))
103 std::vector<GlobalPosition> corners_;
110template<
class EntitySet,
class ctype,
int dimworld>
111inline std::vector<std::size_t>
114 bool isCartesianGrid =
false)
117 std::vector<std::size_t> entities;
126template<
class EntitySet,
class ctype,
int dimworld>
130 std::vector<std::size_t>& entities,
131 bool isCartesianGrid =
false)
141 else if (tree.
isLeaf(bBox, node))
143 const std::size_t entityIdx = bBox.child1;
146 entities.push_back(entityIdx);
149 const auto geometry = tree.
entitySet().entity(entityIdx).geometry();
152 entities.push_back(entityIdx);
168template<
class Geometry,
class EntitySet>
169inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
174 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
175 "Can only intersect geometry and bounding box tree of same world dimension");
178 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
180 static constexpr int dimworld = Geometry::coorddimension;
183 std::array<ctype, 2*Geometry::coorddimension> bBox;
184 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
187 auto corner = geometry.corner(0);
188 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
189 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
192 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
194 corner = geometry.corner(cornerIdx);
195 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
199 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
200 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
209 return intersections;
216template<
class Geometry,
class EntitySet>
219 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
222 typename Geometry::ctype,
223 typename EntitySet::ctype>>& intersections)
226 static constexpr int dimworld = Geometry::coorddimension;
234 if (tree.
isLeaf(bBoxNode, nodeIdx))
237 const auto eIdxA = 0;
238 const auto eIdxB = bBoxNode.child1;
240 const auto geometryTree = tree.
entitySet().entity(eIdxB).geometry();
241 using GeometryTree = std::decay_t<
decltype(geometryTree)>;
244 using Intersection =
typename IntersectionAlgorithm::Intersection;
245 Intersection intersection;
247 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
249 static constexpr int dimIntersection = Policy::dimIntersection;
250 if (dimIntersection >= 2)
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]));
257 intersections.emplace_back(eIdxA, eIdxB, intersection);
273template<
class EntitySet0,
class EntitySet1>
274inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
279 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
280 "Can only intersect bounding box trees of same world dimension");
283 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
291 return intersections;
298template<
class EntitySet0,
class EntitySet1>
301 std::size_t nodeA, std::size_t nodeB,
303 typename EntitySet0::ctype,
304 typename EntitySet1::ctype>>& intersections)
311 static constexpr int dimworld = EntitySet0::dimensionworld;
317 const bool isLeafA = treeA.
isLeaf(bBoxA, nodeA);
318 const bool isLeafB = treeB.
isLeaf(bBoxB, nodeB);
321 if (isLeafA && isLeafB)
323 const auto eIdxA = bBoxA.child1;
324 const auto eIdxB = bBoxB.child1;
326 const auto geometryA = treeA.
entitySet().entity(eIdxA).geometry();
327 const auto geometryB = treeB.
entitySet().entity(eIdxB).geometry();
329 using GeometryA = std::decay_t<
decltype(geometryA)>;
330 using GeometryB = std::decay_t<
decltype(geometryB)>;
333 using Intersection =
typename IntersectionAlgorithm::Intersection;
335 if (Intersection intersection; IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
337 static constexpr int dimIntersection = Policy::dimIntersection;
341 if (dimIntersection >= 2)
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]));
348 intersections.emplace_back(eIdxA, eIdxB, intersection);
368 else if (nodeA > nodeB)
388template<
class ctype,
int dimworld>
390 const Dune::FieldVector<ctype, dimworld>& min,
391 const Dune::FieldVector<ctype, dimworld>& max,
392 const std::array<
int, std::size_t(dimworld)>& cells)
394 std::size_t index = 0;
395 for (
int i = 0; i < dimworld; ++i)
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);
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
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