24#ifndef DUMUX_GEOMETRY_DISTANCE_HH
25#define DUMUX_GEOMETRY_DISTANCE_HH
27#include <dune/common/fvector.hh>
28#include <dune/geometry/quadraturerules.hh>
39template<
class Geometry>
40static inline typename Geometry::ctype
42 const Geometry& geometry,
43 std::size_t integrationOrder = 2)
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();
57static inline typename Point::value_type
60 const auto ab = b - a;
61 const auto t = (p - a)*ab/ab.two_norm2();
64 return (proj - p).two_norm2();
73template<
class Geometry>
74static inline typename Geometry::ctype
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);
88static inline typename Point::value_type
98template<
class Geometry>
99static inline typename Geometry::ctype
108static inline typename Point::value_type
111 const auto ab = b - a;
112 const auto ap = p - a;
113 const auto t = ap*ab;
116 return ap.two_norm2();
118 const auto lengthSq = ab.two_norm2();
120 return (b - p).two_norm2();
123 proj.axpy(t/lengthSq, ab);
124 return (proj - p).two_norm2();
131template<
class Geometry>
132static inline typename Geometry::ctype
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);
146static inline typename Point::value_type
154template<
class Geometry>
155static inline typename Geometry::ctype
165static inline typename Point::value_type
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;
174 const auto ap = p - a;
175 const auto bp = p - b;
176 const auto cp = p - c;
195 const auto tmp =
normal*ap;
204template<
class Geometry>
205static inline typename Geometry::ctype
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);
222static inline typename Point::value_type
230template<
class Geometry>
231static inline typename Geometry::ctype
240template<
class Geometry>
241static inline typename Geometry::ctype
244 if (geometry.corners() == 3)
246 else if (geometry.corners() == 4)
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);
258 DUNE_THROW(Dune::NotImplemented,
"Polygon with " << geometry.corners() <<
" corners not supported");
266template<
class Geometry>
267static inline typename Geometry::ctype
275template<
class Geometry>
276static inline typename Geometry::ctype
278 const typename Geometry::GlobalCoordinate& b,
279 const Geometry& geometry,
280 std::size_t integrationOrder = 2)
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();
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(); }
302template<
class ctype,
int dimWorld>
304 const Dune::FieldVector<ctype, dimWorld>& b)
305{
return (a-b).two_norm2(); }
311template<
class Geo1,
class Geo2,
312 int dimWorld = Geo1::coorddimension,
313 int dim1 = Geo1::mydimension,
int dim2 = Geo2::mydimension>
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)
319 DUNE_THROW(Dune::NotImplemented,
"Geometry distance computation not implemented for dimworld = "
320 << dimWorld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
325template<
class Geo1,
class Geo2,
int dimWorld>
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)
334template<
class Geo1,
class Geo2,
int dimWorld>
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)
343template<
class Geo1,
class Geo2,
int dimWorld>
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)
352template<
class Geo1,
class Geo2,
int dimWorld>
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)
361template<
class Geo1,
class Geo2,
int dimWorld>
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)
375template<
class Geo1,
class Geo2>
383template<
class Geo1,
class Geo2>
384static inline auto distance(
const Geo1& geo1,
const Geo2& geo2)
396template<
class EntitySet,
class ctype,
int dimworld,
397 typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0),
int> = 0>
401 ctype& minSquaredDistance,
416 else if (tree.
isLeaf(bBox, node))
418 const std::size_t entityIdx = bBox.child1;
420 const auto geometry = tree.
entitySet().entity(entityIdx).geometry();
421 if constexpr (EntitySet::Entity::Geometry::mydimension == 2)
423 else if constexpr (EntitySet::Entity::Geometry::mydimension == 1)
426 DUNE_THROW(Dune::NotImplemented,
"squaredDistance to entity with dim>2");
439 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
440 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
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,
454 ctype& minSquaredDistance,
461 if (tree.
isLeaf(bBox, node))
463 const std::size_t entityIdx = bBox.child1;
464 const auto& p = tree.
entitySet().entity(entityIdx).geometry().corner(0);
485 closestEntity(point, tree, bBox.child0, minSquaredDistance, eIdx);
486 closestEntity(point, tree, bBox.child1, minSquaredDistance, eIdx);
505template<
class EntitySet,
class ctype,
int dimworld>
506std::pair<ctype, std::size_t>
closestEntity(
const Dune::FieldVector<ctype, dimworld>& point,
508 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
510 std::size_t eIdx = 0;
513 return { minSquaredDistance, eIdx };
520template<
class EntitySet,
class ctype,
int dimworld>
523 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
532template<
class EntitySet,
class ctype,
int dimworld>
533ctype
distance(
const Dune::FieldVector<ctype, dimworld>& point,
535 ctype minSquaredDistance = std::numeric_limits<ctype>::max())
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