GEOS
Namespaces | Functions | Variables
ComputationalGeometry.hpp File Reference
#include "common/DataTypes.hpp"
#include "common/DataLayouts.hpp"
#include "finiteElement/elementFormulations/H1_Hexahedron_Lagrange1_GaussLegendre2.hpp"
#include "finiteElement/elementFormulations/H1_Pyramid_Lagrange1_Gauss5.hpp"
#include "finiteElement/elementFormulations/H1_Tetrahedron_Lagrange1_Gauss1.hpp"
#include "finiteElement/elementFormulations/H1_Wedge_Lagrange1_Gauss6.hpp"
#include "LvArray/src/output.hpp"
#include "LvArray/src/tensorOps.hpp"

Go to the source code of this file.

Namespaces

 geos
 

Functions

template<typename LINEDIR_TYPE , typename POINT_TYPE , typename NORMAL_TYPE , typename ORIGIN_TYPE , typename INTPOINT_TYPE >
void geos::computationalGeometry::LinePlaneIntersection (LINEDIR_TYPE const &lineDir, POINT_TYPE const &linePoint, NORMAL_TYPE const &planeNormal, ORIGIN_TYPE const &planeOrigin, INTPOINT_TYPE &intersectionPoint)
 Calculate the intersection between a line and a plane. More...
 
template<typename NORMAL_TYPE >
array1d< int > geos::computationalGeometry::orderPointsCCW (arrayView2d< real64 > const &points, NORMAL_TYPE const &normal)
 Reorder a set of points counter-clockwise. More...
 
template<typename NORMAL_TYPE >
real64 geos::computationalGeometry::ComputeSurfaceArea (arrayView2d< real64 const > const &points, NORMAL_TYPE const &&normal)
 Calculate the area of a polygon given the set of points in ccw order defining it. More...
 
template<localIndex DIMENSION, typename POINT_COORDS_TYPE >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 geos::computationalGeometry::computeDiameter (POINT_COORDS_TYPE points, localIndex const &numPoints)
 Calculate the diameter of a set of points in a given dimension. More...
 
template<typename CENTER_TYPE , typename NORMAL_TYPE >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 geos::computationalGeometry::centroid_3DPolygon (arraySlice1d< localIndex const > const pointsIndices, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &points, CENTER_TYPE &&center, NORMAL_TYPE &&normal, real64 const areaTolerance=0.0)
 Calculate the centroid of a convex 3D polygon as well as the normal and the rotation matrix. More...
 
template<typename NORMAL_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::FixNormalOrientation_3D (NORMAL_TYPE &&normal)
 Change the orientation of the input vector to be consistent in a global sense. More...
 
template<typename NORMAL_TYPE , typename MATRIX_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::RotationMatrix_3D (NORMAL_TYPE const &normal, MATRIX_TYPE &&rotationMatrix)
 Calculate the rotation matrix for a face in the 3D space. More...
 
template<typename T >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE int geos::computationalGeometry::sign (T const val)
 Return the sign of a given value as an integer. More...
 
template<typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::isPointInsidePolyhedron (arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &nodeCoordinates, arraySlice1d< localIndex const > const &faceIndices, ArrayOfArraysView< localIndex const > const &facesToNodes, POINT_TYPE const &elemCenter, POINT_TYPE const &point, real64 const areaTolerance=0.0)
 Check if a point is inside a convex polyhedron (3D polygon) More...
 
template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareVertex (POINT_TYPE const ax, POINT_TYPE const ay, POINT_TYPE const az, COORD_TYPE const bx, COORD_TYPE const by, COORD_TYPE const bz)
 Method to perform lexicographic comparison of two nodes based on coordinates. More...
 
template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareEdge (POINT_TYPE const ax, POINT_TYPE const ay, POINT_TYPE const az, COORD_TYPE const e1x, COORD_TYPE const e1y, COORD_TYPE const e1z, COORD_TYPE const e2x, COORD_TYPE const e2y, COORD_TYPE const e2z)
 Method to perform lexicographic comparison of a node and an edge based on coordinates. More...
 
template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareTriangle (POINT_TYPE const ax, POINT_TYPE const ay, POINT_TYPE const az, COORD_TYPE const t1x, COORD_TYPE const t1y, COORD_TYPE const t1z, COORD_TYPE const t2x, COORD_TYPE const t2y, COORD_TYPE const t2z, COORD_TYPE const t3x, COORD_TYPE const t3y, COORD_TYPE const t3z)
 Method to perform lexicographic comparison of a node and a triangle based on coordinates. More...
 
template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findVertexRefElement (arraySlice1d< localIndex const > const &nodeElements, arrayView1d< globalIndex const > const &elementGlobalIndex)
 Method to find the reference element touching a vertex. The element with the lowest global ID is chosen from the list. More...
 
template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findEdgeRefElement (arraySlice1d< localIndex const > const &nodeElements1, arraySlice1d< localIndex const > const &nodeElements2, arrayView1d< globalIndex const > const &elementGlobalIndex)
 Method to find the reference element for an edge. The element with the lowest global ID is chosen from the list. More...
 
template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findTriangleRefElement (arraySlice1d< localIndex const > const &nodeElements1, arraySlice1d< localIndex const > const &nodeElements2, arraySlice1d< localIndex const > const &nodeElements3, arrayView1d< globalIndex const > const &elementGlobalIndex)
 Method to find the reference element for a triangle. The element with the lowest global ID is chosen from the list. More...
 
template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::computeWindingNumber (localIndex element, arrayView2d< COORD_TYPE const, nodes::REFERENCE_POSITION_USD > const &nodeCoordinates, arrayView2d< localIndex const > const &elementsToFaces, ArrayOfArraysView< localIndex const > const &facesToNodes, ArrayOfArraysView< localIndex const > const &nodesToElements, arrayView1d< globalIndex const > const &nodeLocalToGlobal, arrayView1d< globalIndex const > const &elementLocalToGlobal, POINT_TYPE const &elemCenter, POINT_TYPE const &point)
 Computes the winding number of a point with respecto to a mesh element. More...
 
template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::isPointInsideConvexPolyhedronRobust (localIndex element, arrayView2d< COORD_TYPE const, nodes::REFERENCE_POSITION_USD > const &nodeCoordinates, arrayView2d< localIndex const > const &elementsToFaces, ArrayOfArraysView< localIndex const > const &facesToNodes, ArrayOfArraysView< localIndex const > const &nodesToElements, arrayView1d< globalIndex const > const &nodeLocalToGlobal, arrayView1d< globalIndex const > const &elementLocalToGlobal, POINT_TYPE const &elemCenter, POINT_TYPE const &point)
 Check if a point is inside a convex polyhedron (3D polygon), using a robust method to avoid ambiguity when the point lies on an interface. This method is based on the following method: More...
 
template<typename VEC_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::getBoundingBox (localIndex const elemIndex, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &pointIndices, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &pointCoordinates, VEC_TYPE &&boxDims)
 Compute the dimensions of the bounding box containing the element defined here by the coordinates of its vertices. More...
 
template<typename FE_TYPE >
GEOS_HOST_DEVICE real64 geos::computationalGeometry::elementVolume (real64 const (&X)[FE_TYPE::numNodes][3])
 Compute the volume of an element (tetrahedron, pyramid, wedge, hexahedron) More...
 
GEOS_HOST_DEVICE real64 geos::computationalGeometry::hexahedronVolume (real64 const (&X)[8][3])
 Compute the volume of an hexahedron. More...
 
GEOS_HOST_DEVICE real64 geos::computationalGeometry::tetrahedronVolume (real64 const (&X)[4][3])
 Compute the volume of an tetrahedron. More...
 
GEOS_HOST_DEVICE real64 geos::computationalGeometry::wedgeVolume (real64 const (&X)[6][3])
 Compute the volume of a wedge. More...
 
GEOS_HOST_DEVICE real64 geos::computationalGeometry::pyramidVolume (real64 const (&X)[5][3])
 Compute the volume of a pyramid. More...
 
template<integer N>
GEOS_HOST_DEVICE real64 geos::computationalGeometry::prismVolume (real64 const (&X)[2 *N][3])
 Compute the volume of a prism with N-sided polygon base. More...
 

Variables

constexpr real64 geos::computationalGeometry::machinePrecision = LvArray::NumericLimits< real64 >::epsilon
 Machine epsilon for double-precision calculations.
 

Function Documentation

◆ centroid_3DPolygon()

template<typename CENTER_TYPE , typename NORMAL_TYPE >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 geos::computationalGeometry::centroid_3DPolygon ( arraySlice1d< localIndex const > const  pointsIndices,
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &  points,
CENTER_TYPE &&  center,
NORMAL_TYPE &&  normal,
real64 const  areaTolerance = 0.0 
)

Calculate the centroid of a convex 3D polygon as well as the normal and the rotation matrix.

Template Parameters
CENTER_TYPEThe type of center.
NORMAL_TYPEThe type of normal.
Parameters
[in]pointsIndiceslist of index references for the points array in order (CW or CCW) about the polygon loop
[in]points3D point list
[out]center3D center of the given ordered polygon point list
[out]normalnormal to the face
[in]areaTolerancetolerance used in the geometric computations
Returns
area of the convex 3D polygon

if area < - areaTolerance, this function will throw an error, and if (- areaTolerance <= area <= areaTolerance), the area is set to zero

Definition at line 234 of file ComputationalGeometry.hpp.

◆ computeDiameter()

template<localIndex DIMENSION, typename POINT_COORDS_TYPE >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 geos::computationalGeometry::computeDiameter ( POINT_COORDS_TYPE  points,
localIndex const &  numPoints 
)

Calculate the diameter of a set of points in a given dimension.

Template Parameters
DIMENSIONThe dimensionality of the points.
POINT_COORDS_TYPEThe type of the container holding the point coordinates.
Parameters
[in]pointsThe container holding the coordinates of the points.
[in]numPointsThe number of points in the container.
Returns
The diameter of the set of points.

Definition at line 194 of file ComputationalGeometry.hpp.

◆ ComputeSurfaceArea()

template<typename NORMAL_TYPE >
real64 geos::computationalGeometry::ComputeSurfaceArea ( arrayView2d< real64 const > const &  points,
NORMAL_TYPE const &&  normal 
)

Calculate the area of a polygon given the set of points in ccw order defining it.

Template Parameters
NORMAL_TYPEthe type of normal
Parameters
[in]pointscoordinates of the points
[in]normalvector normal to the plane
Returns
the area of the polygon

Definition at line 153 of file ComputationalGeometry.hpp.

◆ computeWindingNumber()

template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::computeWindingNumber ( localIndex  element,
arrayView2d< COORD_TYPE const, nodes::REFERENCE_POSITION_USD > const &  nodeCoordinates,
arrayView2d< localIndex const > const &  elementsToFaces,
ArrayOfArraysView< localIndex const > const &  facesToNodes,
ArrayOfArraysView< localIndex const > const &  nodesToElements,
arrayView1d< globalIndex const > const &  nodeLocalToGlobal,
arrayView1d< globalIndex const > const &  elementLocalToGlobal,
POINT_TYPE const &  elemCenter,
POINT_TYPE const &  point 
)

Computes the winding number of a point with respecto to a mesh element.

Template Parameters
POINT_TYPEtype of point
Parameters
[in]elementthe element to be checked
[in]nodeCoordinatesa global array of nodal coordinates
[in]elementsToFacesmap from elements to faces
[in]facesToNodesmap from faces to nodes
[in]nodesToElementsmap from nodes to elements
[in]nodeLocalToGlobalglobal indices of nodes
[in]elementLocalToGlobalglobal indices of elements
[in]elemCentercoordinates of the element centroid
[in]pointcoordinates of the query point
Returns
the signed winding number, which is positive if and only if the point is inside the mesh element.

Definition at line 659 of file ComputationalGeometry.hpp.

◆ elementVolume()

template<typename FE_TYPE >
GEOS_HOST_DEVICE real64 geos::computationalGeometry::elementVolume ( real64 const (&)  X[FE_TYPE::numNodes][3])
inline

Compute the volume of an element (tetrahedron, pyramid, wedge, hexahedron)

Template Parameters
FE_TYPEthe type of finite element space
Parameters
[in]Xvertices of the element
Returns
the volume of the element

Definition at line 879 of file ComputationalGeometry.hpp.

◆ findEdgeRefElement()

template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findEdgeRefElement ( arraySlice1d< localIndex const > const &  nodeElements1,
arraySlice1d< localIndex const > const &  nodeElements2,
arrayView1d< globalIndex const > const &  elementGlobalIndex 
)

Method to find the reference element for an edge. The element with the lowest global ID is chosen from the list.

Parameters
[in]nodeElements1the list of elements adjacent to the first node
[in]nodeElements2the list of elements adjacent to the second node
[in]elementGlobalIndexthe global IDs for elements
Returns
the element shared by the two nodes, with the minimal global index

Definition at line 577 of file ComputationalGeometry.hpp.

◆ findTriangleRefElement()

template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findTriangleRefElement ( arraySlice1d< localIndex const > const &  nodeElements1,
arraySlice1d< localIndex const > const &  nodeElements2,
arraySlice1d< localIndex const > const &  nodeElements3,
arrayView1d< globalIndex const > const &  elementGlobalIndex 
)

Method to find the reference element for a triangle. The element with the lowest global ID is chosen from the list.

Parameters
[in]nodeElements1the list of elements adjacent to the first node
[in]nodeElements2the list of elements adjacent to the second node
[in]nodeElements3the list of elements adjacent to the third node
[in]elementGlobalIndexthe global IDs for elements
Returns
the element shared by the three nodes, with the minimal global index

Definition at line 613 of file ComputationalGeometry.hpp.

◆ findVertexRefElement()

template<typename ... LIST_TYPE>
GEOS_HOST_DEVICE int geos::computationalGeometry::findVertexRefElement ( arraySlice1d< localIndex const > const &  nodeElements,
arrayView1d< globalIndex const > const &  elementGlobalIndex 
)

Method to find the reference element touching a vertex. The element with the lowest global ID is chosen from the list.

Parameters
[in]nodeElementsthe list of elements adjacent to the vertex
[in]elementGlobalIndexthe global IDs for elements
Returns
element touching the vertex with the least global index

Definition at line 550 of file ComputationalGeometry.hpp.

◆ FixNormalOrientation_3D()

template<typename NORMAL_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::FixNormalOrientation_3D ( NORMAL_TYPE &&  normal)

Change the orientation of the input vector to be consistent in a global sense.

Template Parameters
NORMAL_TYPEtype of normal
Parameters
[in,out]normalnormal to the face

Definition at line 294 of file ComputationalGeometry.hpp.

◆ getBoundingBox()

template<typename VEC_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::getBoundingBox ( localIndex const  elemIndex,
arrayView2d< localIndex const, cells::NODE_MAP_USD > const &  pointIndices,
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &  pointCoordinates,
VEC_TYPE &&  boxDims 
)

Compute the dimensions of the bounding box containing the element defined here by the coordinates of its vertices.

Template Parameters
VEC_TYPEtype of boxDims
Parameters
[in]elemIndexindex of the element in pointIndices.
[in]pointIndicesthe indices of the vertices in pointCoordinates.
[in]pointCoordinatesthe vertices coordinates.
[out]boxDimsThe dimensions of the bounding box.

Definition at line 844 of file ComputationalGeometry.hpp.

◆ hexahedronVolume()

GEOS_HOST_DEVICE real64 geos::computationalGeometry::hexahedronVolume ( real64 const (&)  X[8][3])
inline

Compute the volume of an hexahedron.

Parameters
[in]Xvertices of the hexahedron
Returns
the volume of the hexahedron

Definition at line 896 of file ComputationalGeometry.hpp.

◆ isPointInsideConvexPolyhedronRobust()

template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::isPointInsideConvexPolyhedronRobust ( localIndex  element,
arrayView2d< COORD_TYPE const, nodes::REFERENCE_POSITION_USD > const &  nodeCoordinates,
arrayView2d< localIndex const > const &  elementsToFaces,
ArrayOfArraysView< localIndex const > const &  facesToNodes,
ArrayOfArraysView< localIndex const > const &  nodesToElements,
arrayView1d< globalIndex const > const &  nodeLocalToGlobal,
arrayView1d< globalIndex const > const &  elementLocalToGlobal,
POINT_TYPE const &  elemCenter,
POINT_TYPE const &  point 
)

Check if a point is inside a convex polyhedron (3D polygon), using a robust method to avoid ambiguity when the point lies on an interface. This method is based on the following method:

  • the winding number omega of the point with respect to the cell is used to determine if the point is inside (omega=1) or not (omega=0)
  • corner cases (point lying on a face, edge or vertex of the cell) are detected using a robust method based on lexicographical comparisons
  • these comparisons are made consistent across MPI ranks by consistently arranging items based on global indices (GIDs). In particular:
    • Faces are triangulated using the vertex with the smallest GID as root;
    • Edges and faces are described by vertices in increasing GID order
  • Finally, if the point lies on a (vertex, edge, face), it is assigned to the first shared element with the least global index
    Template Parameters
    POINT_TYPEtype of point
    Parameters
    [in]elementthe element to be checked
    [in]nodeCoordinatesa global array of nodal coordinates
    [in]elementsToFacesmap from elements to faces
    [in]facesToNodesmap from faces to nodes
    [in]nodesToElementsmap from nodes to elements
    [in]nodeLocalToGlobalglobal indices of nodes
    [in]elementLocalToGlobalglobal indices of elements
    [in]elemCentercoordinates of the element centroid
    [in]pointcoordinates of the query point
    Returns
    whether the point is inside

Definition at line 820 of file ComputationalGeometry.hpp.

◆ isPointInsidePolyhedron()

template<typename POINT_TYPE >
GEOS_HOST_DEVICE bool geos::computationalGeometry::isPointInsidePolyhedron ( arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &  nodeCoordinates,
arraySlice1d< localIndex const > const &  faceIndices,
ArrayOfArraysView< localIndex const > const &  facesToNodes,
POINT_TYPE const &  elemCenter,
POINT_TYPE const &  point,
real64 const  areaTolerance = 0.0 
)

Check if a point is inside a convex polyhedron (3D polygon)

Template Parameters
POINT_TYPEtype of point
Parameters
[in]nodeCoordinatesa global array of nodal coordinates
[in]faceIndicesglobal indices of the faces of the cell
[in]facesToNodesmap from face to nodes
[in]elemCentercoordinates of the element center
[in]pointcoordinates of the query point
[in]areaTolerancesame as in centroid_3DPolygon
Returns
whether the point is inside
Note
For faces with n>3 nodes that are non-planar, average normal is used

Definition at line 399 of file ComputationalGeometry.hpp.

◆ lexicographicalCompareEdge()

template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareEdge ( POINT_TYPE const  ax,
POINT_TYPE const  ay,
POINT_TYPE const  az,
COORD_TYPE const  e1x,
COORD_TYPE const  e1y,
COORD_TYPE const  e1z,
COORD_TYPE const  e2x,
COORD_TYPE const  e2y,
COORD_TYPE const  e2z 
)

Method to perform lexicographic comparison of a node and an edge based on coordinates.

Template Parameters
COORD_TYPEtype of coordinate
POINT_TYPEtype of point
Parameters
[in]axx-coordinate of the first vertex
[in]ayy-coordinate of the first vertex
[in]azz-coordinate of the first vertex
[in]e1xx-coordinate of the first edge vertex
[in]e1yy-coordinate of the first edge vertex
[in]e1zz-coordinate of the first edge vertex
[in]e2xx-coordinate of the second edge vertex
[in]e2yy-coordinate of the second edge vertex
[in]e2zz-coordinate of the second edge vertex
Returns
0 if the vertices lies on the edge, +1 if e > a and -1 if a > e

Definition at line 486 of file ComputationalGeometry.hpp.

◆ lexicographicalCompareTriangle()

template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareTriangle ( POINT_TYPE const  ax,
POINT_TYPE const  ay,
POINT_TYPE const  az,
COORD_TYPE const  t1x,
COORD_TYPE const  t1y,
COORD_TYPE const  t1z,
COORD_TYPE const  t2x,
COORD_TYPE const  t2y,
COORD_TYPE const  t2z,
COORD_TYPE const  t3x,
COORD_TYPE const  t3y,
COORD_TYPE const  t3z 
)

Method to perform lexicographic comparison of a node and a triangle based on coordinates.

Template Parameters
COORD_TYPEtype of coordinate
POINT_TYPEtype of point
Parameters
[in]axx-coordinate of the first vertex
[in]ayy-coordinate of the first vertex
[in]azz-coordinate of the first vertex
[in]t1xx-coordinate of the first triangle vertex
[in]t1yy-coordinate of the first triangle vertex
[in]t1zz-coordinate of the first triangle vertex
[in]t2xx-coordinate of the second triangle vertex
[in]t2yy-coordinate of the second triangle vertex
[in]t2zz-coordinate of the second triangle vertex
[in]t3xx-coordinate of the third triangle vertex
[in]t3yy-coordinate of the third triangle vertex
[in]t3zz-coordinate of the third triangle vertex
Returns
0 if the vertices lies on the triangle, +1 if t > a and -1 if t > e

Definition at line 518 of file ComputationalGeometry.hpp.

◆ lexicographicalCompareVertex()

template<typename COORD_TYPE , typename POINT_TYPE >
GEOS_HOST_DEVICE int geos::computationalGeometry::lexicographicalCompareVertex ( POINT_TYPE const  ax,
POINT_TYPE const  ay,
POINT_TYPE const  az,
COORD_TYPE const  bx,
COORD_TYPE const  by,
COORD_TYPE const  bz 
)

Method to perform lexicographic comparison of two nodes based on coordinates.

Template Parameters
COORD_TYPEtype of coordinate
POINT_TYPEtype of point
Parameters
[in]axx-coordinate of the first vertex
[in]ayy-coordinate of the first vertex
[in]azz-coordinate of the first vertex
[in]bxx-coordinate of the second vertex
[in]byy-coordinate of the second vertex
[in]bzz-coordinate of the second vertex
Returns
0 if the vertices coincide, +1 if a > b and -1 if b > a

Definition at line 451 of file ComputationalGeometry.hpp.

◆ LinePlaneIntersection()

template<typename LINEDIR_TYPE , typename POINT_TYPE , typename NORMAL_TYPE , typename ORIGIN_TYPE , typename INTPOINT_TYPE >
void geos::computationalGeometry::LinePlaneIntersection ( LINEDIR_TYPE const &  lineDir,
POINT_TYPE const &  linePoint,
NORMAL_TYPE const &  planeNormal,
ORIGIN_TYPE const &  planeOrigin,
INTPOINT_TYPE &  intersectionPoint 
)

Calculate the intersection between a line and a plane.

Template Parameters
LINEDIR_TYPEthe type of lineDir
POINT_TYPEthe type of linePoint
NORMAL_TYPEthe type of planeNormal
ORIGIN_TYPEthe type of planeOrigin
INTPOINT_TYPEthe type of instersectionPoint
Parameters
[in]lineDirvector defining direction of the line
[in]linePointone point of the line
[in]planeNormalnormal to plane
[in]planeOriginplane origin
[out]intersectionPointthe intersection point

Definition at line 58 of file ComputationalGeometry.hpp.

◆ orderPointsCCW()

template<typename NORMAL_TYPE >
array1d< int > geos::computationalGeometry::orderPointsCCW ( arrayView2d< real64 > const &  points,
NORMAL_TYPE const &  normal 
)

Reorder a set of points counter-clockwise.

Template Parameters
NORMAL_TYPEthe type of normal
Parameters
[in]pointscoordinates of the points
[in]normalvector normal to the plane
Returns
an std::vector containing the original indices of the reordered points.

Definition at line 87 of file ComputationalGeometry.hpp.

◆ prismVolume()

template<integer N>
GEOS_HOST_DEVICE real64 geos::computationalGeometry::prismVolume ( real64 const (&)  X[2 *N][3])
inline

Compute the volume of a prism with N-sided polygon base.

Template Parameters
Nthe number of sides in the polygon base
Parameters
[in]Xvertices of the prism
Returns
the volume of the prism
Note
The volume is computed splitting the prism into wedges. The function can be called only for N > 5. For N = 3 and N = 4 function wedgeVolume and hexahedronVolume, respectively, should be used.

Definition at line 950 of file ComputationalGeometry.hpp.

◆ pyramidVolume()

GEOS_HOST_DEVICE real64 geos::computationalGeometry::pyramidVolume ( real64 const (&)  X[5][3])
inline

Compute the volume of a pyramid.

Parameters
[in]Xvertices of the pyramid
Returns
the volume of the pyramid

Definition at line 932 of file ComputationalGeometry.hpp.

◆ RotationMatrix_3D()

template<typename NORMAL_TYPE , typename MATRIX_TYPE >
GEOS_HOST_DEVICE void geos::computationalGeometry::RotationMatrix_3D ( NORMAL_TYPE const &  normal,
MATRIX_TYPE &&  rotationMatrix 
)

Calculate the rotation matrix for a face in the 3D space.

Template Parameters
NORMAL_TYPEtype of normal
MATRIX_TYPEtype of rotationMatrix
Parameters
[in]normalnormal to the face
[out]rotationMatrixrotation matrix for the face

Definition at line 331 of file ComputationalGeometry.hpp.

◆ sign()

template<typename T >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE int geos::computationalGeometry::sign ( T const  val)

Return the sign of a given value as an integer.

Template Parameters
Ttype of value
Parameters
valthe value in question
Returns
-1, 0 or 1 depending on whether the value is negative, zero or positive

Definition at line 379 of file ComputationalGeometry.hpp.

◆ tetrahedronVolume()

GEOS_HOST_DEVICE real64 geos::computationalGeometry::tetrahedronVolume ( real64 const (&)  X[4][3])
inline

Compute the volume of an tetrahedron.

Parameters
[in]Xvertices of the tetrahedron
Returns
the volume of the tetrahedron

Definition at line 908 of file ComputationalGeometry.hpp.

◆ wedgeVolume()

GEOS_HOST_DEVICE real64 geos::computationalGeometry::wedgeVolume ( real64 const (&)  X[6][3])
inline

Compute the volume of a wedge.

Parameters
[in]Xvertices of the wedge
Returns
the volume of the wedge

Definition at line 920 of file ComputationalGeometry.hpp.