24#ifndef DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
25#define DUMUX_GEOMETRY_INTERSECTING_ENTITIES_HH
33#include <dune/common/fvector.hh>
48template<
int dimworld,
class CoordTypeA,
class CoordTypeB = CoordTypeA>
52 using ctype =
typename Dune::PromotionTraits<CoordTypeA, CoordTypeB>::PromotedType;
56 template<
class Corners>
60 , corners_(c.begin(), c.end())
72 const std::vector<GlobalPosition>&
corners()
const
79 bool cornersMatch(
const std::vector<GlobalPosition>& otherCorners)
const
81 if (otherCorners.size() != corners_.size())
85 ctype eps2 = std::numeric_limits<ctype>::min();
86 for (
int i = 1; i < corners_.size(); ++i)
87 eps2 = max(eps2, (corners_[i] - corners_[0]).two_norm2());
91 eps2 *= 1.5e-7*1.5e-7;
93 for (
int i = 0; i < corners_.size(); ++i)
95 if (std::none_of(otherCorners.begin(),
97 [&] (
const auto& other) { return (corners_[i] - other).two_norm2() < eps2; }))
105 std::vector<GlobalPosition> corners_;
112template<
class EntitySet,
class ctype,
int dimworld>
113inline std::vector<std::size_t>
116 bool isCartesianGrid =
false)
119 std::vector<std::size_t> entities;
128template<
class EntitySet,
class ctype,
int dimworld>
132 std::vector<std::size_t>& entities,
133 bool isCartesianGrid =
false)
143 else if (tree.
isLeaf(bBox, node))
145 const std::size_t entityIdx = bBox.child1;
148 entities.push_back(entityIdx);
151 const auto geometry = tree.
entitySet().entity(entityIdx).geometry();
154 entities.push_back(entityIdx);
170template<
class Geometry,
class EntitySet>
171inline std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>>
176 static_assert(int(Geometry::coorddimension) == int(EntitySet::dimensionworld),
177 "Can only intersect geometry and bounding box tree of same world dimension");
180 std::vector<IntersectionInfo<Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype>> intersections;
182 static constexpr int dimworld = Geometry::coorddimension;
185 std::array<ctype, 2*Geometry::coorddimension> bBox;
186 ctype* xMin = bBox.data(); ctype* xMax = xMin + Geometry::coorddimension;
189 auto corner = geometry.corner(0);
190 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
191 xMin[dimIdx] = xMax[dimIdx] = corner[dimIdx];
194 for (std::size_t cornerIdx = 1; cornerIdx < geometry.corners(); ++cornerIdx)
196 corner = geometry.corner(cornerIdx);
197 for (std::size_t dimIdx = 0; dimIdx < dimworld; ++dimIdx)
201 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
202 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
211 return intersections;
218template<
class Geometry,
class EntitySet>
221 const std::array<typename Geometry::ctype, 2*Geometry::coorddimension>& bBox,
224 typename Geometry::ctype,
225 typename EntitySet::ctype>>& intersections)
228 static constexpr int dimworld = Geometry::coorddimension;
236 if (tree.
isLeaf(bBoxNode, nodeIdx))
239 const auto eIdxA = 0;
240 const auto eIdxB = bBoxNode.child1;
242 const auto geometryTree = tree.
entitySet().entity(eIdxB).geometry();
243 using GeometryTree = std::decay_t<
decltype(geometryTree)>;
246 using Intersection =
typename IntersectionAlgorithm::Intersection;
247 Intersection intersection;
249 if (IntersectionAlgorithm::intersection(geometry, geometryTree, intersection))
251 static constexpr int dimIntersection = Policy::dimIntersection;
252 if (dimIntersection >= 2)
254 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
255 for (
unsigned int i = 0; i < triangulation.size(); ++i)
256 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
259 intersections.emplace_back(eIdxA, eIdxB, intersection);
275template<
class EntitySet0,
class EntitySet1>
276inline std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>>
281 static_assert(int(EntitySet0::dimensionworld) == int(EntitySet1::dimensionworld),
282 "Can only intersect bounding box trees of same world dimension");
285 std::vector<IntersectionInfo<EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype>> intersections;
293 return intersections;
300template<
class EntitySet0,
class EntitySet1>
303 std::size_t nodeA, std::size_t nodeB,
305 typename EntitySet0::ctype,
306 typename EntitySet1::ctype>>& intersections)
313 static constexpr int dimworld = EntitySet0::dimensionworld;
319 const bool isLeafA = treeA.
isLeaf(bBoxA, nodeA);
320 const bool isLeafB = treeB.
isLeaf(bBoxB, nodeB);
323 if (isLeafA && isLeafB)
325 const auto eIdxA = bBoxA.child1;
326 const auto eIdxB = bBoxB.child1;
328 const auto geometryA = treeA.
entitySet().entity(eIdxA).geometry();
329 const auto geometryB = treeB.
entitySet().entity(eIdxB).geometry();
331 using GeometryA = std::decay_t<
decltype(geometryA)>;
332 using GeometryB = std::decay_t<
decltype(geometryB)>;
335 using Intersection =
typename IntersectionAlgorithm::Intersection;
337 if (Intersection intersection; IntersectionAlgorithm::intersection(geometryA, geometryB, intersection))
339 static constexpr int dimIntersection = Policy::dimIntersection;
343 if (dimIntersection >= 2)
345 const auto triangulation = triangulate<dimIntersection, dimworld>(intersection);
346 for (
unsigned int i = 0; i < triangulation.size(); ++i)
347 intersections.emplace_back(eIdxA, eIdxB, std::move(triangulation[i]));
350 intersections.emplace_back(eIdxA, eIdxB, intersection);
370 else if (nodeA > nodeB)
390template<
class ctype,
int dimworld>
392 const Dune::FieldVector<ctype, dimworld>& min,
393 const Dune::FieldVector<ctype, dimworld>& max,
394 const std::array<
int, std::size_t(dimworld)>& cells)
396 std::size_t index = 0;
397 for (
int i = 0; i < dimworld; ++i)
399 using std::clamp;
using std::floor;
400 ctype dimOffset = clamp<ctype>(floor((point[i]-min[i])*cells[i]/(max[i]-min[i])), 0.0, cells[i]-1);
401 for (
int j = 0; j < i; ++j)
402 dimOffset *= cells[j];
403 index +=
static_cast<std::size_t
>(dimOffset);
Define some often used mathematical functions.
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.
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:391
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:40
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:114
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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:306
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:116
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:68
const ctype * getBoundingBoxCoordinates(std::size_t nodeIdx) const
Get an existing bounding box for a given node.
Definition: boundingboxtree.hh:147
bool isLeaf(const BoundingBoxNode &node, std::size_t nodeIdx) const
Definition: boundingboxtree.hh:156
const EntitySet & entitySet() const
the entity set this tree was built with
Definition: boundingboxtree.hh:135
std::size_t numBoundingBoxes() const
Get the number of bounding boxes currently in the tree.
Definition: boundingboxtree.hh:151
const BoundingBoxNode & getBoundingBoxNode(std::size_t nodeIdx) const
Interface to be used by other bounding box trees.
Definition: boundingboxtree.hh:143
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:219
An intersection object resulting from the intersection of two primitives in an entity set.
Definition: intersectingentities.hh:50
static constexpr int dimensionworld
Definition: intersectingentities.hh:53
bool cornersMatch(const std::vector< GlobalPosition > &otherCorners) const
Check if the corners of this intersection match with the given corners.
Definition: intersectingentities.hh:79
const std::vector< GlobalPosition > & corners() const
Get the corners of the intersection geometry.
Definition: intersectingentities.hh:72
IntersectionInfo(std::size_t a, std::size_t b, Corners &&c)
Definition: intersectingentities.hh:57
Dune::FieldVector< ctype, dimworld > GlobalPosition
Definition: intersectingentities.hh:54
std::size_t second() const
Get the index of the intersecting entity belonging to the other grid.
Definition: intersectingentities.hh:68
std::size_t first() const
Get the index of the intersecting entity belonging to this grid.
Definition: intersectingentities.hh:64
typename Dune::PromotionTraits< CoordTypeA, CoordTypeB >::PromotedType ctype
Definition: intersectingentities.hh:52