version 3.8
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-FileCopyrightInfo: 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
296namespace Detail {
297
298// helper struct to compute distance between two geometries, specialized below
299template<class Geo1, class Geo2,
300 int dimWorld = Geo1::coorddimension,
301 int dim1 = Geo1::mydimension, int dim2 = Geo2::mydimension>
303{
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)
306 {
307 DUNE_THROW(Dune::NotImplemented, "Geometry distance computation not implemented for dimworld = "
308 << dimWorld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
309 }
310};
311
312// distance point-point
313template<class Geo1, class Geo2, int dimWorld>
314struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 0>
315{
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)
318 { return Dumux::squaredDistance(geo1.corner(0), geo2.corner(0)); }
319};
320
321// distance segment-point
322template<class Geo1, class Geo2, int dimWorld>
323struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 1, 0>
324{
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)
327 { return squaredDistancePointSegment(geo2.corner(0), geo1); }
328};
329
330// distance point-segment
331template<class Geo1, class Geo2, int dimWorld>
332struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 1>
333{
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)
336 { return squaredDistancePointSegment(geo1.corner(0), geo2); }
337};
338
339// distance point-polygon
340template<class Geo1, class Geo2, int dimWorld>
341struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 2>
342{
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)
345 { return squaredDistancePointPolygon(geo1.corner(0), geo2); }
346};
347
348// distance polygon-point
349template<class Geo1, class Geo2, int dimWorld>
350struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 2, 0>
351{
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)
354 { return squaredDistancePointPolygon(geo2.corner(0), geo1); }
355};
356
357} // end namespace Detail
358
363template<class Geo1, class Geo2>
364static inline auto squaredDistance(const Geo1& geo1, const Geo2& geo2)
366
371template<class Geo1, class Geo2>
372static inline auto distance(const Geo1& geo1, const Geo2& geo2)
373{ using std::sqrt; return sqrt(squaredDistance(geo1, geo2)); }
374
375
377namespace Detail {
378
384template<class EntitySet, class ctype, int dimworld,
385 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0), int> = 0>
386void closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
387 const BoundingBoxTree<EntitySet>& tree,
388 std::size_t node,
389 ctype& minSquaredDistance,
390 std::size_t& eIdx)
391{
392 // Get the bounding box for the current node
393 const auto& bBox = tree.getBoundingBoxNode(node);
394
395 // If bounding box is outside radius, then don't search further
397 point, tree.getBoundingBoxCoordinates(node)
398 );
399
400 // do not continue if the AABB is further away than the current minimum
401 if (squaredDistance > minSquaredDistance) return;
402
403 // If the bounding box is a leaf, update distance and index with primitive test
404 else if (tree.isLeaf(bBox, node))
405 {
406 const std::size_t entityIdx = bBox.child1;
407 const auto squaredDistance = [&]{
408 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
409 if constexpr (EntitySet::Entity::Geometry::mydimension == 2)
410 return squaredDistancePointPolygon(point, geometry);
411 else if constexpr (EntitySet::Entity::Geometry::mydimension == 1)
412 return squaredDistancePointSegment(point, geometry);
413 else
414 DUNE_THROW(Dune::NotImplemented, "squaredDistance to entity with dim>2");
415 }();
416
417 if (squaredDistance < minSquaredDistance)
418 {
419 eIdx = entityIdx;
420 minSquaredDistance = squaredDistance;
421 }
422 }
423
424 // Check the children nodes recursively
425 else
426 {
427 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
428 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
429 }
430}
431
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,
440 const BoundingBoxTree<EntitySet>& tree,
441 std::size_t node,
442 ctype& minSquaredDistance,
443 std::size_t& eIdx)
444{
445 // Get the bounding box for the current node
446 const auto& bBox = tree.getBoundingBoxNode(node);
447
448 // If the bounding box is a leaf, update distance and index with primitive test
449 if (tree.isLeaf(bBox, node))
450 {
451 const std::size_t entityIdx = bBox.child1;
452 const auto& p = tree.entitySet().entity(entityIdx).geometry().corner(0);
453 const auto squaredDistance = (p-point).two_norm2();
454
455 if (squaredDistance < minSquaredDistance)
456 {
457 eIdx = entityIdx;
458 minSquaredDistance = squaredDistance;
459 }
460 }
461
462 // Check the children nodes recursively
463 else
464 {
465 // If bounding box is outside radius, then don't search further
467 point, tree.getBoundingBoxCoordinates(node)
468 );
469
470 // do not continue if the AABB is further away than the current minimum
471 if (squaredDistance > minSquaredDistance) return;
472
473 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
474 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
475 }
476}
477
478} // end namespace Detail
479
493template<class EntitySet, class ctype, int dimworld>
494std::pair<ctype, std::size_t> closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
495 const BoundingBoxTree<EntitySet>& tree,
496 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
497{
498 std::size_t eIdx = 0;
499 Detail::closestEntity(point, tree, tree.numBoundingBoxes() - 1, minSquaredDistance, eIdx);
500 using std::sqrt;
501 return { minSquaredDistance, eIdx };
502}
503
508template<class EntitySet, class ctype, int dimworld>
509ctype squaredDistance(const Dune::FieldVector<ctype, dimworld>& point,
510 const BoundingBoxTree<EntitySet>& tree,
511 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
512{
513 return closestEntity(point, tree, minSquaredDistance).first;
514}
515
520template<class EntitySet, class ctype, int dimworld>
521ctype distance(const Dune::FieldVector<ctype, dimworld>& point,
522 const BoundingBoxTree<EntitySet>& tree,
523 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
524{
525 using std::sqrt;
526 return sqrt(squaredDistance(point, tree, minSquaredDistance));
527}
528
529} // end namespace Dumux
530
531#endif
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.
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: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