12#ifndef DUMUX_GEOMETRY_DISTANCE_HH
13#define DUMUX_GEOMETRY_DISTANCE_HH
15#include <dune/common/fvector.hh>
16#include <dune/geometry/quadraturerules.hh>
27template<
class Geometry>
28static inline typename Geometry::ctype
30 const Geometry& geometry,
31 std::size_t integrationOrder = 2)
33 typename Geometry::ctype avgDist = 0.0;
34 const auto& quad = Dune::QuadratureRules<typename Geometry::ctype, Geometry::mydimension>::rule(geometry.type(), integrationOrder);
35 for (
const auto& qp : quad)
36 avgDist += (geometry.global(qp.position())-p).two_norm()*qp.weight()*geometry.integrationElement(qp.position());
37 return avgDist/geometry.volume();
45static inline typename Point::value_type
48 const auto ab = b - a;
49 const auto t = (p - a)*ab/ab.two_norm2();
52 return (proj - p).two_norm2();
61template<
class Geometry>
62static inline typename Geometry::ctype
65 static_assert(Geometry::mydimension == 1,
"Geometry has to be a line");
66 const auto& a = geometry.corner(0);
67 const auto& b = geometry.corner(1);
76static inline typename Point::value_type
86template<
class Geometry>
87static inline typename Geometry::ctype
96static inline typename Point::value_type
99 const auto ab = b - a;
100 const auto ap = p - a;
101 const auto t = ap*ab;
104 return ap.two_norm2();
106 const auto lengthSq = ab.two_norm2();
108 return (b - p).two_norm2();
111 proj.axpy(t/lengthSq, ab);
112 return (proj - p).two_norm2();
119template<
class Geometry>
120static inline typename Geometry::ctype
123 static_assert(Geometry::mydimension == 1,
"Geometry has to be a segment");
124 const auto& a = geometry.corner(0);
125 const auto& b = geometry.corner(1);
134static inline typename Point::value_type
142template<
class Geometry>
143static inline typename Geometry::ctype
153static inline typename Point::value_type
156 static_assert(Point::dimension == 3,
"Only works in 3D");
157 const auto ab = b - a;
158 const auto bc = c - b;
159 const auto ca = a - c;
162 const auto ap = p - a;
163 const auto bp = p - b;
164 const auto cp = p - c;
183 const auto tmp =
normal*ap;
192template<
class Geometry>
193static inline typename Geometry::ctype
196 static_assert(Geometry::coorddimension == 3,
"Only works in 3D");
197 static_assert(Geometry::mydimension == 2,
"Geometry has to be a triangle");
198 assert(geometry.corners() == 3);
199 const auto& a = geometry.corner(0);
200 const auto& b = geometry.corner(1);
201 const auto& c = geometry.corner(2);
210static inline typename Point::value_type
218template<
class Geometry>
219static inline typename Geometry::ctype
228template<
class Geometry>
229static inline typename Geometry::ctype
232 if (geometry.corners() == 3)
234 else if (geometry.corners() == 4)
236 const auto& a = geometry.corner(0);
237 const auto& b = geometry.corner(1);
238 const auto& c = geometry.corner(2);
239 const auto& d = geometry.corner(3);
246 DUNE_THROW(Dune::NotImplemented,
"Polygon with " << geometry.corners() <<
" corners not supported");
254template<
class Geometry>
255static inline typename Geometry::ctype
263template<
class Geometry>
264static inline typename Geometry::ctype
266 const typename Geometry::GlobalCoordinate& b,
267 const Geometry& geometry,
268 std::size_t integrationOrder = 2)
270 typename Geometry::ctype avgDist = 0.0;
271 const auto& quad = Dune::QuadratureRules<typename Geometry::ctype, Geometry::mydimension>::rule(geometry.type(), integrationOrder);
272 for (
const auto& qp : quad)
273 avgDist +=
distancePointSegment(geometry.global(qp.position()), a, b)*qp.weight()*geometry.integrationElement(qp.position());
274 return avgDist/geometry.volume();
281template<
class ctype,
int dimWorld>
282static inline ctype
distance(
const Dune::FieldVector<ctype, dimWorld>& a,
283 const Dune::FieldVector<ctype, dimWorld>& b)
284{
return (a-b).two_norm(); }
290template<
class ctype,
int dimWorld>
292 const Dune::FieldVector<ctype, dimWorld>& b)
293{
return (a-b).two_norm2(); }
299template<
class Geo1,
class Geo2,
300 int dimWorld = Geo1::coorddimension,
301 int dim1 = Geo1::mydimension,
int dim2 = Geo2::mydimension>
304 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
305 static auto distance(
const Geo1& geo1,
const Geo2& geo2)
307 DUNE_THROW(Dune::NotImplemented,
"Geometry distance computation not implemented for dimworld = "
308 << dimWorld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
313template<
class Geo1,
class Geo2,
int dimWorld>
316 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
317 static auto distance(
const Geo1& geo1,
const Geo2& geo2)
322template<
class Geo1,
class Geo2,
int dimWorld>
325 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
326 static auto distance(
const Geo1& geo1,
const Geo2& geo2)
331template<
class Geo1,
class Geo2,
int dimWorld>
334 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
335 static auto distance(
const Geo1& geo1,
const Geo2& geo2)
340template<
class Geo1,
class Geo2,
int dimWorld>
343 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
344 static inline auto distance(
const Geo1& geo1,
const Geo2& geo2)
349template<
class Geo1,
class Geo2,
int dimWorld>
352 static_assert(Geo1::coorddimension == Geo2::coorddimension,
"Geometries have to have the same coordinate dimensions");
353 static inline auto distance(
const Geo1& geo1,
const Geo2& geo2)
363template<
class Geo1,
class Geo2>
371template<
class Geo1,
class Geo2>
372static inline auto distance(
const Geo1& geo1,
const Geo2& geo2)
384template<
class EntitySet,
class ctype,
int dimworld,
385 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0),
int> = 0>
389 ctype& minSquaredDistance,
404 else if (tree.
isLeaf(bBox, node))
406 const std::size_t entityIdx = bBox.child1;
408 const auto geometry = tree.
entitySet().entity(entityIdx).geometry();
409 if constexpr (EntitySet::Entity::Geometry::mydimension == 2)
411 else if constexpr (EntitySet::Entity::Geometry::mydimension == 1)
414 DUNE_THROW(Dune::NotImplemented,
"squaredDistance to entity with dim>2");
427 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
428 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
437template<
class EntitySet,
class ctype,
int dimworld,
438 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension == 0),
int> = 0>
439void closestEntity(
const Dune::FieldVector<ctype, dimworld>& point,
442 ctype& minSquaredDistance,
449 if (tree.
isLeaf(bBox, node))
451 const std::size_t entityIdx = bBox.child1;
452 const auto& p = tree.
entitySet().entity(entityIdx).geometry().corner(0);
473 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
474 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
493template<
class EntitySet,
class ctype,
int dimworld>
494std::pair<ctype, std::size_t>
closestEntity(
const Dune::FieldVector<ctype, dimworld>& point,
496 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
498 std::size_t eIdx = 0;
501 return { minSquaredDistance, eIdx };
508template<
class EntitySet,
class ctype,
int dimworld>
511 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
520template<
class EntitySet,
class ctype,
int dimworld>
521ctype
distance(
const Dune::FieldVector<ctype, dimworld>& point,
523 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
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
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:642
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:629
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
static Geometry::ctype distancePointPolygon(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
Compute the shortest distance from a point to a given polygon geometry.
Definition: distance.hh:256
static Point::value_type distancePointSegment(const Point &p, const Point &a, const Point &b)
Compute the distance from a point to the segment connecting the points a and b.
Definition: distance.hh:135
static Point::value_type squaredDistancePointTriangle(const Point &p, const Point &a, const Point &b, const Point &c)
Compute the shortest squared distance from a point to the triangle connecting the points a,...
Definition: distance.hh:154
static Point::value_type squaredDistancePointLine(const Point &p, const Point &a, const Point &b)
Compute the squared distance from a point to a line through the points a and b.
Definition: distance.hh:46
static Point::value_type distancePointTriangle(const Point &p, const Point &a, const Point &b, const Point &c)
Compute the shortest distance from a point to the triangle connecting the points a,...
Definition: distance.hh:211
static Geometry::ctype squaredDistancePointPolygon(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
Compute the shortest squared distance from a point to a given polygon geometry.
Definition: distance.hh:230
static Point::value_type squaredDistancePointSegment(const Point &p, const Point &a, const Point &b)
Compute the squared distance from a point to the segment connecting the points a and b.
Definition: distance.hh:97
static ctype squaredDistance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest squared distance between two points.
Definition: distance.hh:291
std::pair< ctype, std::size_t > closestEntity(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, ctype minSquaredDistance=std::numeric_limits< ctype >::max())
Compute the closest entity in an AABB tree to a point (index and shortest squared distance)
Definition: distance.hh:494
static Point::value_type distancePointLine(const Point &p, const Point &a, const Point &b)
Compute the distance from a point to a line through the points a and b.
Definition: distance.hh:77
static Geometry::ctype averageDistanceSegmentGeometry(const typename Geometry::GlobalCoordinate &a, const typename Geometry::GlobalCoordinate &b, const Geometry &geometry, std::size_t integrationOrder=2)
Compute the average distance from a segment to a geometry by integration.
Definition: distance.hh:265
static Geometry::ctype averageDistancePointGeometry(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry, std::size_t integrationOrder=2)
Compute the average distance from a point to a geometry by integration.
Definition: distance.hh:29
void closestEntity(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, std::size_t node, ctype &minSquaredDistance, std::size_t &eIdx)
Compute the closest entity in an AABB tree (index and shortest squared distance) recursively.
Definition: distance.hh:386
Define some often used mathematical functions.
ctype squaredDistancePointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Compute squared distance between point and bounding box.
Definition: boundingboxtree.hh:410
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:317
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:335
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:344
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:326
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:353
Definition: distance.hh:303
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:305