Woolz Image Processing
Version 1.8.3
|
Geometric utility functions. More...
Data Structures | |
struct | _WlzGeomPolyListItem2D |
Item for polygon vertex list. See WlzGeomPolyTriangulate2D(). More... | |
Typedefs | |
typedef struct _WlzGeomPolyListItem2D | WlzGeomPolyListItem2D |
Functions | |
double | cbrt (double c) |
WlzDVertex2 | WlzGeomTriangleCen2D (WlzDVertex2 v0, WlzDVertex2 v1, WlzDVertex2 v2) |
Computes the position of the centroid of the given triangle in 2D. More... | |
WlzDVertex3 | WlzGeomTriangleCen3D (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2) |
Computes the position of the centroid of the given triangle in 3D. More... | |
WlzDVertex3 | WlzGeomTetrahedronCen3D (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, WlzDVertex3 v3) |
Computes the position of the centroid of the given tetrahedron in 3D. More... | |
int | WlzGeomTriangleCircumcentre (WlzDVertex2 *ccVx, WlzDVertex2 vx0, WlzDVertex2 vx1, WlzDVertex2 vx2) |
Computes the circumcentre of the given triangle. More... | |
int | WlzGeomVxInTriangle2D (WlzDVertex2 p0, WlzDVertex2 p1, WlzDVertex2 p2, WlzDVertex2 pP) |
Tests whether the given vertex lies within the given triangle using a barycentric coordinates test. More... | |
int | WlzGeomVxInTriangle3D (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, WlzDVertex3 vQ, double vPMax) |
First finds the closest point on the plane of the triangle to the given point. Then if the distance from the point to the plane is less than the given tolerance vvalue tests to set if the given vertex lies within the given triangle using a barycentric coordinates test (see WlzGeomVxInTriangle2D()). More... | |
int | WlzGeomVxInTetrahedron (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, WlzDVertex3 v3, WlzDVertex3 vP) |
Tests whether the given vertex lies within the given tetrahedron using a barycentric coordinates test. More... | |
WlzLong | WlzGeomTriangleSnArea2I (WlzIVertex2 vx0, WlzIVertex2 vx1, WlzIVertex2 vx2) |
Computes twice the signed area of the given triangle. More... | |
double | WlzGeomTriangleSnArea2 (WlzDVertex2 vx0, WlzDVertex2 vx1, WlzDVertex2 vx2) |
Computes twice the signed area of the given triangle. More... | |
double | WlzGeomTetraSnVolume6 (WlzDVertex3 vx0, WlzDVertex3 vx1, WlzDVertex3 vx2, WlzDVertex3 vx3) |
Computes six times the signed volume of the given tetrahedron. More... | |
int | WlzGeomTetVolZeroI (WlzIVertex3 v0, WlzIVertex3 v1, WlzIVertex3 v2, WlzIVertex3 v3) |
Checks whether the volume of the tetrahedron with the given interger vertices has zero volume. More... | |
int | WlzGeomTetVolZeroD (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, WlzDVertex3 v3) |
Checks whether the volume of the tetrahedron with the given double vertices has zero volume. More... | |
WlzLong | WlzGeomTriangleArea2Sq3I (WlzIVertex3 vx0, WlzIVertex3 vx1, WlzIVertex3 vx2) |
Computes twice the square of the area of the given 3D triangle. More... | |
double | WlzGeomTriangleArea2Sq3 (WlzDVertex3 vx0, WlzDVertex3 vx1, WlzDVertex3 vx2) |
Computes twice the square of the area of the given 3D triangle. More... | |
int | WlzGeomInTriangleCircumcircle (WlzDVertex2 vx0, WlzDVertex2 vx1, WlzDVertex2 vx2, WlzDVertex2 gVx) |
Tests to see if the given vertex is inside the circumcircle of the given triangle. More... | |
int | WlzGeomLineSegmentsIntersect (WlzDVertex2 p0, WlzDVertex2 p1, WlzDVertex2 q0, WlzDVertex2 q1, WlzDVertex2 *dstN) |
Tests to see if the two given line segments intersect. More... | |
int | WlzGeomCmpAngle (WlzDVertex2 p0, WlzDVertex2 p1) |
Given two end connected 2D line segments this function compares the CCW angle of the segments. More... | |
int | WlzGeomVtxEqual2D (WlzDVertex2 pos0, WlzDVertex2 pos1, double tolSq) |
Checks to see if two verticies are the same within some tollerance. More... | |
int | WlzGeomVtxEqual3D (WlzDVertex3 pos0, WlzDVertex3 pos1, double tol) |
Checks to see if two verticies are the same within some tollerance. More... | |
void | WlzGeomVtxSortRadial (int nV, WlzDVertex3 *vP, int *idxBuf, WlzDVertex2 *wP, WlzDVertex3 rV) |
Sorts the given 3D verticies, which lie in a plane perpendicular to the radial vector, in order of their angle the radial vector. More... | |
WlzDVertex3 | WlzGeomTriangleNormal (WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2) |
Computes the unit normal vector perpendicular to the triangle \(v_0, v_1, v_2\). More... | |
int | WlzGeomTriangleAABBIntersect2D (WlzDVertex2 t0, WlzDVertex2 t1, WlzDVertex2 t2, WlzDVertex2 b0, WlzDVertex2 b1, int tst) |
Tests for an intersection between the given triangle and the axis aligned bounding box using the Separating Axis Theorem (SAT). More... | |
int | WlzGeomTetrahedronAABBIntersect3D (WlzDVertex3 t0, WlzDVertex3 t1, WlzDVertex3 t2, WlzDVertex3 t3, WlzDVertex3 b0, WlzDVertex3 b1, int tst) |
Tests for an intersection between the given tetrahedron and the axis aligned bounding box using the Separating Axis Theorem (SAT). More... | |
int | WlzGeomPlaneAABBIntersect (double a, double b, double c, double d, WlzDBox3 box) |
Tests for an intersection between the plane defined by the equation: \(ax + by + cz + d = 0\) and the given axis aligned bounding box. More... | |
int | WlzGeomPlaneLineIntersect (double a, double b, double c, double d, WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 *dstIsn) |
Tests for an intersection between the plane defined by the equation: \(ax + by + cz + d = 0\) and the line segment with end points \(p_0\) and \(p_1\). More... | |
int | WlzGeomPlaneTriangleIntersect (double a, double b, double c, double d, WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 p2, WlzDVertex3 *dstIsn0, WlzDVertex3 *dstIsn1) |
Tests for an intersection between a plane and a triangle. More... | |
double | WlzGeomEllipseVxDistSq (WlzDVertex2 centre, WlzDVertex2 sAx, WlzDVertex2 gPnt) |
Given an ellipse defined by it's centre \(\mathbf{c}\) and it's semi axes \(\mathbf{a}\). This function computes the square of the ratio of the distances from the centre of the ellipse to the given point and from the centre of the ellipse in the direction of the given point to the ellipse. More... | |
unsigned int | WlzGeomHashVtx3D (WlzDVertex3 pos, double tol) |
Computes a hash value from a given 3D double precision position. More... | |
unsigned int | WlzGeomHashVtx2D (WlzDVertex2 pos, double tol) |
Computes a hash value from a given 2D double precision position. More... | |
int | WlzGeomCmpVtx3D (WlzDVertex3 pos0, WlzDVertex3 pos1, double tol) |
Compares the coordinates of the given 3D double precision vertices to find a signed value for sorting. More... | |
int | WlzGeomCmpVtx2D (WlzDVertex2 pos0, WlzDVertex2 pos1, double tol) |
Compares the coordinates of the given 2D double precision vertices to find a signed value for sorting. More... | |
WlzDVertex2 | WlzGeomUnitVector2D (WlzDVertex2 vec) |
Computes the 2D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\). More... | |
WlzDVertex3 | WlzGeomUnitVector3D (WlzDVertex3 vec) |
Computes the 3D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\). More... | |
WlzDVertex2 | WlzGeomUnitVector2D2 (WlzDVertex2 v1, WlzDVertex2 v0) |
Computes the unit 2D vector with the direction given by \(\mathbf{p}_1 - \mathbf{p}_0\). If the two given vertices are coincident then a zero vector is returned instead of a unit vector. More... | |
WlzDVertex3 | WlzGeomUnitVector3D2 (WlzDVertex3 v1, WlzDVertex3 v0) |
Computes the unit 3D vector with the direction given by \(\mathbf{p}_1 - \mathbf{p}_0\). If the two given vertices are coincident then a zero vector is returned instead of a unit vector. More... | |
int | WlzGeomVertexInDiamCircle (WlzDVertex2 lPos0, WlzDVertex2 lPos1, WlzDVertex2 pos) |
Determines whether a point is inside the diametral circle of a line segment. If two vectors \(\mathbf{v_0}\) and \(\mathbf{v_1}\) are directed from the line segment end points to the point, then the angle between the vectors is \(<\) \(90^\circ\) if the point is inside the diametral circle, \(0\) if it lies on the circle and \(>\) \(90^{\circ}\) if it lies outside the circle. This is easily tested by \(\mathbf{v_0} \cdot \mathbf{v_1} < 0\). More... | |
int | WlzGeomItrSpiralRing (int step) |
Computes the ring of a spiral. If two rings differ by more than one then at least one itteration outwards on the spiral has been performed between the rings. More... | |
int | WlzGeomItrSpiral2I (int step, int *pX, int *pY) |
Iterates the given positions coordinates through a 2D expanding integer spiral. More... | |
int | WlzGeomItrSpiralShell (int step) |
Computes the shell of a spiral. If two shells differ by more than one then at least one itteration outwards on the spiral has been performed between the shells. More... | |
int | WlzGeomItrSpiral3I (int step, int *pX, int *pY, int *pZ) |
Iterates the given positions coordinates through a 3D expanding integer spiral. More... | |
double | WlzGeomDistSq2D (WlzDVertex2 v0, WlzDVertex2 v1) |
Computes square of the Euclidean distance between the given two vertices. More... | |
double | WlzGeomDistSq3D (WlzDVertex3 v0, WlzDVertex3 v1) |
Computes square of the Euclidean distance between the given two vertices. More... | |
double | WlzGeomDist2D (WlzDVertex2 v0, WlzDVertex2 v1) |
Computes the Euclidean distance between the given two vertices. More... | |
double | WlzGeomDist3D (WlzDVertex3 v0, WlzDVertex3 v1) |
Computes the Euclidean distance between the given two vertices. More... | |
int | WlzGeomTriangleAffineSolve (double *xTr, double *yTr, double dd, WlzDVertex2 *sVx, WlzDVertex2 *dVx, double thresh) |
If the unsigned area of the triangle is very small then the only the transform translation coefficients are computed with the other coefficients being set to zero. If the unsigned area of the triangle is not very small then a system of linear equations is solved for the coefficients of the 2D affine transform from the source triangle to the destination triangle. More... | |
int | WlzGeomTetraAffineSolveLU (double *tr, WlzDVertex3 *sVx, WlzDVertex3 *dVx) |
Computes the affine transform coefficients from the source to target tetrahedron. More... | |
int | WlzGeomTetraAffineSolve (double *tr, WlzDVertex3 *sVx, WlzDVertex3 *dVx, double thresh) |
Computes the affine transform coefficients from the source to target tetrahedron. More... | |
WlzDVertex2 | WlzGeomObjLineSegIntersect2D (WlzObject *obj, WlzDVertex2 p0, WlzDVertex2 p1, double tol, int inside, int method, int *dstStat) |
Given a Woolz object and two vertices, finds the position along a line segment between the two vertices which is just inside/outside the boundary of the object. The destination pointer is used to return the status of the vertices, using the following code: 0 - One of the given verticies was inside and the other outside, 1 - Both the given verticies were inside, 2 - Both the given verticies were outside. This function assumes that the line segment only crosses the object's boundary once. More... | |
WlzDVertex3 | WlzGeomObjLineSegIntersect3D (WlzObject *obj, WlzDVertex3 p0, WlzDVertex3 p1, double tol, int inside, int method, int *dstStat) |
Given a Woolz object and two vertices, finds the position along a line segment between the two vertices which is just inside/outside the boundary of the object. The destination pointer is used to return the status of the vertices, using the following code: 0 - One of the given verticies was inside and the other outside, 1 - Both the given verticies were inside, 2 - Both the given verticies were outside. This function assumes that the line segment only crosses the object's boundary once. More... | |
double | WlzGeomTetraInSphereDiam (WlzDVertex3 vx0, WlzDVertex3 vx1, WlzDVertex3 vx2, WlzDVertex3 vx3) |
Given the coordinates of the four vertices of a tetrahedron the function computes the maximum diameter of an inscribed sphere. More... | |
double | WlzGeomTetraInSphereRegDiam (double side) |
Given the side length of a regular tetrahedron this function computes the maximum diameter of an inscribed sphere. More... | |
double | WlzGeomPolar2D (WlzDVertex2 org, WlzDVertex2 dst, double *dstRad) |
Computes the angle of a ray from the origin to the destination vertex along with the length of the ray. Angles are: More... | |
double | WlzGeomCos3V (WlzDVertex2 v0, WlzDVertex2 v1, WlzDVertex2 v2) |
Computes the cosine of angle between line segments (v0, v1) and (v1, v2). If any of these vertices are coincident then zero is returned. More... | |
int | WlzGeomVtxOnLineSegment2D (WlzDVertex2 tst, WlzDVertex2 seg0, WlzDVertex2 seg1, double tol) |
Tests whether the given test vertex is on the given line segment. If all three vertices are coincident then the test vertex is considered to line on the line segment. More... | |
int | WlzGeomVtxOnLineSegment3D (WlzDVertex3 pX, WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 *dstN) |
Tests whether the given test vertex is on the given line segment. If all three vertices are coincident then the test vertex is considered to be coincident with an end point on the line segment. More... | |
double | WlzGeomArcLength2D (WlzDVertex2 a, WlzDVertex2 b, WlzDVertex2 c) |
Computes the arc length from a to b traveling CCW on a circle with centre c. More... | |
int | WlzGeomRectFromWideLine (WlzDVertex2 s, WlzDVertex2 t, double w, WlzDVertex2 *v) |
Computes the coordinates of vertices that may be used to draw a wide rectangle. More... | |
WlzDVertex3 | WlzGeomLinePlaneIntersection (WlzDVertex3 v, WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 p2, WlzDVertex3 p3, int *dstPar) |
Computes the intersection of a line with a plane. More... | |
int | WlzGeomLineTriangleIntersect3D (WlzDVertex3 org, WlzDVertex3 dir, WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, int *dstPar, double *dstT, double *dstU, double *dstV) |
Tests whether a line directed from a given origin intersects a triangle in 3D space. This function is based on the algorithm: Tomas Moller and Ben Trumbore, "Fast, Minimum Storage Ray/Triangle Intersection", Journal of Graphics Tools, 1997(2), pp 25–30. More... | |
int | WlzGeomLineLineSegmentIntersect3D (WlzDVertex3 r0, WlzDVertex3 rD, WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 *dstN) |
Tests to see if the two given line segment is intersected by the given line using the ALG_DBL_TOLLERANCE tollerance value. The line is a line which passes through the given point to infinity (on both sides) with the given direction. More... | |
int | WlzGeomVtxOnLine3D (WlzDVertex3 p0, WlzDVertex3 r0, WlzDVertex3 rD) |
Tests whether a vertex is a line. More... | |
int | WlzGeomBaryCoordsTri2D (WlzDVertex2 p0, WlzDVertex2 p1, WlzDVertex2 p2, WlzDVertex2 pX, double *lambda) |
Given the coordinates of the vertices of a 2D triangle and a position within the triangle, this function computes the barycentric coordinates ( \(\lambda_0\), \(\lambda_0\), \(\lambda_2\) of the position. More... | |
double | WlzGeomInterpolateTri2D (WlzDVertex2 p0, WlzDVertex2 p1, WlzDVertex2 p2, double v0, double v1, double v2, WlzDVertex2 pX) |
Given the coordinates of the vertices of a 2D triangle and a set of values at each of these vertices, this function interpolates the value at the given position which is either inside or on an edge of the triangle. This is implimented using barycentric coordinates. Once the barycentric coordinates ( \(\lambda_0\), \(\lambda_0\), \(\lambda_2\)) have been computed then the interpolated value is given by: \(v = \sum_{i=0}^{2}{v_i \lambda_i}\). If the determinant is zero in solving for the barycentric coordinates then the interpolated value is just the mean of the given values. More... | |
double | WlzGeomInterpolatePoly2D (int n, WlzDVertex2 *p, double *v, double *w, WlzDVertex2 q, WlzErrorNum *dstErr) |
Given the vertex coordinates of an irregular possible non-convex 2D polygon ordered counter-clockwise and a set of values at each of these vertices, this function interpolates the value at the given position which must be inside the convex hull of the polygon, on an edge of the convex hull of the polygon or coincident with one of the vertices of it's convex hull. This function first calls WlzConvHullClarkson2D() to compute the convex hull and then WlzGeomInterpolateConvexPoly2D() to perform the interpolation. If the polygonis known to be convex then WlzGeomInterpolateConvexPoly2D() should be called directly since temporary workspaces are allocated. More... | |
double | WlzGeomInterpolateConvexPoly2D (int n, WlzDVertex2 *p, double *v, double *w, WlzDVertex2 q) |
Given the vertex coordinates of an irregular convex 2D polygon ordered counter-clockwise and a set of values at each of these vertices, this function interpolates the value at the given position which must be inside the polygon, on an edge of the polygon or coincident with one of it's vertices. This is implimented using general barycentric coordinates, see the paper: "Generalized Barycentric Coordinates on
Irregular Polygons" Mark Mayer, etal, Journal of Graphics Tools 2002. All parameters of this function must be valid. More... | |
int | WlzGeomBaryCoordsTet3D (WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 p2, WlzDVertex3 p3, WlzDVertex3 pX, double *lambda) |
Given the coordinates of the vertices of a 3D tetrahedron and a position within the tetrahedron, this function computes the barycentric coordinates ( \(\lambda_0\), \(\lambda_0\), \(\lambda_2\), \(\lambda_3\)) of the position. More... | |
double | WlzGeomInterpolateTet3D (WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 p2, WlzDVertex3 p3, double v0, double v1, double v2, double v3, WlzDVertex3 pX) |
Given the coordinates of the vertices of a 3D tetrahedron and a set of values at each of these vertices, this function interpolates the value at the given position which is either inside or on a face/edge of the tetrahedron. This is implimented using barycentric coordinates. Once the barycentric coordinates ( \(\lambda_0\), \(\lambda_0\), \(\lambda_2\), \(\lambda_3\)) have been computed then the interpolated value is given by: \(v = \sum_{i=0}^{3}{v_i \lambda_i}\). If the determinant is zero in solving for the barycentric coordinates then the interpolated value is just the mean of the given values. More... | |
void | WlzGeomMap3DTriangleTo2D (WlzDVertex3 p0, WlzDVertex3 p1, WlzDVertex3 p2, WlzDVertex2 *dstQ1, WlzDVertex2 *dstQ2) |
Given the three vertices of a triangle in 3D computes the 2D coordinates within the plane of the thriangle. More... | |
int | WlzGeomTriangleAABBIntersect3D (WlzDVertex3 t0, WlzDVertex3 t1, WlzDVertex3 t2, WlzDVertex3 b0, WlzDVertex3 b1, int tst) |
Tests for an intersection between the given triangle and the axis aligned bounding box in 3D using the Separating Axis Theorem (SAT). More... | |
int | WlzGeomTriangleTriangleIntersect2D (WlzDVertex2 s0, WlzDVertex2 s1, WlzDVertex2 s2, WlzDVertex2 t0, WlzDVertex2 t1, WlzDVertex2 t2) |
Tests for an intersection between the two triangles in 2D using the by testing for intersections between the line segments of the triangles and then for either triangle being contained within the other. More... | |
int | WlzGeomTriangleTriangleIntersect3D (WlzDVertex3 s0, WlzDVertex3 s1, WlzDVertex3 s2, WlzDVertex3 t0, WlzDVertex3 t1, WlzDVertex3 t2) |
Tests for an intersection between the two triangles in 3D using the Separating Axis Theorem (SAT). More... | |
WlzErrorNum | WlzGeomCurvature (int nC, double *dstC, WlzDVertex3 nrm, int nV, WlzDVertex3 *vtx) |
Computes the principle curvatures of a parabolic surface fitted to the given vertices at the first of these vertices. More... | |
WlzErrorNum | WlzGeometryLSqOPlane (WlzDVertex3 *dstNrm, WlzDVertex3 *dstCen, int nVtx, WlzDVertex3 *vtx) |
Computes the plane which is the least squares best fit to the given vertices, where this is with respect to the orthogonal distance from the vertices to the plane. More... | |
double | WlzGeomTriangleVtxDistSq3D (WlzDVertex3 *dstPT, int *dstZT, int *dstIT, double *dstL0, double *dstL1, WlzDVertex3 vT, WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2) |
Computes the minimum distance from the test vertex to the triangle. This algorithm is based on "Distance Between Point
and Triangle in 3D", David Eberly, Geometric Tools, 1999. In this algorithm the distance is computed for a parameterised triangle in 3D for which the following regions R[0-6] ara e defined: More... | |
double | WlzGeomTriangleVtxDistSq2D (WlzDVertex2 *dstU, int *dstEI, WlzDVertex2 vT, WlzDVertex2 v0, WlzDVertex2 v1, WlzDVertex2 v2) |
Determines the (squared) distance from the test vertex to the closest point on a triangle edge. The triangle vertices should be ordered as for WlzGeomVxInTriangle2D(). More... | |
double | WlzGeomTetrahedronVtxDistSq3D (WlzDVertex3 *dstU, int *dstFI, WlzDVertex3 vT, WlzDVertex3 v0, WlzDVertex3 v1, WlzDVertex3 v2, WlzDVertex3 v3) |
Determines the (squared) distance from the test vertex to the closest point on a tetrahedron face. The indexing of the faces for the computed closest face of the tetrahedron vertices is as follows: \[ f0 = \{v0, v1, v2\} f1 = \{v0, v3, v1\} f2 = \{v0, v2, v3\} f3 = \{v2, v1, v3\} \] This is the same as the vertex ordering used for WlzCMeshElm3D. The method used is first to test whether the vertex is outside, on or inside the tetrahedron and then compute the closest point on each face and the squared distance from the test vertex to that point. The distances and points of intersection are computed using WlzGeomTriangleVtxDistSq3D(). More... | |
WlzErrorNum | WlzGeomPolyTriangulate2D (int nPVtx, WlzDVertex2 *pVtx, int *dstNTri, int **dstTri) |
Given a sorted list or polygon vertex positions. The polygon is triangulated by ear clipping and the resulting triangulation is returned. More... | |
Geometric utility functions.
This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
typedef struct _WlzGeomPolyListItem2D WlzGeomPolyListItem2D |
double cbrt | ( | double | c | ) |
Referenced by WlzGeomItrSpiral3I(), WlzGeomItrSpiralShell(), and WlzGeoModelGridWSpSet3D().