GemaCoreLib
The GeMA Core library
gmGeometryUtils.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
24 #ifndef _GEMA_GEOMETRY_UTILS_H_
25 #define _GEMA_GEOMETRY_UTILS_H_
26 
27 #include "gmVector.h"
28 #include "gmMatrix.h"
29 
30 class GmCellGeometry;
31 
35 namespace GmGeometryUtils
36 {
37  GMC_API_EXPORT double edgeLength (const GmVector& p1, const GmVector& p2);
38  GMC_API_EXPORT double triangleArea(const GmVector& p1, const GmVector& p2, const GmVector& p3);
39  GMC_API_EXPORT double quadArea (const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4);
40  GMC_API_EXPORT double hexahedronVolume(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
41  const GmVector& p5, const GmVector& p6, const GmVector& p7, const GmVector& p8);
42  GMC_API_EXPORT double tetrahedronVolume(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4);
43  GMC_API_EXPORT double wedgeVolume(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
44  const GmVector& p5, const GmVector& p6);
45 
46  GMC_API_EXPORT double pyramidVolume(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
47  const GmVector& p5);
48 
49  GMC_API_EXPORT void edgeCentroid (const GmVector& p1, const GmVector& p2, GmVector& coord);
50  GMC_API_EXPORT void triangleCentroid(const GmVector& p1, const GmVector& p2, const GmVector& p3, GmVector& coord);
51  GMC_API_EXPORT void quadCentroid (const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4, GmVector& coord);
52  GMC_API_EXPORT void hexahedronCentroid(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
53  const GmVector& p5, const GmVector& p6, const GmVector& p7, const GmVector& p8, GmVector& coord);
54  GMC_API_EXPORT void tetrahedronCentroid(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4, GmVector& coord);
55  GMC_API_EXPORT void wedgeCentroid(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
56  const GmVector& p5, const GmVector& p6, GmVector& coord);
57  GMC_API_EXPORT void pyramidCentroid(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4,
58  const GmVector& p5, GmVector& coord);
59 
60  GMC_API_EXPORT bool tri3DIsColinear(const GmMatrix& X, double eps = 1e-3);
61  GMC_API_EXPORT bool polygon3DIsPlanar(const GmMatrix& X, double tol);
62  GMC_API_EXPORT bool polygon2DIsCCW(const GmMatrix& X);
64  GMC_API_EXPORT bool polygonIsConvex (const GmMatrix& X);
65  GMC_API_EXPORT bool pointInConvexPolygon (const GmMatrix& X, const GmVector& p);
68  GMC_API_EXPORT bool pointInPolygon (const GmMatrix& X, const GmCellGeometry& g, const GmVector& p);
69  GMC_API_EXPORT bool pointInTetrahedron (const GmMatrix& X, const GmVector& p);
70  GMC_API_EXPORT bool pointInHexahedron (const GmMatrix& X, const GmVector& p);
71  GMC_API_EXPORT double pointToConvexPolygon (const GmMatrix& X, const GmVector& p, GmVector& proj);
72  GMC_API_EXPORT double pointToTetrahedron (const GmMatrix& X, const GmVector& p, GmVector& proj);
73 
74  GMC_API_EXPORT GmMatrix shuffle(const GmMatrix& X, const int order[], int size);
75 
76  GMC_API_EXPORT bool hex27ValidFacesCenter(const GmMatrix& X, double tol);
77  GMC_API_EXPORT bool hex20ValidEdges(const GmMatrix& X, double proportion_max = 0.5);
78  GMC_API_EXPORT bool hex8PlanarFaces(const GmMatrix& X, double tol);
79  GMC_API_EXPORT bool hex8IsWellDefined(const GmMatrix& X);
80  GMC_API_EXPORT bool tet10ValidEdges(const GmMatrix& X, double proportion_max = 0.5);
81  GMC_API_EXPORT bool tet4IsWellDefined(const GmVector& p1, const GmVector& p2, const GmVector& p3, const GmVector& p4);
82  GMC_API_EXPORT bool pyra5IsWellDefined(const GmMatrix& X, double tol);
83  GMC_API_EXPORT bool pyra13ValidEdges(const GmMatrix& X, double proportion_max = 0.5);
84  GMC_API_EXPORT bool wedge6IsWellDefined(const GmMatrix& X, double tol);
85  GMC_API_EXPORT bool wedge15ValidEdges(const GmMatrix& X, double proportion_max = 0.5);
86 
87  GMC_API_EXPORT bool isCoplanar3D(const GmVector& a, const GmVector& b, const GmVector& c, const GmVector& d, double tol);
90  GMC_API_EXPORT QVector<QVector<int>> edgeFaceList(const QVector<QVector<int>>& faceEdgeList, int edgeCount);
91 }
92 
93 #endif
94 
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.
Definition: gmGeometryUtils.cpp:269
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),...
Definition: gmGeometryUtils.cpp:645
bool polygon3DIsPlanar(const GmMatrix &X, double tol)
Returns true if the polygon delimited by the cell vertices is planar.
Definition: gmGeometryUtils.cpp:587
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 ...
Definition: gmGeometryUtils.cpp:933
void triangleCentroid(const GmVector &p1, const GmVector &p2, const GmVector &p3, GmVector &coord)
Returns the cartesian coordinate of the triangle centroid.
Definition: gmGeometryUtils.cpp:293
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.
Definition: gmGeometryUtils.cpp:419
GmMatrix shuffle(const GmMatrix &X, const int order[], int size)
Returns quad8 vertices in circular order.
Definition: gmGeometryUtils.cpp:612
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.
Definition: gmGeometryUtils.cpp:196
Groups utilitary routines for working with geometry.
Definition: gmGeometryUtils.cpp:32
static bool pointInTetrahedron(const GmMatrix &X, const GmVector &p, QVector< bool > &ss)
Helper function for pointInTetrahedron(X, p) also returning the "same side" vector.
Definition: gmGeometryUtils.cpp:784
Declaration of the GmMatrix class.
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)
Definition: gmGeometryUtils.cpp:325
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) a...
Definition: gmGeometryUtils.cpp:1044
GmMatrix get2DProjectionMatRow(const GmMatrix &X)
Returns a 3x2 (transposed) rigid transform matrix taking v1-origin and v2-origin 3D vectors to the 2D...
Definition: gmGeometryUtils.cpp:1282
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 co...
Definition: gmGeometryUtils.cpp:1095
A class used to return static metadata information about a cell geometry, along with some methods for...
Definition: gmCellGeometry.h:54
void edgeCentroid(const GmVector &p1, const GmVector &p2, GmVector &coord)
Returns the cartesian coordinate of the edge centroid.
Definition: gmGeometryUtils.cpp:283
bool polygonIsConvex(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is convex, or false otherwise.
Definition: gmGeometryUtils.cpp:559
GmMatrix get2DProjectionMatCol(const GmMatrix &X)
Returns a 2x3 rigid transform matrix taking v1-origin and v2-origin 3D vectors to the 2D xy plane.
Definition: gmGeometryUtils.cpp:1261
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.
Definition: gmGeometryUtils.cpp:240
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.
Definition: gmGeometryUtils.cpp:203
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.
Definition: gmGeometryUtils.cpp:1306
bool pointInTesselatedPolygon_Q8(const GmMatrix &X, const GmVector &p)
Returns true if the point p is inside the tesselated Q8 polygon or false otherwise.
Definition: gmGeometryUtils.cpp:700
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.
Definition: gmGeometryUtils.cpp:254
bool polygon2DIsCCW(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is oriented CCW.
Definition: gmGeometryUtils.cpp:527
void tetrahedronCentroid(const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, GmVector &coord)
Returns the cartesian coordinate of the tetrahedron centroid.
Definition: gmGeometryUtils.cpp:373
void quadCentroid(const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, GmVector &coord)
Returns the cartesian coordinate of the quad centroid.
Definition: gmGeometryUtils.cpp:304
Declaration of the GmVector class.
#define GMC_API_EXPORT
Macro for controling if the class is being exported (GEMA_CORE_LIB defined) or imported (GEMA_CORE_LI...
Definition: gmCoreConfig.h:35
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),...
Definition: gmGeometryUtils.cpp:165
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.
Definition: gmGeometryUtils.cpp:383
double edgeLength(const GmVector &p1, const GmVector &p2)
Calculates the length of the edge determined by points p1 and p2.
Definition: gmGeometryUtils.cpp:157
bool pointInTesselatedPolygon_T6(const GmMatrix &X, const GmVector &p)
Returns true if the point p is inside the tesselated T6 polygon or false otherwise.
Definition: gmGeometryUtils.cpp:666
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.
Definition: gmGeometryUtils.cpp:958
bool polygon2DIsSelfIntersecting(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is self-intersecting, or false otherwise.
Definition: gmGeometryUtils.cpp:498
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
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.
Definition: gmGeometryUtils.cpp:745