20 #ifndef GEOS_MESH_UTILITIES_COMPUTATIONALGEOMETRY_HPP_
21 #define GEOS_MESH_UTILITIES_COMPUTATIONALGEOMETRY_HPP_
29 #include "LvArray/src/output.hpp"
30 #include "LvArray/src/tensorOps.hpp"
34 namespace computationalGeometry
53 template<
typename LINEDIR_TYPE,
57 typename INTPOINT_TYPE >
59 POINT_TYPE
const & linePoint,
60 NORMAL_TYPE
const & planeNormal,
61 ORIGIN_TYPE
const & planeOrigin,
62 INTPOINT_TYPE & intersectionPoint )
70 real64 dummy[ 3 ] = LVARRAY_TENSOROPS_INIT_LOCAL_3( planeOrigin );
71 LvArray::tensorOps::subtract< 3 >( dummy, linePoint );
72 real64 const d = LvArray::tensorOps::AiBi< 3 >( dummy, planeNormal ) /
73 LvArray::tensorOps::AiBi< 3 >( lineDir, planeNormal );
75 LvArray::tensorOps::copy< 3 >( intersectionPoint, linePoint );
76 LvArray::tensorOps::scaledAdd< 3 >( intersectionPoint, lineDir, d );
86 template<
typename NORMAL_TYPE >
88 NORMAL_TYPE
const & normal )
99 LvArray::tensorOps::fill< 3 >( centroid, 0 );
102 LvArray::tensorOps::add< 3 >( centroid, points[ a ] );
106 LvArray::tensorOps::scale< 3 >( centroid, 1.0 / numPoints );
108 real64 v0[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( centroid );
109 LvArray::tensorOps::subtract< 3 >( v0, points[ 0 ] );
110 LvArray::tensorOps::normalize< 3 >( v0 );
116 real64 v[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( centroid );
117 LvArray::tensorOps::subtract< 3 >( v, points[ a ] );
118 real64 const dot = LvArray::tensorOps::AiBi< 3 >( v, v0 );
121 LvArray::tensorOps::crossProduct( crossProduct, v, v0 );
122 real64 const det = LvArray::tensorOps::AiBi< 3 >( normal, crossProduct );
124 angle[ a ] = std::atan2( det, dot );
128 std::sort( indices.begin(), indices.end(), [&](
int i,
int j ) { return angle[ i ] < angle[ j ]; } );
134 LvArray::tensorOps::copy< 3 >( orderedPoints[ a ], points[ indices[ a ] ] );
139 LvArray::tensorOps::copy< 3 >( points[a], orderedPoints[a] );
152 template<
typename NORMAL_TYPE >
154 NORMAL_TYPE
const && normal )
160 for(
localIndex a = 0; a < points.size( 0 ); a++ )
162 LvArray::tensorOps::copy< 3 >( orderedPoints[a], points[a] );
167 for(
localIndex a = 0; a < points.size( 0 ) - 2; ++a )
169 real64 v1[ 3 ] = LVARRAY_TENSOROPS_INIT_LOCAL_3( orderedPoints[ a + 1 ] );
170 real64 v2[ 3 ] = LVARRAY_TENSOROPS_INIT_LOCAL_3( orderedPoints[ a + 2 ] );
172 LvArray::tensorOps::subtract< 3 >( v1, orderedPoints[ 0 ] );
173 LvArray::tensorOps::subtract< 3 >( v2, orderedPoints[ 0 ] );
175 real64 triangleNormal[ 3 ];
176 LvArray::tensorOps::crossProduct( triangleNormal, v1, v2 );
177 surfaceArea += LvArray::tensorOps::l2Norm< 3 >( triangleNormal );
180 return surfaceArea * 0.5;
191 template< localIndex DIMENSION,
typename POINT_COORDS_TYPE >
198 for(
localIndex numPoint = 0; numPoint < numPoints; ++numPoint )
200 for(
localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
202 real64 candidateDiameter = 0.0;
205 real64 coordDiff = points[numPoint][i] - points[numOthPoint][i];
206 candidateDiameter += coordDiff * coordDiff;
208 if( diameter < candidateDiameter )
210 diameter = candidateDiameter;
214 return LvArray::math::sqrt< real64 >( diameter );
231 template<
typename CENTER_TYPE,
typename NORMAL_TYPE >
236 CENTER_TYPE && center,
237 NORMAL_TYPE && normal,
238 real64 const areaTolerance = 0.0 )
241 LvArray::tensorOps::fill< 3 >( center, 0 );
242 LvArray::tensorOps::fill< 3 >( normal, 0 );
244 localIndex const numberOfPoints = pointsIndices.size();
248 real64 current[ 3 ], next[ 3 ], crossProduct[ 3 ];
250 LvArray::tensorOps::copy< 3 >( next, points[ pointsIndices[ numberOfPoints - 1 ] ] );
254 LvArray::tensorOps::copy< 3 >( current, next );
255 LvArray::tensorOps::copy< 3 >( next, points[ pointsIndices[ a ] ] );
257 LvArray::tensorOps::crossProduct( crossProduct, current, next );
259 LvArray::tensorOps::add< 3 >( normal, crossProduct );
260 LvArray::tensorOps::add< 3 >( center, next );
263 area = LvArray::tensorOps::l2Norm< 3 >( normal );
264 LvArray::tensorOps::scale< 3 >( center, 1.0 / numberOfPoints );
266 if( area > areaTolerance )
268 LvArray::tensorOps::normalize< 3 >( normal );
271 else if( area < -areaTolerance )
275 GEOS_LOG_RANK(
"Points: " << points[ pointsIndices[ a ] ] <<
" " << pointsIndices[ a ] );
277 GEOS_ERROR(
"Negative area found : " << area );
292 template<
typename NORMAL_TYPE >
300 if( normal[ 2 ] <= -orientationTolerance )
302 LvArray::tensorOps::scale< 3 >( normal, -1.0 );
304 else if( std::fabs( normal[ 2 ] ) < orientationTolerance )
307 if( normal[ 1 ] <= -orientationTolerance )
309 LvArray::tensorOps::scale< 3 >( normal, -1.0 );
311 else if( fabs( normal[ 1 ] ) < orientationTolerance )
314 if( normal[ 0 ] <= -orientationTolerance )
316 LvArray::tensorOps::scale< 3 >( normal, -1.0 );
329 template<
typename NORMAL_TYPE,
typename MATRIX_TYPE >
332 MATRIX_TYPE && rotationMatrix )
334 real64 m1[ 3 ] = { normal[ 2 ], 0.0, -normal[ 0 ] };
335 real64 m2[ 3 ] = { 0.0, normal[ 2 ], -normal[ 1 ] };
336 real64 const norm_m1 = LvArray::tensorOps::l2Norm< 3 >( m1 );
337 real64 const norm_m2 = LvArray::tensorOps::l2Norm< 3 >( m2 );
343 LvArray::tensorOps::crossProduct( m2, normal, m1 );
344 LvArray::tensorOps::normalize< 3 >( m2 );
345 LvArray::tensorOps::normalize< 3 >( m1 );
349 LvArray::tensorOps::crossProduct( m1, normal, m2 );
350 LvArray::tensorOps::scale< 3 >( m1, -1 );
351 LvArray::tensorOps::normalize< 3 >( m1 );
352 LvArray::tensorOps::normalize< 3 >( m2 );
356 rotationMatrix[ 0 ][ 0 ] = normal[ 0 ];
357 rotationMatrix[ 1 ][ 0 ] = normal[ 1 ];
358 rotationMatrix[ 2 ][ 0 ] = normal[ 2 ];
359 rotationMatrix[ 0 ][ 1 ] = m1[ 0 ];
360 rotationMatrix[ 1 ][ 1 ] = m1[ 1 ];
361 rotationMatrix[ 2 ][ 1 ] = m1[ 2 ];
362 rotationMatrix[ 0 ][ 2 ] = m2[ 0 ];
363 rotationMatrix[ 1 ][ 2 ] = m2[ 1 ];
364 rotationMatrix[ 2 ][ 2 ] = m2[ 2 ];
367 "Rotation matrix with determinant different from +1.0" );
376 template<
typename T >
381 return (T( 0 ) < val) - (val < T( 0 ));
397 template<
typename POINT_TYPE >
402 POINT_TYPE
const & elemCenter,
403 POINT_TYPE
const & point,
404 real64 const areaTolerance = 0.0 )
406 localIndex const numFaces = faceIndices.size();
407 R1Tensor faceCenter, faceNormal, cellToFaceVec;
413 centroid_3DPolygon( facesToNodes[faceIndex], nodeCoordinates, faceCenter, faceNormal, areaTolerance );
416 LvArray::tensorOps::copy< 3 >( cellToFaceVec, faceCenter );
417 LvArray::tensorOps::subtract< 3 >( cellToFaceVec, elemCenter );
418 if( LvArray::tensorOps::AiBi< 3 >( cellToFaceVec, faceNormal ) < 0.0 )
420 LvArray::tensorOps::scale< 3 >( faceNormal, -1 );
424 LvArray::tensorOps::subtract< 3 >( faceCenter, point );
425 int const s =
sign( LvArray::tensorOps::AiBi< 3 >( faceNormal, faceCenter ) );
449 template<
typename COORD_TYPE,
typename POINT_TYPE >
452 COORD_TYPE
const bx, COORD_TYPE
const by, COORD_TYPE
const bz )
484 template<
typename COORD_TYPE,
typename POINT_TYPE >
487 COORD_TYPE
const e1x, COORD_TYPE
const e1y, COORD_TYPE
const e1z,
488 COORD_TYPE
const e2x, COORD_TYPE
const e2y, COORD_TYPE
const e2z )
491 ( e1z - az ) * ( e2x - ax ),
492 ( e1z - az ) * ( e2y - ay ),
493 ( e1x - ax ) * ( e2y - ay ),
494 ( e1x - ax ) * ( e2z - az ),
495 ( e1y - ay ) * ( e2z - az ) );
516 template<
typename COORD_TYPE,
typename POINT_TYPE >
519 COORD_TYPE
const t1x, COORD_TYPE
const t1y, COORD_TYPE
const t1z,
520 COORD_TYPE
const t2x, COORD_TYPE
const t2y, COORD_TYPE
const t2z,
521 COORD_TYPE
const t3x, COORD_TYPE
const t3y, COORD_TYPE
const t3z )
523 COORD_TYPE v1x = t1x - ax;
524 COORD_TYPE v1y = t1y - ay;
525 COORD_TYPE v1z = t1z - az;
526 COORD_TYPE v2x = t2x - ax;
527 COORD_TYPE v2y = t2y - ay;
528 COORD_TYPE v2z = t2z - az;
529 COORD_TYPE v3x = t3x - ax;
530 COORD_TYPE v3y = t3y - ay;
531 COORD_TYPE v3z = t3z - az;
532 COORD_TYPE
sign = ( v1x * v2y - v1y * v2x ) * v3z +
533 ( v2x * v3y - v2y * v3x ) * v1z +
534 ( v3x * v1y - v3y * v1x ) * v2z;
548 template<
typename ... LIST_TYPE >
554 globalIndex minElementGID = LvArray::NumericLimits< globalIndex >::max;
555 for(
int i = 0; i < nodeElements.size(); i++ )
558 if( elementGlobalIndex[ e ] < minElementGID )
560 minElementGID = elementGlobalIndex[ e ];
575 template<
typename ... LIST_TYPE >
582 globalIndex minElementGID = LvArray::NumericLimits< globalIndex >::max;
583 for(
int i = 0; i < nodeElements1.size(); i++ )
586 for(
int j = 0; j < nodeElements2.size(); j++ )
591 if( elementGlobalIndex[ e1 ] < minElementGID )
593 minElementGID = elementGlobalIndex[ e1 ];
611 template<
typename ... LIST_TYPE >
619 globalIndex minElementGID = LvArray::NumericLimits< globalIndex >::max;
620 for(
int i = 0; i < nodeElements1.size(); i++ )
623 for(
int j = 0; j < nodeElements2.size(); j++ )
626 for(
int k = 0; k < nodeElements3.size(); k++ )
629 if( e1 == e2 && e2 == e3 )
631 if( elementGlobalIndex[ e1 ] < minElementGID )
633 minElementGID = elementGlobalIndex[ e1 ];
657 template<
typename COORD_TYPE,
typename POINT_TYPE >
666 POINT_TYPE
const & elemCenter,
667 POINT_TYPE
const & point )
670 localIndex const numFaces = faceIndices.size();
677 globalIndex minGlobalId = LvArray::NumericLimits< globalIndex >::max;
679 localIndex numFaceVertices = facesToNodes[faceIndex].size();
680 for(
localIndex v = 0; v < numFaceVertices; v++ )
682 localIndex vIndex = facesToNodes( faceIndex, v );
683 globalIndex globalId = nodeLocalToGlobal[ vIndex ];
684 if( globalId < minGlobalId )
686 minGlobalId = globalId;
692 for(
localIndex v = 0; v < numFaceVertices; v++ )
694 vi[ 1 ] = facesToNodes( faceIndex, v );
695 vi[ 2 ] = facesToNodes( faceIndex, (v + 1) % numFaceVertices );
696 if( vi[ 1 ] != minVertex && vi[ 2 ] != minVertex )
699 if( nodeLocalToGlobal[ vi[ 1 ] ] > nodeLocalToGlobal[ vi[ 2 ] ] )
705 COORD_TYPE v1x = nodeCoordinates( vi[ 0 ], 0 );
706 COORD_TYPE v1y = nodeCoordinates( vi[ 0 ], 1 );
707 COORD_TYPE v1z = nodeCoordinates( vi[ 0 ], 2 );
708 COORD_TYPE v2x = nodeCoordinates( vi[ 1 ], 0 );
709 COORD_TYPE v2y = nodeCoordinates( vi[ 1 ], 1 );
710 COORD_TYPE v2z = nodeCoordinates( vi[ 1 ], 2 );
711 COORD_TYPE v3x = nodeCoordinates( vi[ 2 ], 0 );
712 COORD_TYPE v3y = nodeCoordinates( vi[ 2 ], 1 );
713 COORD_TYPE v3z = nodeCoordinates( vi[ 2 ], 2 );
715 R1Tensor vv1 = { v2x - v1x, v2y - v1y, v2z - v1z };
716 R1Tensor vv2 = { v3x - v1x, v3y - v1y, v3z - v1z };
717 R1Tensor dist = { elemCenter[ 0 ] - ( v1x + v2x + v3x )/3.0,
718 elemCenter[ 1 ] - ( v1y + v2y + v3y )/3.0,
719 elemCenter[ 2 ] - ( v1z + v2z + v3z )/3.0 };
721 LvArray::tensorOps::crossProduct( norm, vv1, vv2 );
723 int sign = LvArray::tensorOps::AiBi< 3 >( norm, dist ) > 0 ? -1 : +1;
749 return findEdgeRefElement( nodesToElements[ vi[ 0 ] ], nodesToElements[ vi[ 1 ] ], elementLocalToGlobal ) == element;
751 facecmp +=
sign * edgecmp;
760 return findEdgeRefElement( nodesToElements[ vi[ 1 ] ], nodesToElements[ vi[ 2 ] ], elementLocalToGlobal ) == element;
762 facecmp +=
sign * edgecmp;
771 return findEdgeRefElement( nodesToElements[ vi[ 0 ] ], nodesToElements[ vi[ 2 ] ], elementLocalToGlobal ) == element;
773 facecmp +=
sign * edgecmp;
785 return findTriangleRefElement( nodesToElements[ vi[ 0 ] ], nodesToElements[ vi[ 1 ] ], nodesToElements[ vi[ 2 ] ], elementLocalToGlobal ) == element;
787 omega +=
sign * facecmp;
818 template<
typename COORD_TYPE,
typename POINT_TYPE >
827 POINT_TYPE
const & elemCenter,
828 POINT_TYPE
const & point )
830 return computeWindingNumber( element, nodeCoordinates, elementsToFaces, facesToNodes, nodesToElements, nodeLocalToGlobal, elementLocalToGlobal, elemCenter, point ) > 0;
842 template<
typename VEC_TYPE >
847 VEC_TYPE && boxDims )
850 R1Tensor minCoords = { LvArray::NumericLimits< real64 >::max,
851 LvArray::NumericLimits< real64 >::max,
852 LvArray::NumericLimits< real64 >::max };
855 LvArray::tensorOps::fill< 3 >( boxDims, LvArray::NumericLimits< real64 >::lowest );
858 for(
localIndex a = 0; a < pointIndices.size( 1 ); ++a )
860 localIndex const id = pointIndices( elemIndex, a );
863 minCoords[ d ] = fmin( minCoords[ d ], pointCoordinates(
id, d ) );
864 boxDims[ d ] = fmax( boxDims[ d ], pointCoordinates(
id, d ) );
868 LvArray::tensorOps::subtract< 3 >( boxDims, minCoords );
877 template<
typename FE_TYPE >
882 for(
localIndex q=0; q<FE_TYPE::numQuadraturePoints; ++q )
884 result = result + FE_TYPE::transformedQuadratureWeight( q, X );
898 return elementVolume< finiteElement::H1_Hexahedron_Lagrange1_GaussLegendre2 >( X );
910 return elementVolume< finiteElement::H1_Tetrahedron_Lagrange1_Gauss1 >( X );
922 return elementVolume< finiteElement::H1_Wedge_Lagrange1_Gauss6 >( X );
934 return elementVolume< finiteElement::H1_Pyramid_Lagrange1_Gauss5 >( X );
947 template<
integer N >
952 static_assert( N > 4,
953 "Function prismVolume can be called for a prism with N-sided polygon base where N > 5." );
960 for(
integer a = 0; a < N; ++a )
962 LvArray::tensorOps::add< 3 >( XGBot, X[a] );
964 for(
integer a = N; a < 2 * N; ++a )
966 LvArray::tensorOps::add< 3 >( XGTop, X[a] );
968 LvArray::tensorOps::scale< 3 >( XGBot, 1.0 / N );
969 LvArray::tensorOps::scale< 3 >( XGTop, 1.0 / N );
972 for(
int a = 0; a < N - 1; ++a )
975 LvArray::tensorOps::copy< 3 >( XWedge[0], X[a] );
976 LvArray::tensorOps::copy< 3 >( XWedge[1], X[a+N] );
977 LvArray::tensorOps::copy< 3 >( XWedge[2], X[a+1] );
978 LvArray::tensorOps::copy< 3 >( XWedge[3], X[a+1+N] );
979 LvArray::tensorOps::copy< 3 >( XWedge[4], XGBot );
980 LvArray::tensorOps::copy< 3 >( XWedge[5], XGTop );
981 result = result + computationalGeometry::elementVolume< finiteElement::H1_Wedge_Lagrange1_Gauss6 >( XWedge );
983 LvArray::tensorOps::copy< 3 >( XWedge[0], X[N-1] );
984 LvArray::tensorOps::copy< 3 >( XWedge[1], X[2*N-1] );
985 LvArray::tensorOps::copy< 3 >( XWedge[2], X[0] );
986 LvArray::tensorOps::copy< 3 >( XWedge[3], X[N] );
987 LvArray::tensorOps::copy< 3 >( XWedge[4], XGBot );
988 LvArray::tensorOps::copy< 3 >( XWedge[5], XGTop );
989 result = result + computationalGeometry::elementVolume< finiteElement::H1_Wedge_Lagrange1_Gauss6 >( XWedge );
GEOS_HOST_DEVICE int 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.
GEOS_HOST_DEVICE bool 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...
GEOS_HOST_DEVICE void FixNormalOrientation_3D(NORMAL_TYPE &&normal)
Change the orientation of the input vector to be consistent in a global sense.
GEOS_HOST_DEVICE real64 prismVolume(real64 const (&X)[2 *N][3])
Compute the volume of a prism with N-sided polygon base.
GEOS_HOST_DEVICE bool 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)
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 centroid_3DPolygon(arraySlice1d< localIndex const > const pointsIndices, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const &points, CENTER_TYPE &¢er, 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.
GEOS_HOST_DEVICE void RotationMatrix_3D(NORMAL_TYPE const &normal, MATRIX_TYPE &&rotationMatrix)
Calculate the rotation matrix for a face in the 3D space.
GEOS_HOST_DEVICE bool 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.
GEOS_HOST_DEVICE int 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 chos...
constexpr real64 machinePrecision
Machine epsilon for double-precision calculations.
GEOS_HOST_DEVICE real64 elementVolume(real64 const (&X)[FE_TYPE::numNodes][3])
Compute the volume of an element (tetrahedron, pyramid, wedge, hexahedron)
GEOS_HOST_DEVICE int 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 ...
GEOS_HOST_DEVICE real64 hexahedronVolume(real64 const (&X)[8][3])
Compute the volume of an hexahedron.
array1d< int > orderPointsCCW(arrayView2d< real64 > const &points, NORMAL_TYPE const &normal)
Reorder a set of points counter-clockwise.
GEOS_HOST_DEVICE int 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.
GEOS_HOST_DEVICE real64 tetrahedronVolume(real64 const (&X)[4][3])
Compute the volume of an tetrahedron.
void 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.
GEOS_HOST_DEVICE int 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 fro...
real64 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.
GEOS_HOST_DEVICE void 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 ...
GEOS_HOST_DEVICE real64 wedgeVolume(real64 const (&X)[6][3])
Compute the volume of a wedge.
GEOS_HOST_DEVICE real64 pyramidVolume(real64 const (&X)[5][3])
Compute the volume of a pyramid.
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 computeDiameter(POINT_COORDS_TYPE points, localIndex const &numPoints)
Calculate the diameter of a set of points in a given dimension.
GEOS_HOST_DEVICE GEOS_FORCE_INLINE int sign(T const val)
Return the sign of a given value as an integer.
GEOS_HOST_DEVICE int 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.
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
#define GEOS_LOG_RANK(msg)
Log a message to the rank output stream.
#define GEOS_ERROR_IF_LT(lhs, rhs)
Raise a hard error if one value compares less than the other.
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
#define GEOS_ERROR_IF(EXP, msg)
Conditionally raise a hard error and terminate the program.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
LvArray::ArrayOfArraysView< T, INDEX_TYPE const, CONST_SIZES, LvArray::ChaiBuffer > ArrayOfArraysView
View of array of variable-sized arrays. See LvArray::ArrayOfArraysView for details.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
std::int32_t integer
Signed integer type.
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Array< T, 1 > array1d
Alias for 1D array.