3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_DISTANCE_HH
25#define DUMUX_GEOMETRY_DISTANCE_HH
26
27#include <dune/common/fvector.hh>
28#include <dune/geometry/quadraturerules.hh>
29
30#include <dumux/common/math.hh>
32
33namespace Dumux {
34
39template<class Geometry>
40static inline typename Geometry::ctype
41averageDistancePointGeometry(const typename Geometry::GlobalCoordinate& p,
42 const Geometry& geometry,
43 std::size_t integrationOrder = 2)
44{
45 typename Geometry::ctype avgDist = 0.0;
46 const auto& quad = Dune::QuadratureRules<typename Geometry::ctype, Geometry::mydimension>::rule(geometry.type(), integrationOrder);
47 for (const auto& qp : quad)
48 avgDist += (geometry.global(qp.position())-p).two_norm()*qp.weight()*geometry.integrationElement(qp.position());
49 return avgDist/geometry.volume();
50}
51
56template<class Point>
57static inline typename Point::value_type
58squaredDistancePointLine(const Point& p, const Point& a, const Point& b)
59{
60 const auto ab = b - a;
61 const auto t = (p - a)*ab/ab.two_norm2();
62 auto proj = a;
63 proj.axpy(t, ab);
64 return (proj - p).two_norm2();
65}
66
73template<class Geometry>
74static inline typename Geometry::ctype
75squaredDistancePointLine(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
76{
77 static_assert(Geometry::mydimension == 1, "Geometry has to be a line");
78 const auto& a = geometry.corner(0);
79 const auto& b = geometry.corner(1);
80 return squaredDistancePointLine(p, a, b);
81}
82
87template<class Point>
88static inline typename Point::value_type
89distancePointLine(const Point& p, const Point& a, const Point& b)
90{ using std::sqrt; return sqrt(squaredDistancePointLine(p, a, b)); }
91
98template<class Geometry>
99static inline typename Geometry::ctype
100distancePointLine(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
101{ using std::sqrt; return sqrt(squaredDistancePointLine(p, geometry)); }
102
107template<class Point>
108static inline typename Point::value_type
109squaredDistancePointSegment(const Point& p, const Point& a, const Point& b)
110{
111 const auto ab = b - a;
112 const auto ap = p - a;
113 const auto t = ap*ab;
114
115 if (t <= 0.0)
116 return ap.two_norm2();
117
118 const auto lengthSq = ab.two_norm2();
119 if (t >= lengthSq)
120 return (b - p).two_norm2();
121
122 auto proj = a;
123 proj.axpy(t/lengthSq, ab);
124 return (proj - p).two_norm2();
125}
126
131template<class Geometry>
132static inline typename Geometry::ctype
133squaredDistancePointSegment(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
134{
135 static_assert(Geometry::mydimension == 1, "Geometry has to be a segment");
136 const auto& a = geometry.corner(0);
137 const auto& b = geometry.corner(1);
138 return squaredDistancePointSegment(p, a, b);
139}
140
145template<class Point>
146static inline typename Point::value_type
147distancePointSegment(const Point& p, const Point& a, const Point& b)
148{ using std::sqrt; return sqrt(squaredDistancePointSegment(p, a, b)); }
149
154template<class Geometry>
155static inline typename Geometry::ctype
156distancePointSegment(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
157{ using std::sqrt; return sqrt(squaredDistancePointSegment(p, geometry)); }
158
164template<class Point>
165static inline typename Point::value_type
166squaredDistancePointTriangle(const Point& p, const Point& a, const Point& b, const Point& c)
167{
168 static_assert(Point::dimension == 3, "Only works in 3D");
169 const auto ab = b - a;
170 const auto bc = c - b;
171 const auto ca = a - c;
172 const auto normal = crossProduct(ab, ca);
173
174 const auto ap = p - a;
175 const auto bp = p - b;
176 const auto cp = p - c;
177
178 const auto sum = sign(crossProduct(ab, normal)*ap)
179 + sign(crossProduct(bc, normal)*bp)
180 + sign(crossProduct(ca, normal)*cp);
181
182 // if there is no orthogonal projection
183 // (point is outside the infinite prism implied by the triangle and its surface normal)
184 // compute distance to the edges (codim-1 facets)
185 if (sum < 2.0)
186 {
187 using std::min;
188 return min({squaredDistancePointSegment(p, a, b),
191 }
192 // compute distance via orthogonal projection
193 else
194 {
195 const auto tmp = normal*ap;
196 return tmp*tmp / (normal*normal);
197 }
198}
199
204template<class Geometry>
205static inline typename Geometry::ctype
206squaredDistancePointTriangle(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
207{
208 static_assert(Geometry::coorddimension == 3, "Only works in 3D");
209 static_assert(Geometry::mydimension == 2, "Geometry has to be a triangle");
210 assert(geometry.corners() == 3);
211 const auto& a = geometry.corner(0);
212 const auto& b = geometry.corner(1);
213 const auto& c = geometry.corner(2);
214 return squaredDistancePointTriangle(p, a, b, c);
215}
216
221template<class Point>
222static inline typename Point::value_type
223distancePointTriangle(const Point& p, const Point& a, const Point& b, const Point& c)
224{ using std::sqrt; return sqrt(squaredDistancePointTriangle(p, a, b, c)); }
225
230template<class Geometry>
231static inline typename Geometry::ctype
232distancePointTriangle(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
233{ using std::sqrt; return sqrt(squaredDistancePointTriangle(p, geometry)); }
234
240template<class Geometry>
241static inline typename Geometry::ctype
242squaredDistancePointPolygon(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
243{
244 if (geometry.corners() == 3)
245 return squaredDistancePointTriangle(p, geometry);
246 else if (geometry.corners() == 4)
247 {
248 const auto& a = geometry.corner(0);
249 const auto& b = geometry.corner(1);
250 const auto& c = geometry.corner(2);
251 const auto& d = geometry.corner(3);
252
253 using std::min;
254 return min(squaredDistancePointTriangle(p, a, b, d),
255 squaredDistancePointTriangle(p, a, d, c));
256 }
257 else
258 DUNE_THROW(Dune::NotImplemented, "Polygon with " << geometry.corners() << " corners not supported");
259}
260
266template<class Geometry>
267static inline typename Geometry::ctype
268distancePointPolygon(const typename Geometry::GlobalCoordinate& p, const Geometry& geometry)
269{ using std::sqrt; return sqrt(squaredDistancePointPolygon(p, geometry)); }
270
275template<class Geometry>
276static inline typename Geometry::ctype
277averageDistanceSegmentGeometry(const typename Geometry::GlobalCoordinate& a,
278 const typename Geometry::GlobalCoordinate& b,
279 const Geometry& geometry,
280 std::size_t integrationOrder = 2)
281{
282 typename Geometry::ctype avgDist = 0.0;
283 const auto& quad = Dune::QuadratureRules<typename Geometry::ctype, Geometry::mydimension>::rule(geometry.type(), integrationOrder);
284 for (const auto& qp : quad)
285 avgDist += distancePointSegment(geometry.global(qp.position()), a, b)*qp.weight()*geometry.integrationElement(qp.position());
286 return avgDist/geometry.volume();
287}
288
293template<class ctype, int dimWorld>
294static inline ctype distance(const Dune::FieldVector<ctype, dimWorld>& a,
295 const Dune::FieldVector<ctype, dimWorld>& b)
296{ return (a-b).two_norm(); }
297
302template<class ctype, int dimWorld>
303static inline ctype squaredDistance(const Dune::FieldVector<ctype, dimWorld>& a,
304 const Dune::FieldVector<ctype, dimWorld>& b)
305{ return (a-b).two_norm2(); }
306
307
308namespace Detail {
309
310// helper struct to compute distance between two geometries, specialized below
311template<class Geo1, class Geo2,
312 int dimWorld = Geo1::coorddimension,
313 int dim1 = Geo1::mydimension, int dim2 = Geo2::mydimension>
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 {
319 DUNE_THROW(Dune::NotImplemented, "Geometry distance computation not implemented for dimworld = "
320 << dimWorld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
321 }
322};
323
324// distance point-point
325template<class Geo1, class Geo2, int dimWorld>
326struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 0>
327{
328 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
329 static auto distance(const Geo1& geo1, const Geo2& geo2)
330 { return Dumux::squaredDistance(geo1.corner(0), geo2.corner(0)); }
331};
332
333// distance segment-point
334template<class Geo1, class Geo2, int dimWorld>
335struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 1, 0>
336{
337 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
338 static auto distance(const Geo1& geo1, const Geo2& geo2)
339 { return squaredDistancePointSegment(geo2.corner(0), geo1); }
340};
341
342// distance point-segment
343template<class Geo1, class Geo2, int dimWorld>
344struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 1>
345{
346 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
347 static auto distance(const Geo1& geo1, const Geo2& geo2)
348 { return squaredDistancePointSegment(geo1.corner(0), geo2); }
349};
350
351// distance point-polygon
352template<class Geo1, class Geo2, int dimWorld>
353struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 0, 2>
354{
355 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
356 static inline auto distance(const Geo1& geo1, const Geo2& geo2)
357 { return squaredDistancePointPolygon(geo1.corner(0), geo2); }
358};
359
360// distance polygon-point
361template<class Geo1, class Geo2, int dimWorld>
362struct GeometrySquaredDistance<Geo1, Geo2, dimWorld, 2, 0>
363{
364 static_assert(Geo1::coorddimension == Geo2::coorddimension, "Geometries have to have the same coordinate dimensions");
365 static inline auto distance(const Geo1& geo1, const Geo2& geo2)
366 { return squaredDistancePointPolygon(geo2.corner(0), geo1); }
367};
368
369} // end namespace Detail
370
375template<class Geo1, class Geo2>
376static inline auto squaredDistance(const Geo1& geo1, const Geo2& geo2)
378
383template<class Geo1, class Geo2>
384static inline auto distance(const Geo1& geo1, const Geo2& geo2)
385{ using std::sqrt; return sqrt(squaredDistance(geo1, geo2)); }
386
387
389namespace Detail {
390
396template<class EntitySet, class ctype, int dimworld,
397 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0), int> = 0>
398void closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
399 const BoundingBoxTree<EntitySet>& tree,
400 std::size_t node,
401 ctype& minSquaredDistance,
402 std::size_t& eIdx)
403{
404 // Get the bounding box for the current node
405 const auto& bBox = tree.getBoundingBoxNode(node);
406
407 // If bounding box is outside radius, then don't search further
409 point, tree.getBoundingBoxCoordinates(node)
410 );
411
412 // do not continue if the AABB is further away than the current minimum
413 if (squaredDistance > minSquaredDistance) return;
414
415 // If the bounding box is a leaf, update distance and index with primitive test
416 else if (tree.isLeaf(bBox, node))
417 {
418 const std::size_t entityIdx = bBox.child1;
419 const auto squaredDistance = [&]{
420 const auto geometry = tree.entitySet().entity(entityIdx).geometry();
421 if constexpr (EntitySet::Entity::Geometry::mydimension == 2)
422 return squaredDistancePointPolygon(point, geometry);
423 else if constexpr (EntitySet::Entity::Geometry::mydimension == 1)
424 return squaredDistancePointSegment(point, geometry);
425 else
426 DUNE_THROW(Dune::NotImplemented, "squaredDistance to entity with dim>2");
427 }();
428
429 if (squaredDistance < minSquaredDistance)
430 {
431 eIdx = entityIdx;
432 minSquaredDistance = squaredDistance;
433 }
434 }
435
436 // Check the children nodes recursively
437 else
438 {
439 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
440 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
441 }
442}
443
449template<class EntitySet, class ctype, int dimworld,
450 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension == 0), int> = 0>
451void closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
452 const BoundingBoxTree<EntitySet>& tree,
453 std::size_t node,
454 ctype& minSquaredDistance,
455 std::size_t& eIdx)
456{
457 // Get the bounding box for the current node
458 const auto& bBox = tree.getBoundingBoxNode(node);
459
460 // If the bounding box is a leaf, update distance and index with primitive test
461 if (tree.isLeaf(bBox, node))
462 {
463 const std::size_t entityIdx = bBox.child1;
464 const auto& p = tree.entitySet().entity(entityIdx).geometry().corner(0);
465 const auto squaredDistance = (p-point).two_norm2();
466
467 if (squaredDistance < minSquaredDistance)
468 {
469 eIdx = entityIdx;
470 minSquaredDistance = squaredDistance;
471 }
472 }
473
474 // Check the children nodes recursively
475 else
476 {
477 // If bounding box is outside radius, then don't search further
479 point, tree.getBoundingBoxCoordinates(node)
480 );
481
482 // do not continue if the AABB is further away than the current minimum
483 if (squaredDistance > minSquaredDistance) return;
484
485 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
486 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
487 }
488}
489
490} // end namespace Detail
491
505template<class EntitySet, class ctype, int dimworld>
506std::pair<ctype, std::size_t> closestEntity(const Dune::FieldVector<ctype, dimworld>& point,
507 const BoundingBoxTree<EntitySet>& tree,
508 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
509{
510 std::size_t eIdx = 0;
511 Detail::closestEntity(point, tree, tree.numBoundingBoxes() - 1, minSquaredDistance, eIdx);
512 using std::sqrt;
513 return { minSquaredDistance, eIdx };
514}
515
520template<class EntitySet, class ctype, int dimworld>
521ctype squaredDistance(const Dune::FieldVector<ctype, dimworld>& point,
522 const BoundingBoxTree<EntitySet>& tree,
523 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
524{
525 return closestEntity(point, tree, minSquaredDistance).first;
526}
527
532template<class EntitySet, class ctype, int dimworld>
533ctype distance(const Dune::FieldVector<ctype, dimworld>& point,
534 const BoundingBoxTree<EntitySet>& tree,
535 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
536{
537 using std::sqrt;
538 return sqrt(squaredDistance(point, tree, minSquaredDistance));
539}
540
541} // end namespace Dumux
542
543#endif
Define some often used mathematical functions.
An axis-aligned bounding box volume hierarchy for dune grids.
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:294
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:38
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:268
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:147
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:166
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:58
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:223
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:242
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:109
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:303
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:506
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:89
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:277
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:41
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:398
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:654
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
ctype squaredDistancePointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Compute squared distance between point and bounding box.
Definition: boundingboxtree.hh:422
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
Definition: distance.hh:315
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:329
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:338
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:347
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:356
static auto distance(const Geo1 &geo1, const Geo2 &geo2)
Definition: distance.hh:365