version 3.11-dev
distance.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_GEOMETRY_DISTANCE_HH
13#define DUMUX_GEOMETRY_DISTANCE_HH
14
15#include <dune/common/fvector.hh>
16#include <dune/geometry/quadraturerules.hh>
17
18#include <dumux/common/math.hh>
20
21namespace Dumux {
22
27template<class Geometry>
28static inline typename Geometry::ctype
29averageDistancePointGeometry(const typename Geometry::GlobalCoordinate& p,
30 const Geometry& geometry,
31 std::size_t integrationOrder = 2)
32{
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();
38}
39
44template<class Point>
45static inline typename Point::value_type
46squaredDistancePointLine(const Point& p, const Point& a, const Point& b)
47{
48 const auto ab = b - a;
49 const auto t = (p - a)*ab/ab.two_norm2();
50 auto proj = a;
51 proj.axpy(t, ab);
52 return (proj - p).two_norm2();
53}
54
61template<class Geometry>
62static inline typename Geometry::ctype
63squaredDistancePointLine(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
64{
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);
68 return squaredDistancePointLine(p, a, b);
69}
70
75template<class Point>
76static inline typename Point::value_type
77distancePointLine(const Point& p, const Point& a, const Point& b)
78{ using std::sqrt; return sqrt(squaredDistancePointLine(p, a, b)); }
79
86template<class Geometry>
87static inline typename Geometry::ctype
88distancePointLine(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
89{ using std::sqrt; return sqrt(squaredDistancePointLine(p, geometry)); }
90
95template<class Point>
96static inline typename Point::value_type
97squaredDistancePointSegment(const Point& p, const Point& a, const Point& b)
98{
99 const auto ab = b - a;
100 const auto ap = p - a;
101 const auto t = ap*ab;
102
103 if (t <= 0.0)
104 return ap.two_norm2();
105
106 const auto lengthSq = ab.two_norm2();
107 if (t >= lengthSq)
108 return (b - p).two_norm2();
109
110 auto proj = a;
111 proj.axpy(t/lengthSq, ab);
112 return (proj - p).two_norm2();
113}
114
119template<class Geometry>
120static inline typename Geometry::ctype
121squaredDistancePointSegment(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
122{
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);
126 return squaredDistancePointSegment(p, a, b);
127}
128
133template<class Point>
134static inline typename Point::value_type
135distancePointSegment(const Point& p, const Point& a, const Point& b)
136{ using std::sqrt; return sqrt(squaredDistancePointSegment(p, a, b)); }
137
142template<class Geometry>
143static inline typename Geometry::ctype
144distancePointSegment(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
145{ using std::sqrt; return sqrt(squaredDistancePointSegment(p, geometry)); }
146
152template<class Point>
153static inline typename Point::value_type
154squaredDistancePointTriangle(const Point& p, const Point& a, const Point& b, const Point& c)
155{
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;
160 const auto normal = crossProduct(ab, ca);
161
162 const auto ap = p - a;
163 const auto bp = p - b;
164 const auto cp = p - c;
165
166 const auto sum = sign(crossProduct(ab, normal)*ap)
167 + sign(crossProduct(bc, normal)*bp)
168 + sign(crossProduct(ca, normal)*cp);
169
170 // if there is no orthogonal projection
171 // (point is outside the infinite prism implied by the triangle and its surface normal)
172 // compute distance to the edges (codim-1 facets)
173 if (sum < 2.0)
174 {
175 using std::min;
176 return min({squaredDistancePointSegment(p, a, b),
179 }
180 // compute distance via orthogonal projection
181 else
182 {
183 const auto tmp = normal*ap;
184 return tmp*tmp / (normal*normal);
185 }
186}
187
192template<class Geometry>
193static inline typename Geometry::ctype
194squaredDistancePointTriangle(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
195{
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);
202 return squaredDistancePointTriangle(p, a, b, c);
203}
204
209template<class Point>
210static inline typename Point::value_type
211distancePointTriangle(const Point& p, const Point& a, const Point& b, const Point& c)
212{ using std::sqrt; return sqrt(squaredDistancePointTriangle(p, a, b, c)); }
213
218template<class Geometry>
219static inline typename Geometry::ctype
220distancePointTriangle(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
221{ using std::sqrt; return sqrt(squaredDistancePointTriangle(p, geometry)); }
222
228template<class Geometry>
229static inline typename Geometry::ctype
230squaredDistancePointPolygon(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
231{
232 if (geometry.corners() == 3)
233 return squaredDistancePointTriangle(p, geometry);
234 else if (geometry.corners() == 4)
235 {
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);
240
241 using std::min;
242 return min(squaredDistancePointTriangle(p, a, b, d),
243 squaredDistancePointTriangle(p, a, d, c));
244 }
245 else
246 DUNE_THROW(Dune::NotImplemented, "Polygon with " << geometry.corners() << " corners not supported");
247}
248
254template<class Geometry>
255static inline typename Geometry::ctype
256distancePointPolygon(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
257{ using std::sqrt; return sqrt(squaredDistancePointPolygon(p, geometry)); }
258
263template<class Geometry>
264static inline typename Geometry::ctype
265averageDistanceSegmentGeometry(const typename Geometry::GlobalCoordinate& a,
266 const typename Geometry::GlobalCoordinate& b,
267 const Geometry& geometry,
268 std::size_t integrationOrder = 2)
269{
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();
275}
276
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(); }
285
290template<class ctype, int dimWorld>
291static inline ctype squaredDistance(const Dune::FieldVector<ctype, dimWorld>& a,
292 const Dune::FieldVector<ctype, dimWorld>& b)
293{ return (a-b).two_norm2(); }
294
295
296#ifndef DOXYGEN
297namespace Detail {
298
299// helper struct to compute distance between two geometries, specialized below
300template<class Geo1, class Geo2,
301 int dimWorld = Geo1::coorddimension,
302 int dim1 = Geo1::mydimension, int dim2 = Geo2::mydimension>
303struct GeometrySquaredDistance
304{
305 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
306 static auto distance(const Geo1& geo1, const Geo2& geo2)
307 {
308 DUNE_THROW(Dune::NotImplemented, "Geometry distance computation not implemented for dimworld = "
309 << dimWorld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
310 }
311};
312
313// distance point-point
314template<class Geo1, class Geo2, int dimWorld>
315struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 0>
316{
317 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
318 static auto distance(const Geo1& geo1, const Geo2& geo2)
319 { return Dumux::squaredDistance(geo1.corner(0), geo2.corner(0)); }
320};
321
322// distance segment-point
323template<class Geo1, class Geo2, int dimWorld>
324struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 1, 0>
325{
326 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
327 static auto distance(const Geo1& geo1, const Geo2& geo2)
328 { return squaredDistancePointSegment(geo2.corner(0), geo1); }
329};
330
331// distance point-segment
332template<class Geo1, class Geo2, int dimWorld>
333struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 1>
334{
335 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
336 static auto distance(const Geo1& geo1, const Geo2& geo2)
337 { return squaredDistancePointSegment(geo1.corner(0), geo2); }
338};
339
340// distance point-polygon
341template<class Geo1, class Geo2, int dimWorld>
342struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 2>
343{
344 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
345 static inline auto distance(const Geo1& geo1, const Geo2& geo2)
346 { return squaredDistancePointPolygon(geo1.corner(0), geo2); }
347};
348
349// distance polygon-point
350template<class Geo1, class Geo2, int dimWorld>
351struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 2, 0>
352{
353 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
354 static inline auto distance(const Geo1& geo1, const Geo2& geo2)
355 { return squaredDistancePointPolygon(geo2.corner(0), geo1); }
356};
357
358} // end namespace Detail
359#endif // DOXYGEN
360
365template<class Geo1, class Geo2>
366static inline auto squaredDistance(const Geo1& geo1, const Geo2& geo2)
368
373template<class Geo1, class Geo2>
374static inline auto distance(const Geo1& geo1, const Geo2& geo2)
375{ using std::sqrt; return sqrt(squaredDistance(geo1, geo2)); }
376
377
378#ifndef DOXYGEN
379namespace Detail {
380
386template<class EntitySet, class ctype, int dimworld,
387 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0), int> = 0>
388void closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
389 const BoundingBoxTree<EntitySet>& tree,
390 std::size_t node,
391 ctype& minSquaredDistance,
392 std::size_t& eIdx)
393{
394 // Get the bounding box for the current node
395 const auto& bBox = tree.getBoundingBoxNode(node);
396
397 // If bounding box is outside radius, then don't search further
399 point, tree.getBoundingBoxCoordinates(node)
400 );
401
402 // do not continue if the AABB is further away than the current minimum
403 if (squaredDistance > minSquaredDistance) return;
404
405 // If the bounding box is a leaf, update distance and index with primitive test
406 else if (tree.isLeaf(bBox, node))
407 {
408 const std::size_t entityIdx = bBox.child1;
409 const auto squaredDistance = [&]{
410 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
411 if constexpr (EntitySet::Entity::Geometry::mydimension == 2)
412 return squaredDistancePointPolygon(point, geometry);
413 else if constexpr (EntitySet::Entity::Geometry::mydimension == 1)
414 return squaredDistancePointSegment(point, geometry);
415 else
416 DUNE_THROW(Dune::NotImplemented, "squaredDistance to entity with dim>2");
417 }();
418
419 if (squaredDistance < minSquaredDistance)
420 {
421 eIdx = entityIdx;
422 minSquaredDistance = squaredDistance;
423 }
424 }
425
426 // Check the children nodes recursively
427 else
428 {
429 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
430 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
431 }
432}
433
439template<class EntitySet, class ctype, int dimworld,
440 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension == 0), int> = 0>
441void closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
442 const BoundingBoxTree<EntitySet>& tree,
443 std::size_t node,
444 ctype& minSquaredDistance,
445 std::size_t& eIdx)
446{
447 // Get the bounding box for the current node
448 const auto& bBox = tree.getBoundingBoxNode(node);
449
450 // If the bounding box is a leaf, update distance and index with primitive test
451 if (tree.isLeaf(bBox, node))
452 {
453 const std::size_t entityIdx = bBox.child1;
454 const auto& p = tree.entitySet().entity(entityIdx).geometry().corner(0);
455 const auto squaredDistance = (p-point).two_norm2();
456
457 if (squaredDistance < minSquaredDistance)
458 {
459 eIdx = entityIdx;
460 minSquaredDistance = squaredDistance;
461 }
462 }
463
464 // Check the children nodes recursively
465 else
466 {
467 // If bounding box is outside radius, then don't search further
469 point, tree.getBoundingBoxCoordinates(node)
470 );
471
472 // do not continue if the AABB is further away than the current minimum
473 if (squaredDistance > minSquaredDistance) return;
474
475 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
476 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
477 }
478}
479
480} // end namespace Detail
481#endif // DOXYGEN
482
496template<class EntitySet, class ctype, int dimworld>
497std::pair<ctype, std::size_t> closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
498 const BoundingBoxTree<EntitySet>& tree,
499 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
500{
501 std::size_t eIdx = 0;
502 Detail::closestEntity(point, tree, tree.numBoundingBoxes() - 1, minSquaredDistance, eIdx);
503 using std::sqrt;
504 return { minSquaredDistance, eIdx };
505}
506
511template<class EntitySet, class ctype, int dimworld>
512ctype squaredDistance(const Dune::FieldVector<ctype, dimworld>& point,
513 const BoundingBoxTree<EntitySet>& tree,
514 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
515{
516 return closestEntity(point, tree, minSquaredDistance).first;
517}
518
523template<class EntitySet, class ctype, int dimworld>
524ctype distance(const Dune::FieldVector<ctype, dimworld>& point,
525 const BoundingBoxTree<EntitySet>& tree,
526 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
527{
528 using std::sqrt;
529 return sqrt(squaredDistance(point, tree, minSquaredDistance));
530}
531
532} // end namespace Dumux
533
534#endif
An axis-aligned bounding box volume hierarchy for dune grids.
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:67
std::size_t numBoundingBoxes() const
Get the number of bounding boxes currently in the tree.
Definition: boundingboxtree.hh:150
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:671
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:658
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
ctype distance(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, ctype minSquaredDistance=std::numeric_limits< ctype >::max())
Compute the shortest distance to entities in an AABB tree.
Definition: distance.hh:524
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:497
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
Define some often used mathematical functions.
Definition: adapt.hh:17
ctype squaredDistancePointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Compute squared distance between point and bounding box.
Definition: boundingboxtree.hh:431