3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Files | Classes | Typedefs | Functions
Geometry

Algorithms for geometry computations (intersections, distances, ...). More...

Description

Algorithms for geometry computations (intersections, distances, ...).

Files

file  boundingboxtree.hh
 An axis-aligned bounding box volume hierarchy for dune grids.
 
file  boundingsphere.hh
 A function to compute bounding spheres of points clouds or convex polytopes.
 
file  center.hh
 Compute the center point of a convex polytope geometry or a random-access container of corner points.
 
file  circumsphere.hh
 A function to compute the circumsphere of a given geometry.
 
file  diameter.hh
 A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
 
file  distance.hh
 Helper functions for distance queries.
 
file  distancefield.hh
 
file  geometricentityset.hh
 An interface for a set of geometric entities.
 
file  geometryintersection.hh
 A class for collision detection of two geometries and computation of intersection corners.
 
file  grahamconvexhull.hh
 A function to compute the convex hull of a point cloud.
 
file  intersectingentities.hh
 Algorithms that finds which geometric entities intersect.
 
file  intersectionentityset.hh
 A class representing the intersection entities of two geometric entity sets.
 
file  intersectspointgeometry.hh
 Detect if a point intersects a geometry.
 
file  intersectspointsimplex.hh
 Detect if a point intersects a simplex (including boundary)
 
file  makegeometry.hh
 Create Dune geometries from user-specified points.
 
file  normal.hh
 Helper functions for normals.
 
file  refinementquadraturerule.hh
 A quadrature based on refinement.
 
file  sphere.hh
 A function to compute bounding spheres of points clouds or convex polytopes.
 
file  triangulation.hh
 Functionality to triangulate point clouds.
 
file  volume.hh
 Compute the volume of several common geometry types.
 

Classes

class  Dumux::BoundingBoxTree< GeometricEntitySet >
 An axis-aligned bounding box volume tree implementation. More...
 
class  Dumux::AABBDistanceField< Geometry >
 Class to calculate the closest distance from a point to a given set of geometries describing the domain's boundaries. Internally uses an AABB tree representation of the geometries for logarithmic distance queries. More...
 
class  Dumux::GridViewGeometricEntitySet< GridView, codim, Mapper >
 An interface for a set of geometric entities based on a GridView. More...
 
class  Dumux::GeometriesEntitySet< GeoType >
 An interface for a set of geometric entities. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, dimworld, dim1, dim2 >
 A class for geometry collision detection and intersection calculation The class can be specialized for combinations of dimworld, dim1, dim2, where dimworld is the world dimension embedding a grid of dim1 and a grid of dim2. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 2, 1, 1 >
 A class for segment–segment intersection in 2d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 2, 2, 1 >
 A class for polygon–segment intersection in 2d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 2, 1, 2 >
 A class for segment–polygon intersection in 2d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 2, 2, 2 >
 A class for polygon–polygon intersection in 2d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 3, 1 >
 A class for polyhedron–segment intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 1, 3 >
 A class for segment–polyhedron intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 2, 2 >
 A class for polygon–polygon intersections in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 3, 2 >
 A class for polyhedron–polygon intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 2, 3 >
 A class for polygon–polyhedron intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 3, 3 >
 A class for polyhedron–polyhedron intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 2, 1 >
 A class for polygon–segment intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 1, 2 >
 A class for segment–polygon intersection in 3d space. More...
 
class  Dumux::GeometryIntersection< Geometry1, Geometry2, Policy, 3, 1, 1 >
 A class for segment–segment intersection in 3d space. More...
 
class  Dumux::IntersectionInfo< dimworld, CoordTypeA, CoordTypeB >
 An intersection object resulting from the intersection of two primitives in an entity set. More...
 
class  Dumux::IntersectionEntitySet< DomainEntitySet, TargetEntitySet >
 A class representing the intersection entities two geometric entity sets. More...
 
class  Dumux::RefinementQuadratureRule< ct, mydim >
 A "quadrature" based on virtual refinement. More...
 
class  Dumux::Sphere< Scalar, dim >
 A simple sphere type. More...
 

Typedefs

template<class Geometry >
using Dumux::DistanceField = AABBDistanceField< Geometry >
 Class to calculate the closest distance from a point to a given set of geometries describing the domain's boundaries. More...
 
template<int dim, int dimWorld, class ctype >
using Dumux::Triangulation = std::vector< std::array< Dune::FieldVector< ctype, dimWorld >, dim+1 > >
 The default data type to store triangulations. More...
 

Functions

template<class ctype , int dimworld>
bool Dumux::intersectsPointBoundingBox (const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &min, const Dune::FieldVector< ctype, dimworld > &max)
 Determine if a point intersects an axis-aligned bounding box The bounding box is given by the lower left corner (min) and the upper right corner (max) More...
 
template<class ConvexGeometry >
static Sphere< typename ConvexGeometry::ctype, ConvexGeometry::coorddimension > Dumux::boundingSphere (const ConvexGeometry &geometry)
 Computes a bounding sphere of a convex polytope geometry (Dune::Geometry interface) More...
 
template<class Scalar , int dim>
static Sphere< Scalar, dim > Dumux::boundingSphere (const std::vector< Dune::FieldVector< Scalar, dim > > &points)
 Computes a bounding sphere of points. More...
 
template<class Corners >
Corners::value_type Dumux::center (const Corners &corners)
 The center of a given list of corners. More...
 
template<class Point >
static Sphere< typename Point::value_type, Point::dimension > Dumux::circumSphereOfTriangle (const Point &a, const Point &b, const Point &c)
 Computes the circumsphere of a triangle. More...
 
template<class Geometry , typename std::enable_if_t< Geometry::mydimension==2, int > = 0>
static Sphere< typename Geometry::ctype, Geometry::coorddimension > Dumux::circumSphereOfTriangle (const Geometry &triangle)
 Computes the circumsphere of a triangle. More...
 
template<class Geometry >
Geometry::ctype Dumux::diameter (const Geometry &geo)
 Computes the longest distance between points of a geometry. More...
 
template<class Geometry >
static Geometry::ctype Dumux::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. More...
 
template<class Point >
static Point::value_type Dumux::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. More...
 
template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointLine (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the squared distance from a point to a line given by a geometry with two corners. More...
 
template<class Point >
static Point::value_type Dumux::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. More...
 
template<class Geometry >
static Geometry::ctype Dumux::distancePointLine (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the distance from a point to a line given by a geometry with two corners. More...
 
template<class Point >
static Point::value_type Dumux::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. More...
 
template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointSegment (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the squared distance from a point to a given segment geometry. More...
 
template<class Point >
static Point::value_type Dumux::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. More...
 
template<class Geometry >
static Geometry::ctype Dumux::distancePointSegment (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the distance from a point to a given segment geometry. More...
 
template<class Point >
static Point::value_type Dumux::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, b and c See https://www.iquilezles.org/www/articles/triangledistance/triangledistance.htm. More...
 
template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointTriangle (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the shortest squared distance from a point to a given triangle geometry. More...
 
template<class Point >
static Point::value_type Dumux::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, b and c. More...
 
template<class Geometry >
static Geometry::ctype Dumux::distancePointTriangle (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the shortest distance from a point to a given triangle geometry. More...
 
template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointPolygon (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the shortest squared distance from a point to a given polygon geometry. More...
 
template<class Geometry >
static Geometry::ctype Dumux::distancePointPolygon (const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
 Compute the shortest distance from a point to a given polygon geometry. More...
 
template<class Geometry >
static Geometry::ctype Dumux::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. More...
 
template<class ctype , int dimWorld>
static ctype Dumux::distance (const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
 Compute the shortest distance between two points. More...
 
template<class ctype , int dimWorld>
static ctype Dumux::squaredDistance (const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
 Compute the shortest squared distance between two points. More...
 
template<class Geo1 , class Geo2 >
static auto Dumux::squaredDistance (const Geo1 &geo1, const Geo2 &geo2)
 Compute the shortest distance between two geometries. More...
 
template<class Geo1 , class Geo2 >
static auto Dumux::distance (const Geo1 &geo1, const Geo2 &geo2)
 Compute the shortest distance between two geometries. More...
 
template<class EntitySet , class ctype , int dimworld, typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0), int > = 0>
void Dumux::Detail::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. More...
 
template<class EntitySet , class ctype , int dimworld>
std::pair< ctype, std::size_t > Dumux::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) More...
 
template<class EntitySet , class ctype , int dimworld>
ctype Dumux::squaredDistance (const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, ctype minSquaredDistance=std::numeric_limits< ctype >::max())
 Compute the shortest squared distance to entities in an AABB tree. More...
 
template<class EntitySet , class ctype , int dimworld>
ctype Dumux::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. More...
 
template<class Geo1 , class Geo2 , class ctype , class GetFacetCornerIndices , class ComputeNormalFunction >
bool Dumux::Detail::computeSegmentIntersection (const Geo1 &geo1, const Geo2 &geo2, ctype baseEps, ctype &tfirst, ctype &tlast, const GetFacetCornerIndices &getFacetCornerIndices, const ComputeNormalFunction &computeNormal)
 Algorithm to find segment-like intersections of a polygon/polyhedron with a segment. The result is stored in the form of the local coordinates tfirst and tlast on the segment geo1. More...
 
template<class ctype >
int Dumux::getOrientation (const Dune::FieldVector< ctype, 3 > &a, const Dune::FieldVector< ctype, 3 > &b, const Dune::FieldVector< ctype, 3 > &c, const Dune::FieldVector< ctype, 3 > &normal)
 Returns the orientation of a sequence a-->b-->c in one plane (defined by normal vector) More...
 
template<int dim, class ctype , std::enable_if_t<(dim==2), int > = 0>
std::vector< Dune::FieldVector< ctype, 3 > > Dumux::grahamConvexHullImpl (std::vector< Dune::FieldVector< ctype, 3 > > &points)
 Compute the points making up the convex hull around the given set of unordered points. More...
 
template<int dim, class ctype , std::enable_if_t<(dim==2), int > = 0>
std::vector< Dune::FieldVector< ctype, 2 > > Dumux::grahamConvexHullImpl (const std::vector< Dune::FieldVector< ctype, 2 > > &points)
 Compute the points making up the convex hull around the given set of unordered points. More...
 
template<int dim, class ctype , int dimWorld>
std::vector< Dune::FieldVector< ctype, dimWorld > > Dumux::grahamConvexHull (std::vector< Dune::FieldVector< ctype, dimWorld > > &points)
 Compute the points making up the convex hull around the given set of unordered points. More...
 
template<int dim, class ctype , int dimWorld>
std::vector< Dune::FieldVector< ctype, dimWorld > > Dumux::grahamConvexHull (const std::vector< Dune::FieldVector< ctype, dimWorld > > &points)
 Compute the points making up the convex hull around the given set of unordered points. More...
 
template<class EntitySet , class ctype , int dimworld>
std::vector< std::size_t > Dumux::intersectingEntities (const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
 Compute all intersections between entities and a point. More...
 
template<class EntitySet , class ctype , int dimworld>
void Dumux::intersectingEntities (const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, std::size_t node, std::vector< std::size_t > &entities, bool isCartesianGrid=false)
 Compute intersections with point for all nodes of the bounding box tree recursively. More...
 
template<class Geometry , class EntitySet >
std::vector< IntersectionInfo< Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype > > Dumux::intersectingEntities (const Geometry &geometry, const BoundingBoxTree< EntitySet > &tree)
 Compute all intersections between a geometry and a bounding box tree. More...
 
template<class Geometry , class EntitySet >
void Dumux::intersectingEntities (const Geometry &geometry, const BoundingBoxTree< EntitySet > &tree, const std::array< typename Geometry::ctype, 2 *Geometry::coorddimension > &bBox, std::size_t nodeIdx, std::vector< IntersectionInfo< Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype > > &intersections)
 Compute intersections with point for all nodes of the bounding box tree recursively. More...
 
template<class EntitySet0 , class EntitySet1 >
std::vector< IntersectionInfo< EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype > > Dumux::intersectingEntities (const BoundingBoxTree< EntitySet0 > &treeA, const BoundingBoxTree< EntitySet1 > &treeB)
 Compute all intersections between two bounding box trees. More...
 
template<class EntitySet0 , class EntitySet1 >
void Dumux::intersectingEntities (const BoundingBoxTree< EntitySet0 > &treeA, const BoundingBoxTree< EntitySet1 > &treeB, std::size_t nodeA, std::size_t nodeB, std::vector< IntersectionInfo< EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype > > &intersections)
 Compute all intersections between two all bounding box tree nodes recursively. More...
 
template<class ctype , int dimworld>
std::size_t Dumux::intersectingEntityCartesianGrid (const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &min, const Dune::FieldVector< ctype, dimworld > &max, const std::array< int, std::size_t(dimworld)> &cells)
 Compute the index of the intersecting element of a Cartesian grid with a point The grid is given by the lower left corner (min), the upper right corner (max) and the number of cells in each direction (cells). More...
 
template<class ctype , int dimworld, class Geometry , typename std::enable_if_t<(Geometry::mydimension==3), int > = 0>
bool Dumux::intersectsPointGeometry (const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
 Find out whether a point is inside a three-dimensional geometry. More...
 
template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3), int > = 0>
bool Dumux::intersectsPointSimplex (const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
 Find out whether a point is inside the tetrahedron (p0, p1, p2, p3) (dimworld is 3) More...
 
template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3), int > = 0>
bool Dumux::intersectsPointSimplex (const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2)
 Find out whether a point is inside the triangle (p0, p1, p2) (dimworld is 3) More...
 
template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3||dimworld==2), int > = 0>
bool Dumux::intersectsPointSimplex (const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1)
 Find out whether a point is inside the interval (p0, p1) (dimworld is 2 or 3) More...
 
template<class CoordScalar >
bool Dumux::pointsAreCoplanar (const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points, const CoordScalar scale)
 Checks if four points lie within the same plane. More...
 
template<class CoordScalar >
bool Dumux::pointsAreCoplanar (const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
 Checks if four points lie within the same plane. More...
 
template<class CoordScalar >
std::vector< Dune::FieldVector< CoordScalar, 3 > > Dumux::getReorderedPoints (const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
 Returns a vector of points following the dune ordering. Convenience method that creates a temporary object in case no array of orientations is desired. More...
 
template<class CoordScalar >
std::vector< Dune::FieldVector< CoordScalar, 3 > > Dumux::getReorderedPoints (const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points, std::array< int, 4 > &orientations)
 Returns a vector of points following the dune ordering. More...
 
template<class CoordScalar , bool enableSanityCheck = true>
auto Dumux::makeDuneQuadrilaterial (const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
 Creates a dune quadrilateral geometry given 4 corner points. More...
 
template<class Vector >
Vector Dumux::normal (const Vector &v)
 Create a vector normal to the given one (v is expected to be non-zero) More...
 
template<class Vector >
Vector Dumux::unitNormal (const Vector &v)
 Create a vector normal to the given one (v is expected to be non-zero) More...
 
template<int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>, class RandomAccessContainer , std::enable_if_t< std::is_same_v< Policy, TriangulationPolicy::DelaunayPolicy > &&dim==1, int > = 0>
Triangulation< dim, dimWorld, typename RandomAccessContainer::value_type::value_type > Dumux::triangulate (const RandomAccessContainer &points)
 Triangulate area given points of a convex hull (1d) More...
 
template<int dim, class CornerF >
auto Dumux::convexPolytopeVolume (Dune::GeometryType type, const CornerF &c)
 Compute the volume of several common geometry types. More...
 
template<class Geometry >
auto Dumux::convexPolytopeVolume (const Geometry &geo)
 The volume of a given geometry. More...
 
template<class Geometry >
auto Dumux::volume (const Geometry &geo, unsigned int integrationOrder=4)
 The volume of a given geometry. More...
 
template<class Geometry , class Transformation >
auto Dumux::volume (const Geometry &geo, Transformation transformation, unsigned int integrationOrder=4)
 The volume of a given geometry with an extrusion/transformation policy. More...
 

Typedef Documentation

◆ DistanceField

template<class Geometry >
using Dumux::DistanceField = typedef AABBDistanceField<Geometry>

Class to calculate the closest distance from a point to a given set of geometries describing the domain's boundaries.

Template Parameters
GeometryThe (dune) geometry type.
Note
Defaults to Dumux::AABBDistanceField

◆ Triangulation

template<int dim, int dimWorld, class ctype >
using Dumux::Triangulation = typedef std::vector< std::array<Dune::FieldVector<ctype, dimWorld>, dim+1> >

The default data type to store triangulations.

Note
Stores each simplex separate and without connectivity information
We usually use this for sub-triangulation to allow for quadrature rules and interpolation on intersection geometries This is neither meant to be used to store large amounts of data nor as a mesh-like object

Function Documentation

◆ averageDistancePointGeometry()

template<class Geometry >
static Geometry::ctype Dumux::averageDistancePointGeometry ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry,
std::size_t  integrationOrder = 2 
)
inlinestatic

Compute the average distance from a point to a geometry by integration.

◆ averageDistanceSegmentGeometry()

template<class Geometry >
static Geometry::ctype Dumux::averageDistanceSegmentGeometry ( const typename Geometry::GlobalCoordinate &  a,
const typename Geometry::GlobalCoordinate &  b,
const Geometry &  geometry,
std::size_t  integrationOrder = 2 
)
inlinestatic

Compute the average distance from a segment to a geometry by integration.

◆ boundingSphere() [1/2]

template<class ConvexGeometry >
static Sphere< typename ConvexGeometry::ctype, ConvexGeometry::coorddimension > Dumux::boundingSphere ( const ConvexGeometry &  geometry)
inlinestatic

Computes a bounding sphere of a convex polytope geometry (Dune::Geometry interface)

Note
The bounding sphere is not necessarily the minimum enclosing sphere but computation is fast (AABB method)

◆ boundingSphere() [2/2]

template<class Scalar , int dim>
static Sphere< Scalar, dim > Dumux::boundingSphere ( const std::vector< Dune::FieldVector< Scalar, dim > > &  points)
inlinestatic

Computes a bounding sphere of points.

Note
The bounding sphere is not necessarily the minimum enclosing sphere but computation is fast (AABB method)

◆ center()

template<class Corners >
Corners::value_type Dumux::center ( const Corners &  corners)

The center of a given list of corners.

◆ circumSphereOfTriangle() [1/2]

template<class Geometry , typename std::enable_if_t< Geometry::mydimension==2, int > = 0>
static Sphere< typename Geometry::ctype, Geometry::coorddimension > Dumux::circumSphereOfTriangle ( const Geometry &  triangle)
inlinestatic

Computes the circumsphere of a triangle.

◆ circumSphereOfTriangle() [2/2]

template<class Point >
static Sphere< typename Point::value_type, Point::dimension > Dumux::circumSphereOfTriangle ( const Point &  a,
const Point &  b,
const Point &  c 
)
inlinestatic

Computes the circumsphere of a triangle.

Note
See https://www.ics.uci.edu/~eppstein/junkyard/circumcenter.html

◆ closestEntity() [1/2]

template<class EntitySet , class ctype , int dimworld>
std::pair< ctype, std::size_t > Dumux::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)

Parameters
pointthe point
treethe AABB tree
minSquaredDistanceconservative estimate of the minimum distance
Note
it is important that if an estimate is provided for minSquaredDistance to choose this estimate to be larger than the expected result. If the minSquaredDistance is smaller or equal to the result the returned entity index will be zero and the distance is equal to the estimate. However, this can also be the correct result. When in doubt, use the default parameter value.

◆ closestEntity() [2/2]

template<class EntitySet , class ctype , int dimworld, typename std::enable_if_t<(EntitySet::Entity::Geometry::mydimension > 0), int > = 0>
void Dumux::Detail::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.

Note
specialization for geometries with dimension larger than points
specialization for point geometries (point cloud AABB tree)

◆ computeSegmentIntersection()

template<class Geo1 , class Geo2 , class ctype , class GetFacetCornerIndices , class ComputeNormalFunction >
bool Dumux::Detail::computeSegmentIntersection ( const Geo1 &  geo1,
const Geo2 &  geo2,
ctype  baseEps,
ctype &  tfirst,
ctype &  tlast,
const GetFacetCornerIndices &  getFacetCornerIndices,
const ComputeNormalFunction &  computeNormal 
)

Algorithm to find segment-like intersections of a polygon/polyhedron with a segment. The result is stored in the form of the local coordinates tfirst and tlast on the segment geo1.

Parameters
geo1the first geometry
geo2the second geometry (dimGeo2 < dimGeo1)
baseEpsthe base epsilon used for floating point comparisons
tfirststores the local coordinate of beginning of intersection segment
tlaststores the local coordinate of the end of intersection segment
getFacetCornerIndicesFunction to obtain the facet corner indices on geo1
computeNormalFunction to obtain the normal vector on a facet on geo1

◆ convexPolytopeVolume() [1/2]

template<class Geometry >
auto Dumux::convexPolytopeVolume ( const Geometry &  geo)

The volume of a given geometry.

◆ convexPolytopeVolume() [2/2]

template<int dim, class CornerF >
auto Dumux::convexPolytopeVolume ( Dune::GeometryType  type,
const CornerF &  c 
)

Compute the volume of several common geometry types.

Parameters
typethe geometry type
ca function returning the ith corner (in Dune reference element order) e.g. [&](unsigned int i){ return corners[i]; }, where corners is a random access container storing the corners, and the returned corner is stored in a container (e.g. Dune::FieldVector) that has ::value_type and ::dimension.
Template Parameters
dimthe dimension of the geometry
CornerFthe function type (is deduced)
Returns
volume of the geometry or NaN signalling not implemented
Note
This is only correct for convex polytopes (flat sides)

◆ diameter()

template<class Geometry >
Geometry::ctype Dumux::diameter ( const Geometry &  geo)

Computes the longest distance between points of a geometry.

Note
Useful e.g. to compute the maximum cell diameter of a grid

◆ distance() [1/3]

template<class ctype , int dimWorld>
static ctype Dumux::distance ( const Dune::FieldVector< ctype, dimWorld > &  a,
const Dune::FieldVector< ctype, dimWorld > &  b 
)
inlinestatic

Compute the shortest distance between two points.

◆ distance() [2/3]

template<class EntitySet , class ctype , int dimworld>
ctype Dumux::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.

◆ distance() [3/3]

template<class Geo1 , class Geo2 >
static auto Dumux::distance ( const Geo1 &  geo1,
const Geo2 &  geo2 
)
inlinestatic

Compute the shortest distance between two geometries.

◆ distancePointLine() [1/2]

template<class Point >
static Point::value_type Dumux::distancePointLine ( const Point &  p,
const Point &  a,
const Point &  b 
)
inlinestatic

Compute the distance from a point to a line through the points a and b.

◆ distancePointLine() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::distancePointLine ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the distance from a point to a line given by a geometry with two corners.

Note
We currently lack the a representation of a line geometry. This convenience function assumes a segment geometry (with two corners) is passed which represents a line geometry.

◆ distancePointPolygon()

template<class Geometry >
static Geometry::ctype Dumux::distancePointPolygon ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the shortest distance from a point to a given polygon geometry.

Note
We only support triangles and quadrilaterals so far.

◆ distancePointSegment() [1/2]

template<class Point >
static Point::value_type Dumux::distancePointSegment ( const Point &  p,
const Point &  a,
const Point &  b 
)
inlinestatic

Compute the distance from a point to the segment connecting the points a and b.

◆ distancePointSegment() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::distancePointSegment ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the distance from a point to a given segment geometry.

◆ distancePointTriangle() [1/2]

template<class Point >
static Point::value_type Dumux::distancePointTriangle ( const Point &  p,
const Point &  a,
const Point &  b,
const Point &  c 
)
inlinestatic

Compute the shortest distance from a point to the triangle connecting the points a, b and c.

◆ distancePointTriangle() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::distancePointTriangle ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the shortest distance from a point to a given triangle geometry.

◆ getOrientation()

template<class ctype >
int Dumux::getOrientation ( const Dune::FieldVector< ctype, 3 > &  a,
const Dune::FieldVector< ctype, 3 > &  b,
const Dune::FieldVector< ctype, 3 > &  c,
const Dune::FieldVector< ctype, 3 > &  normal 
)

Returns the orientation of a sequence a-->b-->c in one plane (defined by normal vector)

Returns
-1 if a-->b-->c forms a counter-clockwise turn (given the normal vector) +1 for a clockwise turn, 0 if they are on one line (colinear)

◆ getReorderedPoints() [1/2]

template<class CoordScalar >
std::vector< Dune::FieldVector< CoordScalar, 3 > > Dumux::getReorderedPoints ( const std::vector< Dune::FieldVector< CoordScalar, 3 > > &  points)

Returns a vector of points following the dune ordering. Convenience method that creates a temporary object in case no array of orientations is desired.

Parameters
pointsThe user-specified vector of points (potentially in wrong order).

◆ getReorderedPoints() [2/2]

template<class CoordScalar >
std::vector< Dune::FieldVector< CoordScalar, 3 > > Dumux::getReorderedPoints ( const std::vector< Dune::FieldVector< CoordScalar, 3 > > &  points,
std::array< int, 4 > &  orientations 
)

Returns a vector of points following the dune ordering.

Parameters
pointsThe user-specified vector of points (potentially in wrong order).
orientationsAn array of orientations that can be useful for further processing.

◆ grahamConvexHull() [1/2]

template<int dim, class ctype , int dimWorld>
std::vector< Dune::FieldVector< ctype, dimWorld > > Dumux::grahamConvexHull ( const std::vector< Dune::FieldVector< ctype, dimWorld > > &  points)

Compute the points making up the convex hull around the given set of unordered points.

Note
We assume that all points are coplanar and there are no identical points in the list
This is the overload if we are not allowed to write into the given points vector

◆ grahamConvexHull() [2/2]

template<int dim, class ctype , int dimWorld>
std::vector< Dune::FieldVector< ctype, dimWorld > > Dumux::grahamConvexHull ( std::vector< Dune::FieldVector< ctype, dimWorld > > &  points)

Compute the points making up the convex hull around the given set of unordered points.

Note
We assume that all points are coplanar and there are no identical points in the list

◆ grahamConvexHullImpl() [1/2]

template<int dim, class ctype , std::enable_if_t<(dim==2), int > = 0>
std::vector< Dune::FieldVector< ctype, 2 > > Dumux::grahamConvexHullImpl ( const std::vector< Dune::FieldVector< ctype, 2 > > &  points)

Compute the points making up the convex hull around the given set of unordered points.

Note
This is the specialization for 2d space. Here, we make use of the generic implementation for the case of coplanar points in 3d space (a more efficient implementation could be provided).

◆ grahamConvexHullImpl() [2/2]

template<int dim, class ctype , std::enable_if_t<(dim==2), int > = 0>
std::vector< Dune::FieldVector< ctype, 3 > > Dumux::grahamConvexHullImpl ( std::vector< Dune::FieldVector< ctype, 3 > > &  points)

Compute the points making up the convex hull around the given set of unordered points.

Note
We assume that all points are coplanar and there are no identical points in the list
This algorithm changes the order of the given points a bit as they are unordered anyway this shouldn't matter too much

◆ intersectingEntities() [1/6]

template<class EntitySet0 , class EntitySet1 >
std::vector< IntersectionInfo< EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype > > Dumux::intersectingEntities ( const BoundingBoxTree< EntitySet0 > &  treeA,
const BoundingBoxTree< EntitySet1 > &  treeB 
)
inline

Compute all intersections between two bounding box trees.

◆ intersectingEntities() [2/6]

template<class EntitySet0 , class EntitySet1 >
void Dumux::intersectingEntities ( const BoundingBoxTree< EntitySet0 > &  treeA,
const BoundingBoxTree< EntitySet1 > &  treeB,
std::size_t  nodeA,
std::size_t  nodeB,
std::vector< IntersectionInfo< EntitySet0::dimensionworld, typename EntitySet0::ctype, typename EntitySet1::ctype > > &  intersections 
)

Compute all intersections between two all bounding box tree nodes recursively.

◆ intersectingEntities() [3/6]

template<class EntitySet , class ctype , int dimworld>
std::vector< std::size_t > Dumux::intersectingEntities ( const Dune::FieldVector< ctype, dimworld > &  point,
const BoundingBoxTree< EntitySet > &  tree,
bool  isCartesianGrid = false 
)
inline

Compute all intersections between entities and a point.

◆ intersectingEntities() [4/6]

template<class EntitySet , class ctype , int dimworld>
void Dumux::intersectingEntities ( const Dune::FieldVector< ctype, dimworld > &  point,
const BoundingBoxTree< EntitySet > &  tree,
std::size_t  node,
std::vector< std::size_t > &  entities,
bool  isCartesianGrid = false 
)

Compute intersections with point for all nodes of the bounding box tree recursively.

◆ intersectingEntities() [5/6]

template<class Geometry , class EntitySet >
std::vector< IntersectionInfo< Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype > > Dumux::intersectingEntities ( const Geometry &  geometry,
const BoundingBoxTree< EntitySet > &  tree 
)
inline

Compute all intersections between a geometry and a bounding box tree.

◆ intersectingEntities() [6/6]

template<class Geometry , class EntitySet >
void Dumux::intersectingEntities ( const Geometry &  geometry,
const BoundingBoxTree< EntitySet > &  tree,
const std::array< typename Geometry::ctype, 2 *Geometry::coorddimension > &  bBox,
std::size_t  nodeIdx,
std::vector< IntersectionInfo< Geometry::coorddimension, typename Geometry::ctype, typename EntitySet::ctype > > &  intersections 
)

Compute intersections with point for all nodes of the bounding box tree recursively.

◆ intersectingEntityCartesianGrid()

template<class ctype , int dimworld>
std::size_t Dumux::intersectingEntityCartesianGrid ( const Dune::FieldVector< ctype, dimworld > &  point,
const Dune::FieldVector< ctype, dimworld > &  min,
const Dune::FieldVector< ctype, dimworld > &  max,
const std::array< int, std::size_t(dimworld)> &  cells 
)
inline

Compute the index of the intersecting element of a Cartesian grid with a point The grid is given by the lower left corner (min), the upper right corner (max) and the number of cells in each direction (cells).

Note
If there are several options the lowest matching cell index will be returned
Returns the index i + I*j + I*J*k for the intersecting element i,j,k of a grid with I,J,K cells

◆ intersectsPointBoundingBox()

template<class ctype , int dimworld>
bool Dumux::intersectsPointBoundingBox ( const Dune::FieldVector< ctype, dimworld > &  point,
const Dune::FieldVector< ctype, dimworld > &  min,
const Dune::FieldVector< ctype, dimworld > &  max 
)
inline

Determine if a point intersects an axis-aligned bounding box The bounding box is given by the lower left corner (min) and the upper right corner (max)

◆ intersectsPointGeometry()

template<class ctype , int dimworld, class Geometry , typename std::enable_if_t<(Geometry::mydimension==3), int > = 0>
bool Dumux::intersectsPointGeometry ( const Dune::FieldVector< ctype, dimworld > &  point,
const Geometry &  g 
)

Find out whether a point is inside a three-dimensional geometry.

Find out whether a point is inside a one-dimensional geometry.

Find out whether a point is inside a two-dimensional geometry.

◆ intersectsPointSimplex() [1/3]

template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3||dimworld==2), int > = 0>
bool Dumux::intersectsPointSimplex ( const Dune::FieldVector< ctype, dimworld > &  point,
const Dune::FieldVector< ctype, dimworld > &  p0,
const Dune::FieldVector< ctype, dimworld > &  p1 
)

Find out whether a point is inside the interval (p0, p1) (dimworld is 2 or 3)

Find out whether a point is inside the interval (p0, p1) (dimworld is 1)

Note
We assume the given interval has non-zero length and use it to scale the epsilon

◆ intersectsPointSimplex() [2/3]

template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3), int > = 0>
bool Dumux::intersectsPointSimplex ( const Dune::FieldVector< ctype, dimworld > &  point,
const Dune::FieldVector< ctype, dimworld > &  p0,
const Dune::FieldVector< ctype, dimworld > &  p1,
const Dune::FieldVector< ctype, dimworld > &  p2 
)

Find out whether a point is inside the triangle (p0, p1, p2) (dimworld is 3)

Find out whether a point is inside the triangle (p0, p1, p2) (dimworld is 2)

◆ intersectsPointSimplex() [3/3]

template<class ctype , int dimworld, typename std::enable_if_t<(dimworld==3), int > = 0>
bool Dumux::intersectsPointSimplex ( const Dune::FieldVector< ctype, dimworld > &  point,
const Dune::FieldVector< ctype, dimworld > &  p0,
const Dune::FieldVector< ctype, dimworld > &  p1,
const Dune::FieldVector< ctype, dimworld > &  p2,
const Dune::FieldVector< ctype, dimworld > &  p3 
)

Find out whether a point is inside the tetrahedron (p0, p1, p2, p3) (dimworld is 3)

◆ makeDuneQuadrilaterial()

template<class CoordScalar , bool enableSanityCheck = true>
auto Dumux::makeDuneQuadrilaterial ( const std::vector< Dune::FieldVector< CoordScalar, 3 > > &  points)

Creates a dune quadrilateral geometry given 4 corner points.

Template Parameters
CoordScalarThe CoordScalar type.
enableSanityCheckTurn on/off sanity check and reordering of points
Parameters
pointsThe user-specified vector of points (potentially in wrong order).

◆ normal()

template<class Vector >
Vector Dumux::normal ( const Vector &  v)
inline

Create a vector normal to the given one (v is expected to be non-zero)

Note
This returns some orthogonal vector with arbitrary length

◆ pointsAreCoplanar() [1/2]

template<class CoordScalar >
bool Dumux::pointsAreCoplanar ( const std::vector< Dune::FieldVector< CoordScalar, 3 > > &  points)

Checks if four points lie within the same plane.

◆ pointsAreCoplanar() [2/2]

template<class CoordScalar >
bool Dumux::pointsAreCoplanar ( const std::vector< Dune::FieldVector< CoordScalar, 3 > > &  points,
const CoordScalar  scale 
)

Checks if four points lie within the same plane.

◆ squaredDistance() [1/3]

template<class ctype , int dimWorld>
static ctype Dumux::squaredDistance ( const Dune::FieldVector< ctype, dimWorld > &  a,
const Dune::FieldVector< ctype, dimWorld > &  b 
)
inlinestatic

Compute the shortest squared distance between two points.

◆ squaredDistance() [2/3]

template<class EntitySet , class ctype , int dimworld>
ctype Dumux::squaredDistance ( const Dune::FieldVector< ctype, dimworld > &  point,
const BoundingBoxTree< EntitySet > &  tree,
ctype  minSquaredDistance = std::numeric_limits<ctype>::max() 
)

Compute the shortest squared distance to entities in an AABB tree.

◆ squaredDistance() [3/3]

template<class Geo1 , class Geo2 >
static auto Dumux::squaredDistance ( const Geo1 &  geo1,
const Geo2 &  geo2 
)
inlinestatic

Compute the shortest distance between two geometries.

◆ squaredDistancePointLine() [1/2]

template<class Point >
static Point::value_type Dumux::squaredDistancePointLine ( const Point &  p,
const Point &  a,
const Point &  b 
)
inlinestatic

Compute the squared distance from a point to a line through the points a and b.

◆ squaredDistancePointLine() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointLine ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the squared distance from a point to a line given by a geometry with two corners.

Note
We currently lack the a representation of a line geometry. This convenience function assumes a segment geometry (with two corners) is passed which represents a line geometry.

◆ squaredDistancePointPolygon()

template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointPolygon ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the shortest squared distance from a point to a given polygon geometry.

Note
We only support triangles and quadrilaterals so far.

◆ squaredDistancePointSegment() [1/2]

template<class Point >
static Point::value_type Dumux::squaredDistancePointSegment ( const Point &  p,
const Point &  a,
const Point &  b 
)
inlinestatic

Compute the squared distance from a point to the segment connecting the points a and b.

◆ squaredDistancePointSegment() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointSegment ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the squared distance from a point to a given segment geometry.

◆ squaredDistancePointTriangle() [1/2]

template<class Point >
static Point::value_type Dumux::squaredDistancePointTriangle ( const Point &  p,
const Point &  a,
const Point &  b,
const Point &  c 
)
inlinestatic

Compute the shortest squared distance from a point to the triangle connecting the points a, b and c See https://www.iquilezles.org/www/articles/triangledistance/triangledistance.htm.

◆ squaredDistancePointTriangle() [2/2]

template<class Geometry >
static Geometry::ctype Dumux::squaredDistancePointTriangle ( const typename Geometry::GlobalCoordinate &  p,
const Geometry &  geometry 
)
inlinestatic

Compute the shortest squared distance from a point to a given triangle geometry.

◆ triangulate()

template<int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>, class RandomAccessContainer , std::enable_if_t< std::is_same_v< Policy, TriangulationPolicy::DelaunayPolicy > &&dim==1, int > = 0>
Triangulation< dim, dimWorld, typename RandomAccessContainer::value_type::value_type > Dumux::triangulate ( const RandomAccessContainer &  points)
inline

Triangulate area given points of a convex hull (1d)

Triangulate volume given a point cloud (3d)

Triangulate area given points (2d)

Triangulate area given points of a convex hull (2d)

Template Parameters
dimSpecifies the dimension of the resulting triangulation
dimWorldThe dimension of the coordinate space
PolicySpecifies the algorithm to be used for triangulation
Note
this specialization is for 1d discretization using segments
sorts the points and connects them via segments
Template Parameters
dimSpecifies the dimension of the resulting triangulation
dimWorldThe dimension of the coordinate space
PolicySpecifies the algorithm to be used for triangulation
Note
this specialization is for 2d triangulations using mid point policy
Assumes all points of the convex hull are coplanar
This inserts a mid point and connects all corners with that point to triangles
Assumes points are given as a ordered sequence representing the polyline forming the convex hull
Template Parameters
dimSpecifies the dimension of the resulting triangulation
dimWorldThe dimension of the coordinate space
PolicySpecifies the algorithm to be used for triangulation
Note
this specialization is for 2d triangulations using the convex hull policy That means we first construct the convex hull of the points and then triangulate the convex hull using the midpoint policy
Assumes points are unique and not all colinear (will throw a Dune::InvalidStateException)
Template Parameters
dimSpecifies the dimension of the resulting triangulation
dimWorldThe dimension of the coordinate space
PolicySpecifies the algorithm to be used for triangulation
Note
this specialization is for 3d triangulations using the convex hull policy
Assumes points are unique and not all coplanar
Todo:
sort points and create polyline

◆ unitNormal()

template<class Vector >
Vector Dumux::unitNormal ( const Vector &  v)
inline

Create a vector normal to the given one (v is expected to be non-zero)

Note
This returns some orthogonal vector with unit length

◆ volume() [1/2]

template<class Geometry , class Transformation >
auto Dumux::volume ( const Geometry &  geo,
Transformation  transformation,
unsigned int  integrationOrder = 4 
)

The volume of a given geometry with an extrusion/transformation policy.

Note
depending on the transformation this might not be an accurate quadrature rule anymore

◆ volume() [2/2]

template<class Geometry >
auto Dumux::volume ( const Geometry &  geo,
unsigned int  integrationOrder = 4 
)

The volume of a given geometry.