![]() |
GemaCoreLib
The GeMA Core library
|
Groups utilitary routines for working with geometry. More...
Functions | |
static double | angleScore2D (const GmVector &a, const GmVector &b, const GmVector &c) |
Returns the angle score of the ca and cb vectors in 2D. If score is: More... | |
static double | orientScore2D (const GmVector &a, const GmVector &b, const GmVector &c) |
Returns the orientation score of the a, b and c points in 2D. Returns a positive value if the points a, b and c are arranged in CCW order, a negative value if the points are in CW order, and zero if the points are collinear. A more common (but less symmetric) interpretation is that it returns a positive value if a lies to the left of the directed line bc (see https://people.eecs.berkeley.edu/~jrs/meshpapers/robnotes.pdf). Can be used in 3D but only x and y coordinates will be considered. | |
static double | orientScore3D (const GmVector &a, const GmVector &b, const GmVector &c, const GmVector &d) |
Given four points a, b, c, and d, returns the signed volume of the parallelepiped determined by the vectors ad, bd and cd. It is positive if the points occur in the orientation illustrated in Figure 3 (see https://people.eecs.berkeley.edu/~jrs/meshpapers/robnotes.pdf), negative if they occur in the mirror-image orientation, and zero if the four points are coplanar. You can apply a right-hand rule : orient your right hand with fingers curled to follow the circular sequence bcd. If your thumb points toward a, a positive value is returned. | |
static bool | sameSide3D (const GmVector &p, const GmVector &n, const GmVector &a, const GmVector &b) |
Returns true if points a and b are on the same side of plane defined by point p and normal n. | |
static bool | sameSide3D (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &a, const GmVector &b) |
Returns true if points a and b are on the same side of plane defined by points p1, p2 and p3. | |
static GmMatrix | rotationMatrix (const GmVector &u, double angle) |
double | edgeLength (const GmVector &p1, const GmVector &p2) |
Calculates the length of the edge determined by points p1 and p2. | |
double | triangleArea (const GmVector &p1, const GmVector &p2, const GmVector &p3) |
Calculates the area of the triangle determined by points p1, p2 and p3 (2D or 3D), organized in CCW order. | |
double | quadArea (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4) |
Calculates the area of a quad defined by points p1 to p4 (2D or 3D), defined in CCW order. | |
double | hexahedronVolume (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5, const GmVector &p6, const GmVector &p7, const GmVector &p8) |
Calculates the volume of the given hexahedron. Points follow the standard Hex cell definition. | |
double | tetrahedronVolume (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4) |
Calculates the volume of the given tetrahedron. Points follow the standard Tet cell definition. | |
double | wedgeVolume (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5, const GmVector &p6) |
Calculates the volume of the given wedge. Points follow the standard Wedge cell definition. | |
double | pyramidVolume (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5) |
Calculates the volume of the given pyramid. Points follow the standard pyramid cell definition. | |
void | edgeCentroid (const GmVector &p1, const GmVector &p2, GmVector &coord) |
Returns the cartesian coordinate of the edge centroid. | |
void | triangleCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, GmVector &coord) |
Returns the cartesian coordinate of the triangle centroid. | |
void | quadCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, GmVector &coord) |
Returns the cartesian coordinate of the quad centroid. | |
void | hexahedronCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5, const GmVector &p6, const GmVector &p7, const GmVector &p8, GmVector &coord) |
void | tetrahedronCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, GmVector &coord) |
Returns the cartesian coordinate of the tetrahedron centroid. | |
void | wedgeCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5, const GmVector &p6, GmVector &coord) |
Returns the cartesian coordinate of the wedge centroid. | |
void | pyramidCentroid (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, const GmVector &p5, GmVector &coord) |
Returns the cartesian coordinate of the pyramid centroid. | |
double | intersectionLength (double a0, double a1, double b0, double b1) |
Intersection length between two intervals Positive for intervals intersecting, negative otherwise. | |
bool | lineIntersection (GmVector &A, GmVector &B, GmVector &C, GmVector &D, double eps=1e-7) |
Returns true if the lines AB and CD intersect Assumes A, B, C and D are 2D, doesnt cosider the segments extreme points for intersection. | |
bool | polygon2DIsSelfIntersecting (const GmMatrix &X) |
Returns true if the polygon delimited by the cell vertices is self-intersecting, or false otherwise. More... | |
bool | polygon2DIsCCW (const GmMatrix &X) |
Returns true if the polygon delimited by the cell vertices is oriented CCW. More... | |
bool | polygonIsConvex (const GmMatrix &X) |
Returns true if the polygon delimited by the cell vertices is convex, or false otherwise. More... | |
bool | polygon3DIsPlanar (const GmMatrix &X, double tol) |
Returns true if the polygon delimited by the cell vertices is planar. More... | |
GmMatrix | shuffle (const GmMatrix &X, const int order[], int size) |
Returns quad8 vertices in circular order. More... | |
bool | tri3DIsColinear (const GmMatrix &X, double eps) |
bool | pointInConvexPolygon (const GmMatrix &X, const GmVector &p) |
Returns true if the point p is inside the polygon delimited by the cell vertices (convex hull), or false otherwise. More... | |
bool | pointInTesselatedPolygon_T6 (const GmMatrix &X, const GmVector &p) |
Returns true if the point p is inside the tesselated T6 polygon or false otherwise. More... | |
bool | pointInTesselatedPolygon_Q8 (const GmMatrix &X, const GmVector &p) |
Returns true if the point p is inside the tesselated Q8 polygon or false otherwise. More... | |
bool | pointInPolygon (const GmMatrix &X, const GmCellGeometry &g, const GmVector &p) |
Returns true if the point p is inside the polygon delimited by the cell vertices, or false otherwise. More... | |
static bool | pointInTetrahedron (const GmMatrix &X, const GmVector &p, QVector< bool > &ss) |
Helper function for pointInTetrahedron(X, p) also returning the "same side" vector. | |
bool | pointInTetrahedron (const GmMatrix &X, const GmVector &p) |
Returns true if the point p is inside the tetrahedron, or false otherwise. The X matrix should be the node matrix for a GeMA TET4 element with nodes in columns. | |
bool | pointInHexahedron (const GmMatrix &X, const GmVector &p) |
bool | isTriNormalOutward (const GmVector &p0, const GmVector &p1, const GmVector &p2, const GmVector &ref) |
bool | isQuadNormalOutward (const GmVector &p0, const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &ref) |
bool | hex8PlanarFaces (const GmMatrix &X, double tol) |
bool | hex8IsWellDefined (const GmMatrix &X) |
void | lineProjection (const GmVector &a, const GmVector &b, const GmVector &p, double *t, double *u) |
template<int N> | |
bool | checkQuadraticEdges (const GmMatrix &X, double proportion_max, int edges[N][3]) |
bool | tet10ValidEdges (const GmMatrix &X, double proportion_max) |
bool | hex20ValidEdges (const GmMatrix &X, double proportion_max) |
Returns true if the points defined between the edges of hex20 follows: 1) Their projection to the 3d line containing the edge is within the edge, i.e., their parametrization is in (0,1) 2) The distance to the projection is not greater than a proportion of the edge length. | |
bool | hex27ValidFacesCenter (const GmMatrix &X, double tol) |
Returns true if face center point is collinear with the face. Assumes that the quad faces are planar. | |
bool | pyra5IsWellDefined (const GmMatrix &X, double tol) |
bool | pyra13ValidEdges (const GmMatrix &X, double proportion_max) |
bool | wedge6IsWellDefined (const GmMatrix &X, double tol) |
bool | wedge15ValidEdges (const GmMatrix &X, double proportion_max) |
double | pointToConvexPolygon (const GmMatrix &X, const GmVector &p, GmVector &proj) |
Returns the dquared distance of point p to the polygon delimited by the cell vertices (convex hull) and fills the coordinates of the corresponding projected point on cell. More... | |
double | pointToTetrahedron (const GmMatrix &X, const GmVector &p, GmVector &proj) |
Returns the squared distance of point p to the tetrahedron delimited by its vertices and fills the coordinates of the corresponding projected point on tetrahedron. More... | |
bool | tet4IsWellDefined (const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4) |
bool | isCoplanar3D (const GmVector &a, const GmVector &b, const GmVector &c, const GmVector &d, double tol) |
void | get2DPlaneRotationVectors (const GmVector &origin, const GmVector &v1, const GmVector &v2, GmVector &ix, GmVector &iy) |
GmMatrix | get2DProjectionMatCol (const GmMatrix &X) |
Returns a 2x3 rigid transform matrix taking v1-origin and v2-origin 3D vectors to the 2D xy plane. | |
GmMatrix | get2DProjectionMatRow (const GmMatrix &X) |
Returns a 3x2 (transposed) rigid transform matrix taking v1-origin and v2-origin 3D vectors to the 2D xy plane. | |
QVector< QVector< int > > | edgeFaceList (const QVector< QVector< int >> &faceEdgeList, int edgeCount) |
Returns a edgeFace Vector of faces indexed by edge Its input is the inverse map faceEdge Vector. | |
Variables | |
static int | hex2tet [6][4] |
Returns true if the point p is inside the hexahedron, or false otherwise. The X matrix should be the node matrix for a GeMA HEX8 element with nodes in columns. More... | |
Groups utilitary routines for working with geometry.
|
static |
Returns the angle score of the ca and cb vectors in 2D. If score is:
GMC_API_EXPORT void GmGeometryUtils::hexahedronCentroid | ( | const GmVector & | p1, |
const GmVector & | p2, | ||
const GmVector & | p3, | ||
const GmVector & | p4, | ||
const GmVector & | p5, | ||
const GmVector & | p6, | ||
const GmVector & | p7, | ||
const GmVector & | p8, | ||
GmVector & | coord | ||
) |
Returns the cartesian coordinate of the hexahedron centroid fails for non-convex
GMC_API_EXPORT bool GmGeometryUtils::pointInConvexPolygon | ( | const GmMatrix & | X, |
const GmVector & | p | ||
) |
Returns true if the point p is inside the polygon delimited by the cell vertices (convex hull), or false otherwise.
The X matrix should be 2 x n where n is the number of polygon vertices and it should provide the nodes in a "circular" order. Keep in mind that this is NOT true for the standard GeMA node matrix for QUADRATIC elements.
GMC_API_EXPORT bool GmGeometryUtils::pointInPolygon | ( | const GmMatrix & | X, |
const GmCellGeometry & | g, | ||
const GmVector & | p | ||
) |
Returns true if the point p is inside the polygon delimited by the cell vertices, or false otherwise.
Based on point-in-polygon algorithm, see http://alienryderflex.com/polygon/ for more details. If the node is exactly on the edge of the polygon, then the function may return true or false. Note that division by zero is avoided because the division is protected by the "if" clause which surrounds it.
The X matrix should be the node matrix for a GeMA surfcae element with nodes in columns. The cell geometry object is used to traverse the polygon nodes in a "circular" order, even for quadratic elements.
GMC_API_EXPORT bool GmGeometryUtils::pointInTesselatedPolygon_Q8 | ( | const GmMatrix & | X, |
const GmVector & | p | ||
) |
Returns true if the point p is inside the tesselated Q8 polygon or false otherwise.
The Quad8 surface is approximated by converting the Quad8 into 4 quad4. The X matrix should be the node matrix for a GeMA QUAD8 element with nodes in columns.
GMC_API_EXPORT bool GmGeometryUtils::pointInTesselatedPolygon_T6 | ( | const GmMatrix & | X, |
const GmVector & | p | ||
) |
Returns true if the point p is inside the tesselated T6 polygon or false otherwise.
The Tri6 surface is approximated by converting the TRI6 into 4 triangles. The X matrix should be the node matrix for a GeMA TRI6 element with nodes in columns.
GMC_API_EXPORT double GmGeometryUtils::pointToConvexPolygon | ( | const GmMatrix & | X, |
const GmVector & | p, | ||
GmVector & | proj | ||
) |
Returns the dquared distance of point p to the polygon delimited by the cell vertices (convex hull) and fills the coordinates of the corresponding projected point on cell.
Assume p is outside of the polygon. The X matrix should be 2 x n where n is the number of polygon vertices and it should provide the nodes in a "circular" order. Keep in mind that this is NOT true for the standard GeMA node matrix for QUADRATIC elements.
GMC_API_EXPORT double GmGeometryUtils::pointToTetrahedron | ( | const GmMatrix & | X, |
const GmVector & | p, | ||
GmVector & | proj | ||
) |
Returns the squared distance of point p to the tetrahedron delimited by its vertices and fills the coordinates of the corresponding projected point on tetrahedron.
Assume p is outside of the tetrahedron. The X matrix should be the node matrix for a GeMA TET4 element with nodes in columns.
IMPORTANT (by cmendes): Does this function work? Was it even tested? I don't think so.... :(
GMC_API_EXPORT bool GmGeometryUtils::polygon2DIsCCW | ( | const GmMatrix & | X | ) |
Returns true if the polygon delimited by the cell vertices is oriented CCW.
The X matrix should be 2 x n where n is the number of polygon vertices.
GMC_API_EXPORT bool GmGeometryUtils::polygon2DIsSelfIntersecting | ( | const GmMatrix & | X | ) |
Returns true if the polygon delimited by the cell vertices is self-intersecting, or false otherwise.
The X matrix should be 2 x n where n is the number of polygon vertices and it should provide the nodes in a "circular" order. Keep in mind that this is NOT true for the standard GeMA node matrix for QUADRATIC elements.
GMC_API_EXPORT bool GmGeometryUtils::polygon3DIsPlanar | ( | const GmMatrix & | X, |
double | tol | ||
) |
Returns true if the polygon delimited by the cell vertices is planar.
The X matrix should be 3 x n where n is the number of polygon vertices. The eps is the floating point comparison tolerance, default is 1e-7.
We could probably speed up it by re-using the cross product of each vertice but the function 'isCoplanar3D' is already implemented, tested and uses a formula with smaller fp error.
GMC_API_EXPORT bool GmGeometryUtils::polygonIsConvex | ( | const GmMatrix & | X | ) |
Returns true if the polygon delimited by the cell vertices is convex, or false otherwise.
The X matrix should be 2 x n where n is the number of polygon vertices and it should provide the nodes in a "circular" order. Keep in mind that this is NOT true for the standard GeMA node matrix for QUADRATIC elements.
GMC_API_EXPORT GmMatrix GmGeometryUtils::shuffle | ( | const GmMatrix & | X, |
const int | order[], | ||
int | size | ||
) |
Returns quad8 vertices in circular order.
The X matrix should be n x 8 where n is the number of dimensions. The original ordering of quad8 vertices can be found in gmQuadCellGeometryInfo.cpp
|
static |
Returns true if the point p is inside the hexahedron, or false otherwise. The X matrix should be the node matrix for a GeMA HEX8 element with nodes in columns.