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

Basic geometries in DuMux More...

Description

Basic geometries in DuMux

Files

file  boundingboxtree.hh
 An axis-aligned bounding box volume hierarchy for dune grids.
 
file  boundingboxtreeintersection.hh
 A class storing intersections from intersecting two bounding box trees.
 
file  diameter.hh
 A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
 
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 and a function to triangulate the polygon area spanned by the convex hull.
 
file  intersectingentities.hh
 Algorithms that finds which geometric entites intersect.
 
file  intersectionentityset.hh
 A class representing the intersection entites 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  refinementquadraturerule.hh
 A quadrature based on refinement.
 
file  triangulation.hh
 Functionality to triangulate point clouds.
 

Classes

class  Dumux::BoundingBoxTree< GeometricEntitySet >
 An axis-aligned bounding box volume tree implementation. More...
 
class  Dumux::BoundingBoxTreeIntersection< EntitySet0, EntitySet1 >
 An intersection object resulting from the intersection of two bounding box tree primitives. 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, 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, 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 entites two geometric entity sets. More...
 
class  Dumux::RefinementQuadratureRule< ct, mydim >
 A "quadrature" based on virtual refinement. More...
 

Functions

template<class Geometry >
Geometry::ctype Dumux::diameter (const Geometry &geo)
 Computes the longest distance between points of a geometry. 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, 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 a 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 a triangle (p0, p1, p2, p3) (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 a interval (p0, p1, p2, p3) (dimworld is 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<int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>, class RandomAccessContainer , std::enable_if_t< std::is_same< Policy, TriangulationPolicy::MidPointPolicy >::value &&dim==2, int > = 0>
Triangulation< dim, dimWorld, typename RandomAccessContainer::value_type::value_type > Dumux::triangulate (const RandomAccessContainer &convexHull)
 Triangulate area given points of a convex hull. More...
 

Function Documentation

◆ 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

◆ 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 indentical 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 indentical 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 indentical 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.

◆ 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 a interval (p0, p1, p2, p3) (dimworld is 3)

Find out whether a point is inside a interval (p0, p1, p2, p3) (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 a triangle (p0, p1, p2, p3) (dimworld is 3)

Find out whether a point is inside a triangle (p0, p1, p2, p3) (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 a 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).

◆ 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.

◆ triangulate()

template<int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>, class RandomAccessContainer , std::enable_if_t< std::is_same< Policy, TriangulationPolicy::MidPointPolicy >::value &&dim==2, 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.

Template Parameters
dimSpecifies the dimension of the resulting triangulation (2 or 3)
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
Template Parameters
dimSpecifies the dimension of the resulting triangulation (2 or 3)
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
Todo:
sort points and create polyline
Todo:
sort points and create polyline