Woolz Image Processing  Version 1.8.3
WlzGeometry

Files

file  WlzGeometry.c
 Geometric utility functions.
 

Data Structures

struct  _WlzGeomPolyListItem2D
 Item for polygon vertex list. See WlzGeomPolyTriangulate2D(). More...
 

Functions

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...
 

Detailed Description

Function Documentation

◆ WlzGeomTriangleCen2D()

WlzDVertex2 WlzGeomTriangleCen2D ( WlzDVertex2  v0,
WlzDVertex2  v1,
WlzDVertex2  v2 
)

Computes the position of the centroid of the given triangle in 2D.

Returns
Position of the centroid of the given triangle.
Parameters
v0First vertex of the triangle.
v1Second vertex of the triangle.
v2Third vertex of the triangle.

References WLZ_VTX_2_ADD3, and WLZ_VTX_2_SCALE.

◆ WlzGeomTriangleCen3D()

WlzDVertex3 WlzGeomTriangleCen3D ( WlzDVertex3  v0,
WlzDVertex3  v1,
WlzDVertex3  v2 
)

Computes the position of the centroid of the given triangle in 3D.

Returns
Position of the centroid of the given triangle.
Parameters
v0First vertex of the triangle.
v1Second vertex of the triangle.
v2Third vertex of the triangle.

References WLZ_VTX_3_ADD3, and WLZ_VTX_3_SCALE.

◆ WlzGeomTetrahedronCen3D()

WlzDVertex3 WlzGeomTetrahedronCen3D ( WlzDVertex3  v0,
WlzDVertex3  v1,
WlzDVertex3  v2,
WlzDVertex3  v3 
)

Computes the position of the centroid of the given tetrahedron in 3D.

Returns
Position of the centroid of the given tetrahedron.
Parameters
v0First vertex of the tetrahedron.
v1Second vertex of the tetrahedron.
v2Third vertex of the tetrahedron.
v3Fourth vertex of the tetrahedron.

References WLZ_VTX_3_ADD4, and WLZ_VTX_3_SCALE.

◆ WlzGeomTriangleCircumcentre()

int WlzGeomTriangleCircumcentre ( WlzDVertex2 ccVx,
WlzDVertex2  vx0,
WlzDVertex2  vx1,
WlzDVertex2  vx2 
)

Computes the circumcentre of the given triangle.

Returns
Zero if the circumcentre of the triangle lies at infinity, else non-zero.

Given a triangle \((a_0, a_1), (b_0, b_1), (c_0, c_1)\) then the circumcentre \((p_0, p_1)\) is given by:

\[ p0 = (a_0^2 b_1 - a_0^2 c_1 - b_1^2 a_1 + c_1^2 a_1 + b_0^2 c_1 + a_1^2 b_1 + c_0^2 a_1 - c_1^2 b_1 - c_0^2 b_1 - b_0^2 a_1 + b_1^2 c_1 - a_1^2 c_1) / D \]

\[ p1 = (a_0^2 c_0 + a_1^2 c_0 + b_0^2 a_0 - b_0^2 c_0 + b_1^2 a_0 - b_1^2 c_0 - a_0^2 b_0 - a_1^2 b_0 - c_0^2 a_0 + c_0^2 b_0 - c_1^2 a_0 + c_1^2 b_0) / D \]

Where:

\[ D = 2 (a_1 c_0 + b_1 a_0 - b_1 c_0 - a_1 b_0 - c_1 a_0 + c_1 b_0) \]

This is taken from J. O'Rourke: Computational Geometry in C, p201.

Parameters
ccVxDestination pointer for the circumcentre.
vx0First vertex of triangle.
vx1Second vertex of triangle.
vx2Third vertex of triangle.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

Referenced by WlzMeshElemSplit().

◆ WlzGeomVxInTriangle2D()

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.

Returns
Value indicating the position of the vertex with respect to the triangle: +ve if the vertex is inside the triangle, 0 if the vertex is on an edge of the triangle and -ve if the vertex is outside the triangle.

If a triangle has vertices \(p_0, p_1, p_2\), then any point in the plane containing the triangle can be represented by: \(p = \lambda_0 p_0 + \lambda_1 p_2 + \lambda_2 p_3\) subject to the constraint: \(\lambda_0 + \lambda_1 + \lambda_2 = 1\) \(p\) is outside the triangle at one or more of \(\lambda_0\), \(\lambda_1\) and \(\lambda_2\) is -ve. It is inside if all are +ve and on an edge of the triangle if any are close to zero (ie < 10 * ALG_DBL_TOLLERANCE).

Parameters
p0First vertex of triangle.
p1Second vertex of triangle.
p2Third vertex of triangle.
pPGiven vertex.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.

Referenced by WlzCMeshElmEnclosesPos2D(), and WlzGeomTriangleVtxDistSq2D().

◆ WlzGeomVxInTriangle3D()

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()).

Returns
Value indicating the position of the vertex with respect to the triangle: +ve if the vertex is inside the triangle, 0 if the vertex is on an edge of the triangle and -ve if the vertex is outside the triangle.
Parameters
v0First vertex of triangle.
v1Second vertex of triangle.
v2Third vertex of triangle.
vQGiven query vertex.
vPMaxMaximum plane vertex distance.

References ALG_DBL_TOLLERANCE, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SCALE, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshElmEnclosesPos2D5(), and WlzGeomTetrahedronVtxDistSq3D().

◆ WlzGeomVxInTetrahedron()

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.

Returns
Value indicating the position of the vertex with respect to the tetrahedron: +ve if the vertex is inside the tetrahedron, 0 if the vertex is on an edge of the tetrahedron and -ve if the vertex is outside the tetrahedron.

If a tetrahedron has vertices \(p_0, p_1, p_2, p_3\), then any point in the 3D space containing the tetrahedron can be represented by:

\[p = \lambda_0 p_0 + \lambda_1 p_1 + \lambda_2 p_2 + \lambda_3 p_3\]

subject to the constraint:

\[\lambda_0 + \lambda_1 + \lambda_2 + \lambda_3 = 1\]

\(p\) is outside the tetrahedron if one or more of \(\lambda_0\), \(\lambda_1\), \(\lambda_2\) and \(\lambda_3\) is -ve. It is inside if all are +ve and on an edge of the tetrahedron if any are close to zero (ie fabs(x) < 10 * ALG_DBL_TOLLERANCE).

    The barycentric coordinates are computed by inverting the

The tetrahedron vertices

\[ V = \left[ \begin{array}{cccc} vx_0 & vx_1 & vx_2 & vx_3 \\ vy_0 & vy_1 & vy_2 & vy_3 \\ vz_0 & vz_1 & vz_2 & vz_3 \\ 1 & 1 & 1 & 1 \end{array} \right] \]

the barycentric coordinates

\[ L = \left[ \begin{array}{c} \lambda_0 \\ \lambda_1 \\ \lambda_2 \\ \lambda_3 \end{array} \right] \]

and the point to be queried

\[ P = \left[ \begin{array}{c} px \\ py \\ pz \\ 1 \end{array} \right] \]

can be written

\[V L = P \]

and solved for the the barycentric coordinates using

\[L = V^{-1} P\]

.

Parameters
v0First vertex of tetrahedron.
v1Second vertex of tetrahedron.
v2Third vertex of tetrahedron.
v3Fourth vertex of tetrahedron.
vPGiven vertex.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

Referenced by WlzCMeshElmEnclosesPos3D(), and WlzGeomTetrahedronVtxDistSq3D().

◆ WlzGeomTriangleSnArea2I()

WlzLong WlzGeomTriangleSnArea2I ( WlzIVertex2  vx0,
WlzIVertex2  vx1,
WlzIVertex2  vx2 
)

Computes twice the signed area of the given triangle.

Returns
Twice the signed area of the given triangle.
    Computes twice the signed area of the given triangle.
    The determinant is NOT computed with:

\[(x_0 - x_1)(y_1 - y_2) - (y_0 - y_1)(x_1 - x_2)\]

instead the factorized form is used because it is more robust numericaly.
Parameters
vx0First vertex of triangle.
vx1Second vertex of triangle.
vx2Third vertex of triangle.

References _WlzIVertex2::vtX, and _WlzIVertex2::vtY.

◆ WlzGeomTriangleSnArea2()

double WlzGeomTriangleSnArea2 ( WlzDVertex2  vx0,
WlzDVertex2  vx1,
WlzDVertex2  vx2 
)

Computes twice the signed area of the given triangle.

Returns
Twice the signed area of the given triangle.
    Computes twice the signed area of the given triangle.
    The determinant is NOT computed with:

\[(x_0 - x_1)(y_1 - y_2) - (y_0 - y_1)(x_1 - x_2)\]

instead the factorized form is used because it is more robust numericaly.
Parameters
vx0First vertex of triangle.
vx1Second vertex of triangle.
vx2Third vertex of triangle.

References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

Referenced by WlzCMeshBoundConform2D(), WlzCMeshCmpElmFeat2D(), WlzCMeshElmSnArea22D(), WlzCMeshFixNegativeElms2D(), WlzCMeshSetElm2D(), WlzGMEdgeTInsertRadial(), WlzMeshClosestNod2D(), WlzMeshElemFindVxWalk(), WlzMeshElemReplace1(), WlzMeshElemReplaceNWithN(), WlzMeshElemSplit(), WlzMeshElemVerify(), WlzMeshIDomAdd(), WlzMeshTransformAdapt(), and WlzMeshVxVecAdd().

◆ WlzGeomTetraSnVolume6()

double WlzGeomTetraSnVolume6 ( WlzDVertex3  vx0,
WlzDVertex3  vx1,
WlzDVertex3  vx2,
WlzDVertex3  vx3 
)

Computes six times the signed volume of the given tetrahedron.

Returns
Six times the signed volume of the given tetrahedron.
          The signed volume is computed using simple determinant
          evaluation:

\[ area \times 6 = \left| \begin{array}{cccc} x_0 & y_0 & z_0 & 1 \\ x_1 & y_1 & z_1 & 1 \\ x_2 & y_2 & z_2 & 1 \\ x_3 & y_3 & z_3 & 1 \end{array} \right| \]

which can be written

\[ area \times 6 = x_0 ((y_2 z_3 + y_1 (z_2 - z_3)) - (y_3 z_2 + z_1 (y_2 - y_3))) - y_0 ((x_2 z_3 + x_1 (z_2 - z_3)) - (x_3 z_2 + z_1 (x_2 - x_3))) + z_0 ((x_2 y_3 + x_1 (y_2 - y_3)) - (x_3 y_2 + y_1 (x_2 - x_3))) - x_1 (y_2 z_3 - y_3 z_2) + y_1 (x_2 z_3 - x_3 z_2) - z_1 (x_2 y_3 - x_3 y_2) \]

Simple evaluation of this determinant is not robust.
Parameters
vx0First vertex of tetrahedron.
vx1Second vertex of tetrahedron.
vx2Third vertex of tetrahedron.
vx3Forth vertex of tetrahedron.

References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

Referenced by WlzCMeshBoundConform3D(), WlzCMeshCmpElmFeat3D(), WlzCMeshElmSnVolume63D(), WlzCMeshFixNegativeElms3D(), WlzCMeshSetElm3D(), and WlzGeomTetraInSphereDiam().

◆ WlzGeomTetVolZeroI()

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.

Returns
Non zero if the volume of tetrahedron is zero.
Parameters
v0First vertex of tetrahedron.
v1Second vertex of tetrahedron.
v2Third vertex of tetrahedron.
v3Forth vertex of tetrahedron.

References WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SUB, and WlzGeomTriangleArea2Sq3I().

◆ WlzGeomTetVolZeroD()

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.

Returns
Non zero if the volume of tetrahedron is zero.
Parameters
v0First vertex of tetrahedron.
v1Second vertex of tetrahedron.
v2Third vertex of tetrahedron.
v3Forth vertex of tetrahedron.

References WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SUB, and WlzGeomTriangleArea2Sq3().

◆ WlzGeomTriangleArea2Sq3I()

WlzLong WlzGeomTriangleArea2Sq3I ( WlzIVertex3  vx0,
WlzIVertex3  vx1,
WlzIVertex3  vx2 
)

Computes twice the square of the area of the given 3D triangle.

Returns
Twice the square of the area of the given triangle.
    A nieve approach is used in which the area \form#12 is
    computed using:

\[ 2 A^2 = \left|\left| \: \left| \begin{array}{ccc} \mathbf{i} & \mathbf{j} & \mathbf{k} \\ a_x & a_y & a_z \\ b_x & b_y & b_z \end{array} \right| \: \right|\right|^2 \]

Where \(\mathbf{a} = \mathbf{v_0} - \mathbf{v_1}\) and \(\mathbf{b} = \mathbf{v_2} - \mathbf{v_1}\).
Parameters
vx0First vertex of triangle \(\mathbf{v_0}\).
vx1Second vertex of triangle \(\mathbf{v_1}\).
vx2Third vertex of triangle \(\mathbf{v_2}\).

References WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzGeomTetVolZeroI().

◆ WlzGeomTriangleArea2Sq3()

double WlzGeomTriangleArea2Sq3 ( WlzDVertex3  vx0,
WlzDVertex3  vx1,
WlzDVertex3  vx2 
)

Computes twice the square of the area of the given 3D triangle.

Returns
Twice the square of the area of the given triangle.
    A nieve approach is used in which the area \form#12 is
    computed using:

\[ 2 A^2 = \left|\left| \: \left| \begin{array}{ccc} \mathbf{i} & \mathbf{j} & \mathbf{k} \\ a_x & a_y & a_z \\ b_x & b_y & b_z \end{array} \right| \: \right|\right|^2 \]

Where \(\mathbf{a} = \mathbf{v_0} - \mathbf{v_1}\) and \(\mathbf{b} = \mathbf{v_2} - \mathbf{v_1}\).
Parameters
vx0First vertex of triangle \(\mathbf{v_0}\).
vx1Second vertex of triangle \(\mathbf{v_1}\).
vx2Third vertex of triangle \(\mathbf{v_2}\).

References WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshElmSqArea22D5(), WlzCMeshSetElm2D5(), WlzGeomTetraInSphereDiam(), WlzGeomTetVolZeroD(), WlzGMModelConstructSimplex3N(), and WlzGMModelSpxStats().

◆ WlzGeomInTriangleCircumcircle()

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.

Returns
Non zero if the given vertex is inside the circumcircle of the given triangle.
Parameters
vx0First vertex of triangle.
vx1Second vertex of triangle.
vx2Third vertex of triangle.
gVxGiven vertex to test.

References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

◆ WlzGeomLineSegmentsIntersect()

int WlzGeomLineSegmentsIntersect ( WlzDVertex2  p0,
WlzDVertex2  p1,
WlzDVertex2  q0,
WlzDVertex2  q1,
WlzDVertex2 dstN 
)

Tests to see if the two given line segments intersect.

Returns
Integer value which classifies the intersection of the line segments, with values:
  • 0 no intersection.
  • 1 intersection along line segments or all points are coincident.
  • 2 intersection at end points.
  • 3 intersection at a single point (not end points).

Tests to see if the two given line segments intersect using the ALG_DBL_TOLLERANCE tollerance value. This is taken from J. O'Rourke: Computational Geometry in C, p250, but has ben modified to include the use of ALG_DBL_TOLLERANCE.

Parameters
p01st vertex of 1st line segment.
p12nd vertex 1st line segment.
q01st vertex of 2nd line segment.
q12nd vertex of 2nd line segment.
dstNDestination ptr for intersection vertex, may be NULL. The intersection value is not set if there is no intersection or the intersection is along the line segmants.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

◆ WlzGeomCmpAngle()

int WlzGeomCmpAngle ( WlzDVertex2  p0,
WlzDVertex2  p1 
)

Given two end connected 2D line segments this function compares the CCW angle of the segments.

Returns
Result of comparison: -ve, 0 or +ve. Only the sign is meaningful.

Given two end connected 2D line segments: \((p_0, O)\) and \((p_1, O)\), compares the CCW angle of the segments, where \(O\) is the origin \((0, 0)\).

Parameters
p01st segment endpoint vertex.
p12nd segment endpoint vertex.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

◆ WlzGeomVtxEqual2D()

int WlzGeomVtxEqual2D ( WlzDVertex2  pos0,
WlzDVertex2  pos1,
double  tolSq 
)

Checks to see if two verticies are the same within some tollerance.

Returns
1 if node positions are equal, else 0.
Parameters
pos0First node position.
pos1Second node position.
tolSqSquare of tollerance value.

References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

Referenced by WlzGeomVtxOnLineSegment2D(), and WlzGMModelConstructSimplex2N().

◆ WlzGeomVtxEqual3D()

int WlzGeomVtxEqual3D ( WlzDVertex3  pos0,
WlzDVertex3  pos1,
double  tol 
)

Checks to see if two verticies are the same within some tollerance.

Returns
Non-zero if node positions are equal, else 0.
Parameters
pos0First node position.
pos1Second node position.
tolTollerance value.

References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

◆ WlzGeomVtxSortRadial()

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.

Returns
void
    Sorts the given 3D verticies, which lie in a plane
          perpendicular to the radial vector, in order of their
          angle the radial vector.
          No checks are made of the given parameters validity,
          it's assumed that:
            (nV > 0) &&
            (vP != NULL) && (wP != NULL) && (iP != NULL)
            (|rV| > 0) && (rV.(uV = *vP) == 0)
          Note that it is the indicies that are sorted NOT the
          verticies themselves.
Parameters
nVNumber of 3D verticies.
vPThe 3D verticies.
idxBufBuffer of nV indicies used for sorting the verticies.
wPWorkspace with nV 2D verticies.
rVThe radial vector.

References AlgHeapSortIdx(), WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SCALE.

◆ WlzGeomTriangleNormal()

WlzDVertex3 WlzGeomTriangleNormal ( WlzDVertex3  v0,
WlzDVertex3  v1,
WlzDVertex3  v2 
)

Computes the unit normal vector perpendicular to the triangle \(v_0, v_1, v_2\).

Returns
Normal vector.
Parameters
v0First vertex of triangle.
v1Second vertex of triangle.
v2Third vertex of triangle.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_CROSS, WLZ_VTX_3_LENGTH, WLZ_VTX_3_SCALE, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshComputeNormalsElm(), WlzCMeshIntersectDom2D5(), and WlzGMVertexNormal3D().

◆ WlzGeomTriangleAABBIntersect2D()

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).

Returns
The result of the intersection test: 0 - no intersection, 1 - triangle and box are touch or intersect.

Given an axis aligned bounding box and a triangle this function tests for an intersection using the Separating Axis Theorem (SAT) which states : Two 2D convex domains do not intersect iff there exists a line, called a separating axis, on which projection intervals of the domains do not intersect. The minimal set of axes that need to be considered is formed by the normals to all edges of the (polygonal) domains. For an axis aligned bounding box and a triangle in 2D, this is equivalent to testing for the intersection of the given axis aligned bounding box with the axis aligned bounding box of the triangle and the axes normal to the faces of the triangle. The mathematics are simplified by the box being axis aligned.

The algorithm may return false positives when the domains are very close to touching.

Parameters
t0First vertex of triangle.
t1Second vertex of triangle.
t2Third vertex of triangle.
b0Minimum coordinates of axis aligned bounding box.
b1Maximum coordinates of axis aligned bounding box.
tstDetermines the actual intersection tests used: 0 - AABB / triangle. 1 - AABB / AABB(triangle) only. 2 - AABB / triangle omitting the AABB / AABB(triangle) test this is probably only useful if the AABB / AABB(triangle) are known to intersect.

References ALG_DBL_TOLLERANCE, ALG_MAX3, ALG_MIN3, _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.

◆ WlzGeomTetrahedronAABBIntersect3D()

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).

Returns
The result of the intersection test: 0 - no intersection, 1 - tetrahedron and box touch or intersect.

Given an axis aligned bounding box and a tetrahedron this function tests for an intersection using the Separating Axis Theorem (SAT) which states : Two 3D convex domains do not intersect iff there exists a line, called a separating axis, on which projection intervals of the domains do not intersect. The minimal set of axes that need to be considered is formed by the normals to all faces of the (polyhedral) domains and the cross product of all edge combinations in which one edge is from each polyhedron. For an axis aligned bounding box and a tetrahedron in 3D, this is equivalent to testing for the intersection of the given axis aligned bounding box with the axis aligned bounding box of the tetrahedron, testing for intersections on the axes normal to the faces of the tetrahedron and testing for intersection along the cross product of the axis aligned bounding box - tetrahedron edges. The mathematics are simplified by the box being axis aligned.

The algorithm may return false positives when the domains are very close to touching.

Parameters
t0First vertex of tetrahedron.
t1Second vertex of tetrahedron.
t2Third vertex of tetrahedron.
t3Fourth vertex of tetrahedron.
b0Minimum coordinates of axis aligned bounding box.
b1Maximum coordinates of axis aligned bounding box.
tstDetermines the actual intersection tests used: 0 - AABB / tetrahedron. 1 - AABB / AABB(tetrahedron) only. 2 - AABB / tetrahedron omitting the AABB / AABB(tetrahedron) test this is probably only useful if the AABB / AABB(tetrahedron) are known to intersect.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

◆ WlzGeomPlaneAABBIntersect()

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.

Returns
Non-zero if there is an intersection.
Parameters
aPlane X parameter.
bPlane Y parameter.
cPlane Z parameter.
dOther plane parameter.
boxAxis aligned bounding box.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, _WlzDBox3::xMax, _WlzDBox3::xMin, _WlzDBox3::yMax, _WlzDBox3::yMin, _WlzDBox3::zMax, and _WlzDBox3::zMin.

Referenced by Wlz3DViewIntersectAABB().

◆ WlzGeomPlaneLineIntersect()

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\).

Returns
0 if the plane and line do not intersect, 1 if the line segment intersects the plane at a single point or 2 if the line segment is wholly on the plane.
Parameters
aPlane X parameter.
bPlane Y parameter.
cPlane Z parameter.
dOther plane parameter.
p0First end point of the line segment.
p1Second end point of the line segment.
dstIsnDestination pointer for point of intersection, may be NULL.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

Referenced by WlzGeomPlaneTriangleIntersect().

◆ WlzGeomPlaneTriangleIntersect()

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.

Returns
An intersection code:
  • 0 if the plane and triangle do not intersect;
  • 1 if a single vertex of the triangle is on the plane;
  • 2 if a single edge or a line through the triangle is on the plane;
  • 3 if the triangle is wholly on the plane.

Tests for an intersection between the plane defined by the equation: \(ax + by + cz + d = 0\) and the triangle with end verticies \(p_0\), \(p_1\) and \(p_2\). If the destination pointers for the intersection points are not NULL and the intersection code is 1 then the single point of intersection is returned in dstIsn0. If the destination pointers are not NULL and the intersection code is either 1 or 2 then the twther the single intersection point is returned in dstIsn0 or the two intersection points are returned in dstIsn0 and dstIsn1.

Parameters
aPlane X parameter.
bPlane Y parameter.
cPlane Z parameter.
dOther plane parameter.
p0First triangle vertex.
p1Second triangle vertex.
p2Third triangle vertex.
dstIsn0Destination pointer for first point of intersection, may be NULL.
dstIsn1Destination pointer for second point of intersection, may be NULL.

References ALG_DBL_TOLLERANCE, WLZ_VTX_3_SQRLEN, WLZ_VTX_3_SUB, and WlzGeomPlaneLineIntersect().

Referenced by WlzGetSectionFromGMModel().

◆ WlzGeomEllipseVxDistSq()

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.

Returns
Distance ratio.
    Equation of ellipse is:

\[ {(\frac{x}{a})}^2 + {(\frac{y}{b})}^2 = 1 \]

and a straight line through the origin:

\[ y = m x \]

Solving for \(x\) and \(y\) at the ellipse gives the square of the distance ratio for given point \(\mathbf{p}\) is

\[ d = \frac{({(p_x - c_x)}^2 + {(p_y - c_y)}^2) (a_x^2 m^2 + a_y^2)} {a_x^2 a_y^2 ( 1 + m^2)} \]

with \(m = \frac{p_y - c_y}{p_x - c_x}\).
Parameters
centreCentre of ellipse.
sAxEllipse semi axes, both \(> 0.0\).
gPntGiven point.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.

◆ WlzGeomHashVtx3D()

unsigned int WlzGeomHashVtx3D ( WlzDVertex3  pos,
double  tol 
)

Computes a hash value from a given 3D double precision position.

Returns
Hash value.
Parameters
posGiven position.
tolTolerance, \(x = x \pm tol\).

References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

Referenced by WlzGMModelAddVertexToHT(), WlzGMModelMatchVertexG3D(), and WlzGMModelRemVertex().

◆ WlzGeomHashVtx2D()

unsigned int WlzGeomHashVtx2D ( WlzDVertex2  pos,
double  tol 
)

Computes a hash value from a given 2D double precision position.

Returns
Hash value.
Parameters
posGiven position.
tolTolerance, \(x = x \pm tol\).

References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

Referenced by WlzGMModelMatchVertexG2D().

◆ WlzGeomCmpVtx3D()

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.

Returns
The result of the vertex comparison: -1, 0 or +1.
Parameters
pos0First vertex.
pos1Second vertex.
tolTolerance, \(x = x \pm tol\).

References _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.

Referenced by WlzBasisFnIMQ3DFromCPts(), WlzBasisFnMQ3DFromCPts(), and WlzGMVertexCmpSign3D().

◆ WlzGeomCmpVtx2D()

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.

Returns
The sign of the vertex comparison: -1, 0 or +1.
Parameters
pos0First vertex.
pos1Second vertex.
tolTolerance, \(x = x \pm tol\).

References _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.

Referenced by WlzBasisFnGauss2DFromCPts(), WlzBasisFnIMQ2DFromCPts(), WlzBasisFnMQ2DFromCPts(), WlzBasisFnTPS2DFromCPts(), and WlzGMVertexCmpSign2D().

◆ WlzGeomUnitVector2D()

WlzDVertex2 WlzGeomUnitVector2D ( WlzDVertex2  vec)

Computes the 2D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\).

Returns
Unit vector or zero vector if components are both zero.
Parameters
vecGiven vector, \(\mathbf{v}\).

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SCALE.

◆ WlzGeomUnitVector3D()

WlzDVertex3 WlzGeomUnitVector3D ( WlzDVertex3  vec)

Computes the 3D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\).

Returns
Unit vector or zero vector if components are both zero.
Parameters
vecGiven vector, \(\mathbf{v}\).

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SCALE.

◆ WlzGeomUnitVector2D2()

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.

Returns
Unit vector or zero vector if vertices are coincident.
Parameters
v1Position of vertex, \(\mathbf{p}_1\).
v0Position of vertex, \(\mathbf{p}_0\).

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_LENGTH, WLZ_VTX_2_SCALE, and WLZ_VTX_2_SUB.

◆ WlzGeomUnitVector3D2()

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.

Returns
Unit vector or zero vector if vertices are coincident.
Parameters
v1Position of vertex, \(\mathbf{p}_1\).
v0Position of vertex, \(\mathbf{p}_0\).

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_LENGTH, WLZ_VTX_3_SCALE, and WLZ_VTX_3_SUB.

◆ WlzGeomVertexInDiamCircle()

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\).

Returns
Non-zero if point is inside diametral circle.
Parameters
lPos0Vertex at one end of line segment.
lPos1Vertex at other end of line segment.
posPosition of point.

References WLZ_VTX_2_DOT, and WLZ_VTX_2_SUB.

◆ WlzGeomItrSpiralRing()

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.

Returns
Ring of spiral.
Parameters
stepSpiral step count.

◆ WlzGeomItrSpiral2I()

int WlzGeomItrSpiral2I ( int  step,
int *  pX,
int *  pY 
)

Iterates the given positions coordinates through a 2D expanding integer spiral.

Returns
Incremented spiral step count.
Parameters
stepSpiral step count, must be zero when this function is called for the for first step.
pXDestination pointer for column coordinate.
pYDestination pointer for line coordinate.

◆ WlzGeomItrSpiralShell()

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.

Returns
Shell of spiral.
Parameters
stepSpiral step count.

References cbrt().

◆ WlzGeomItrSpiral3I()

int WlzGeomItrSpiral3I ( int  step,
int *  pX,
int *  pY,
int *  pZ 
)

Iterates the given positions coordinates through a 3D expanding integer spiral.

Returns
Incremented spiral step count.
Parameters
stepSpiral step count, must be zero when this function is called for the for first step.
pXDestination pointer for column coordinate.
pYDestination pointer for line coordinate.
pZDestination pointer for plane coordinate.

References cbrt().

◆ WlzGeomDistSq2D()

double WlzGeomDistSq2D ( WlzDVertex2  v0,
WlzDVertex2  v1 
)

Computes square of the Euclidean distance between the given two vertices.

Returns
Euclidean distance between the given vertices.
Parameters
v0First of the given vertices.
v1Second of the given vertices.

References WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.

Referenced by WlzCMeshElmMaxSqEdgLen2D(), and WlzCMeshUpdateMaxSqEdgLen2D().

◆ WlzGeomDistSq3D()

double WlzGeomDistSq3D ( WlzDVertex3  v0,
WlzDVertex3  v1 
)

Computes square of the Euclidean distance between the given two vertices.

Returns
Euclidean distance between the given vertices.
Parameters
v0First of the given vertices.
v1Second of the given vertices.

References WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshElmMaxSqEdgLen2D5(), WlzCMeshElmMaxSqEdgLen3D(), WlzCMeshUpdateMaxSqEdgLen2D5(), and WlzCMeshUpdateMaxSqEdgLen3D().

◆ WlzGeomDist2D()

double WlzGeomDist2D ( WlzDVertex2  v0,
WlzDVertex2  v1 
)

Computes the Euclidean distance between the given two vertices.

Returns
Euclidean distance between the given vertices.
Parameters
v0First of the given vertices.
v1Second of the given vertices.

References WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SUB.

◆ WlzGeomDist3D()

double WlzGeomDist3D ( WlzDVertex3  v0,
WlzDVertex3  v1 
)

Computes the Euclidean distance between the given two vertices.

Returns
Euclidean distance between the given vertices.
Parameters
v0First of the given vertices.
v1Second of the given vertices.

References WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SUB.

◆ WlzGeomTriangleAffineSolve()

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.

Returns
Non zero if the area of the triangle is very small.
Parameters
xTrTransform coordinates for x.
yTrTransform coordinates for y.
ddTwice the area of the source triangle.
sVxSource triangle vertices.
dVxDestination triangle vertices.
threshThreshold value for twice the area.

References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

◆ WlzGeomTetraAffineSolveLU()

int WlzGeomTetraAffineSolveLU ( double *  tr,
WlzDVertex3 sVx,
WlzDVertex3 dVx 
)

Computes the affine transform coefficients from the source to target tetrahedron.

Returns
Non zero if the volume of the tetrahedron is very small.
    This is done by solving for general 3D affine transform
    matrix which transforms the source tetrahedrons vertices
    to those of the destination tetrahedron
    If the destination tetrahedron vertices are given by

\[ D = \left( \begin{array}{cccc} d_{x0} & d_{x0} & d_{x1} & d_{x2} \\ d_{y0} & d_{y0} & d_{y1} & d_{y2} \\ d_{z0} & d_{z0} & d_{z1} & d_{z2} \\ 1 & 1 & 1 & 1 \end{array} \right) \]

and the source vertices by

\[ S = \left( \begin{array}{cccc} s_{x0} & s_{x0} & s_{x1} & s_{x2} \\ s_{y0} & s_{y0} & s_{y1} & s_{y2} \\ s_{z0} & s_{z0} & s_{z1} & s_{z2} \\ 1 & 1 & 1 & 1 \end{array} \right) \]

then the transform being sought satisfies

\[ D = T S \]

Solving for \(T\)

\[ T = D S^{-1} \]

LU decomposition is used to solve for \(T\). For the degenerate cases the transformation is given as a simple translation. This function is slower than WlzGeomTetraAffineSolve() by about a factor of 8, but it may be more robust.
Parameters
trTransform matrix with 4x4 contiguous coefficients which are equivalent to the base storage of the matrix in a WlzAffineTransform.
sVxSource tetrahedron vertices.
dVxDestination tetrahedron vertices.

References AlgMatrixLUInvertRaw4(), _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

◆ WlzGeomTetraAffineSolve()

int WlzGeomTetraAffineSolve ( double *  tr,
WlzDVertex3 sVx,
WlzDVertex3 dVx,
double  thresh 
)

Computes the affine transform coefficients from the source to target tetrahedron.

Returns
Non zero if the volume of the tetrahedron is very small.
    This is done by solving for general 3D affine transform
    matrix which transforms the source tetrahedrons vertices
    to those of the destination tetrahedron
    If the destination tetrahedron vertices are given by

\[ D = \left( \begin{array}{cccc} d_{x0} & d_{x0} & d_{x1} & d_{x2} \\ d_{y0} & d_{y0} & d_{y1} & d_{y2} \\ d_{z0} & d_{z0} & d_{z1} & d_{z2} \\ 1 & 1 & 1 & 1 \end{array} \right) \]

and the source vertices by

\[ S = \left( \begin{array}{cccc} s_{x0} & s_{x0} & s_{x1} & s_{x2} \\ s_{y0} & s_{y0} & s_{y1} & s_{y2} \\ s_{z0} & s_{z0} & s_{z1} & s_{z2} \\ 1 & 1 & 1 & 1 \end{array} \right) \]

then the transform being sought satisfies

\[ D = T S \]

Solving for \(T\)

\[ T = D S^{-1} \]

Setting

\[ U = S^{-1} \]

with

\[ U = \left( \begin{array}{cccc} u_{11} & u_{12} & u_{13} & u_{14} \\ u_{21} & u_{22} & u_{23} & u_{24} \\ u_{31} & u_{32} & u_{33} & u_{34} \\ u_{41} & u_{42} & u_{43} & u_{44} \end{array} \right) \]

and

\[ d = \det(S) \]

For efficiency the followingcollect common sub-expressions are used

\begin{eqnarray*} cz2z3 = s_{z2} - s_{z3} \\ cz1z3 = s_{z1} - s_{z3} \\ cz1z2 = s_{z1} - s_{z2} \\ cy2y3 = s_{y2} - s_{y3} \\ cy1y3 = s_{y1} - s_{y3} \\ cy1y2 = s_{y1} - s_{y2} \\ cy2z3y3z2 = s_{y2} s_{z3} - s_{y3} s_{z2} \\ cy1z3y3z1 = s_{y1} s_{z3} - s_{y3} s_{z1} \\ cy1z2y2z1 = s_{y1} s_{z2} - s_{y2} s_{z1} \\ cz0z3 = s_{z0} - s_{z3} \\ cz0z2 = s_{z0} - s_{z2} \\ cy0y3 = s_{y0} - s_{y3} \\ cy0y2 = s_{y0} - s_{y2} \\ cy0z3y3z0 = s_{y0} s_{z3} - s_{y3} s_{z0} \\ cy0z2y2z0 = s_{y0} s_{z2} - s_{y2} s_{z0} \\ cz0z1 = s_{z0} - s_{z1} \\ cy0y1 = s_{y0} - s_{y1} \\ cy0z1y1z0 = s_{y0} s_{z1} - s_{y1} s_{z0} \end{eqnarray*}

giving

\[ d = s_{x0}( s_{y1} cz2z3 - s_{y2} cz1z3 + s_{y3} (cz1z2)) + s_{x1}(-s_{y0} cz2z3 + s_{y2} cz0z3 - s_{y3} (cz0z2)) + s_{x2}( s_{y0} cz1z3 - s_{y1} cz0z3 + s_{y3} (cz0z1)) + s_{x3}(-s_{y0} cz1z2 + s_{y1} cz0z2 - s_{y2} (cz0z1)) \]

and for the non-degenerate case, when \(d \not= 0\)

\begin{eqnarray*} u11 = \frac{1}{d}( s_{y1} cz2z3 - s_{y2} cz1z3 + s_{y3} cz1z2) \\ u12 = \frac{1}{d}(-s_{x1} cz2z3 + s_{x2} cz1z3 - s_{x3} cz1z2) \\ u13 = \frac{1}{d}( s_{x1} cy2y3 - s_{x2} cy1y3 + s_{x3} cy1y2) \\ u14 = \frac{1}{d}(-s_{x1} cy2z3y3z2 + s_{x2} cy1z3y3z1 - s_{x3} cy1z2y2z1) \\ u21 = \frac{1}{d}(-s_{y0} cz2z3 + s_{y2} cz0z3 - s_{y3} cz0z2) \\ u22 = \frac{1}{d}( s_{x0} cz2z3 - s_{x2} cz0z3 + s_{x3} cz0z2) \\ u23 = \frac{1}{d}(-s_{x0} cy2y3 + s_{x2} cy0y3 - s_{x3} cy0y2) \\ u24 = \frac{1}{d}( s_{x0} cy2z3y3z2 - s_{x2} cy0z3y3z0 + s_{x3} cy0z2y2z0) \\ u31 = \frac{1}{d}( s_{y0} cz1z3 - s_{y1} cz0z3 + s_{y3} cz0z1) \\ u32 = \frac{1}{d}(-s_{x0} cz1z3 + s_{x1} cz0z3 - s_{x3} cz0z1) \\ u33 = \frac{1}{d}( s_{x0} cy1y3 - s_{x1} cy0y3 + s_{x3} cy0y1) \\ u34 = \frac{1}{d}(-s_{x0} cy1z3y3z1 + s_{x1} cy0z3y3z0 - s_{x3} cy0z1y1z0) \\ u41 = \frac{1}{d}(-s_{y0} cz1z2 + s_{y1} cz0z2 - s_{y2} cz0z1) \\ u42 = \frac{1}{d}( s_{x0} cz1z2 - s_{x1} cz0z2 + s_{x2} cz0z1) \\ u43 = \frac{1}{d}(-s_{x0} cy1y2 + s_{x1} cy0y2 - s_{x2} cy0y1) \\ u44 = \frac{1}{d}( s_{x0} cy1z2y2z1 - s_{x1} cy0z2y2z0 + s_{x2} cy0z1y1z0) \end{eqnarray*}

and so ( \(T = D U\))

\begin{eqnarray*} t11 = d_{x3} u41 + d_{x2} u31 + d_{x1} u21 + d_{x0} u11 \\ t12 = d_{x3} u42 + d_{x2} u32 + d_{x1} u22 + d_{x0} u12 \\ t13 = d_{x3} u43 + d_{x2} u33 + d_{x1} u23 + d_{x0} u13 \\ t14 = d_{x3} u44 + d_{x2} u34 + d_{x1} u24 + d_{x0} u14 \\ t21 = d_{y3} u41 + d_{y2} u31 + d_{y1} u21 + d_{y0} u11 \\ t22 = d_{y3} u42 + d_{y2} u32 + d_{y1} u22 + d_{y0} u12 \\ t23 = d_{y3} u43 + d_{y2} u33 + d_{y1} u23 + d_{y0} u13 \\ t24 = d_{y3} u44 + d_{y2} u34 + d_{y1} u24 + d_{y0} u14 \\ t31 = d_{z3} u41 + d_{z2} u31 + d_{z1} u21 + d_{z0} u11 \\ t32 = d_{z3} u42 + d_{z2} u32 + d_{z1} u22 + d_{z0} u12 \\ t33 = d_{z3} u43 + d_{z2} u33 + d_{z1} u23 + d_{z0} u13 \\ t34 = d_{z3} u44 + d_{z2} u34 + d_{z1} u24 + d_{z0} u14 \\ t41 = u41 + u31 + u21 + u11 (0) \\ t42 = u42 + u32 + u22 + u12 (0) \\ t43 = u43 + u33 + u23 + u13 (0) \\ t44 = u44 + u34 + u24 + u14 (1) \end{eqnarray*}

For the degenerate cases in which \(d \approx 0\) then the transformation is given as a simple translation.
Parameters
trTransform matrix with 4x4 contiguous coefficients which are equivalent to the base storage of the matrix in a WlzAffineTransform.
sVxSource tetrahedron vertices.
dVxDestination tetrahedron vertices.
threshThreshold value which is the lower limit of 6 x the source tetrahedron's volume.

References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.

◆ WlzGeomObjLineSegIntersect2D()

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.

Returns
Position of intersection.
          Given the line segment \form#247, \form#248 any position
          along the segment can be given by a parameter \form#176
          (range [0-1]), where \form#284.
Parameters
objGiven object. Object must have a valid domain.
p0First vertex.
p1Second vertex.
tolAcceptable placement error.
insideNon-zero if the returned position should be inside or on the boundary, if zero it will be outside or on the boundary.
methodMethod for finding intersection: 0 - bisection, 1 - increment. If increment is used each point along the line segment will be tested until termination, this may be very slow if tol is small, with possibly 1/tol incremental steps.
dstStatDestination pointer for status, may be NULL.

References ALG_MAX, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_ADD, WLZ_VTX_2_LENGTH, WLZ_VTX_2_SCALE, WLZ_VTX_2_SQRLEN, WLZ_VTX_2_SUB, and WlzInsideDomain().

Referenced by WlzCMeshBoundConform2D().

◆ WlzGeomObjLineSegIntersect3D()

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.

Returns
Position of intersection.
          Given the line segment \form#247, \form#248 any position
          along the segment can be given by a parameter \form#176
          (range [0-1]), where \form#284.
Parameters
objGiven object. Object must have a valid domain.
p0First vertex.
p1Second vertex.
tolAcceptable placement error.
insideNon-zero if the returned position should be inside or on the boundary, if zero it will be outside or on the boundary.
methodMethod for finding intersection: 0 - bisection, 1 - increment. If increment is used each point along the line segment will be tested until termination, this may be very slow if tol is small, with possibly 1/tol incremental steps.
dstStatDestination pointer for status, may be NULL.

References ALG_MAX3, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_ADD, WLZ_VTX_3_LENGTH, WLZ_VTX_3_SCALE, WLZ_VTX_3_SQRLEN, WLZ_VTX_3_SUB, and WlzInsideDomain().

Referenced by WlzCMeshBoundConform3D().

◆ WlzGeomTetraInSphereDiam()

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.

Returns
Maximum diameter sphere inscribed within the tetrahedron.
            Diameter of the (maximum) inscribed sphere \form#4 is given
            by:

\[ d = \frac{6V}{A} \]

where \(V\) is the volume of the tetrahedron and \(A\) is it's surface area. "Encyclopedia of Mathematics" ISBN 1402006098 (http://eom.springer.de).

Parameters
vx0First vertex of tetrahedron.
vx1Second vertex of tetrahedron.
vx2Third vertex of tetrahedron.
vx3Forth vertex of tetrahedron.

References ALG_DBL_TOLLERANCE, WlzGeomTetraSnVolume6(), and WlzGeomTriangleArea2Sq3().

◆ WlzGeomTetraInSphereRegDiam()

double WlzGeomTetraInSphereRegDiam ( double  side)

Given the side length of a regular tetrahedron this function computes the maximum diameter of an inscribed sphere.

Returns
Maximum diameter sphere inscribed within the regular tetrahedron.
            Diameter of the (maximum) inscribed sphere \form#4 is given
            by:

\[ d = \frac{6V}{A} \]

where \(V\) is the volume of the tetrahedron and \(A\) is it's surface area. "Encyclopedia of Mathematics" ISBN 1402006098 (http://eom.springer.de). Standard formulae are used for the area and volume.

Parameters
sideLength of an edge.

References ALG_DBL_TOLLERANCE, ALG_M_SQRT2, and ALG_M_SQRT3.

◆ WlzGeomPolar2D()

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:

Returns
Angle, range \([0 - 2\pi]\).
                    (y)
                     ^  pi/2
                     |
                     |
               pi    |
              <------+-------> (x)
                     |       0
                     |
                3pi/2|
                     V
Parameters
orgPosition of the origin.
dstPosition of the destination.
dstRadDestination pointer for the length of the ray, may be NULL.

References ALG_M_PI, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SUB.

◆ WlzGeomCos3V()

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.

Returns
Cosine of angle between line segments (v0, v1) and (v1, v2) with value in the range [0-1].
Parameters
v0First vertex.
v1Second vertex (the common one).
v2Third vertex.

References ALG_DBL_TOLLERANCE, WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.

◆ WlzGeomVtxOnLineSegment2D()

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.

Returns
Non-zero value only if test vertex is on the line segment.
Parameters
tstTest vertex.
seg0First vertex of line segment.
seg1Second vertex of line segment.
tolTollerance.

References _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_SUB, WlzGeomVtxEqual2D(), _WlzDBox2::xMax, _WlzDBox2::xMin, _WlzDBox2::yMax, and _WlzDBox2::yMin.

Referenced by WlzGeomTriangleVtxDistSq2D().

◆ WlzGeomVtxOnLineSegment3D()

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.

Returns
Integer code corresponding to position of the test vertex with respect to the line segment:
  • 0 no intersection.
  • 1 intersection at an end point.
  • 2 intersection at a single point (not end points).
    Consider a line segment from a vertex at \form#287
    to another vertex at \form#288 , with a third test
    vertex at \form#289 , the shortest path from
\(\mathbf{p_x}\) to the line segment will be perpendicular to the line segment. Let the position of the intersection of this perpendicular with the line segment be at \(\mathbf{p}\) , with distance \(d\) from \(\mathbf{p_x}\) , then:

\[ \mathbf{p} = \mathbf{p_0} + s (\mathbf{p_1} - \mathbf{p_0}, \]

\[ d^2 = \| (\mathbf{p_0} - \mathbf{p_x}) + (\mathbf{p} - \mathbf{p_0}) \|^2, \]

\[ s = \frac{ (\mathbf{p_0} - \mathbf{p_x}) \cdot (\mathbf{p_1} - \mathbf{p_0})} {\| \mathbf{p_1} - \mathbf{p_0} \|^2}, \]

and after substitution:

\[ d^2 = \frac{\| \mathbf{p_0} - \mathbf{p_x} \|^2 \| \mathbf{p_1} - \mathbf{p_0} \|^2 - \left[(\mathbf{p_1} - \mathbf{p_0})\right]^2 } {\| \mathbf{p_1} - \mathbf{p_0} \|^2}, \]

Observing that the numerator is a vector quad product

\[ \mathbf{A} \times \mathbf{B} = \|\mathbf{A}\|^2\|\mathbf{A}\|^2 - (\mathbf{A} \cdot \mathbf{Ba})^2 \]

gives

\[ d^2 = \frac{\| \mathbf{p_1} - \mathbf{p_0} \times \mathbf{p_0} - \mathbf{p_x} \|^2 } {\| \mathbf{p_1} - \mathbf{p_0} \|^2}, \]

Obviously if \(\mathbf{p_x}\) is on the line segment then \(d^2 = 0\). ALG_DBL_TOLLERANCE as a tolerance for squared distances in this function.
Parameters
pXTest vertex.
p0First vertex of line segment.
p1Second vertex of line segment.
dstNDestination pointer for the intersection, may be NULL.

References ALG_DBL_TOLLERANCE, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SCALE_ADD, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzGeomLineLineSegmentIntersect3D().

◆ WlzGeomArcLength2D()

double WlzGeomArcLength2D ( WlzDVertex2  a,
WlzDVertex2  b,
WlzDVertex2  c 
)

Computes the arc length from a to b traveling CCW on a circle with centre c.

Returns
The arc length.
Parameters
aStart point.
bEnd point.
cCirecle centre.

References ALG_DBL_TOLLERANCE, ALG_M_PI, ALG_M_PI_2, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.

◆ WlzGeomRectFromWideLine()

int WlzGeomRectFromWideLine ( WlzDVertex2  s,
WlzDVertex2  t,
double  w,
WlzDVertex2 v 
)

Computes the coordinates of vertices that may be used to draw a wide rectangle.

Returns
Non-zero if \(S\) and \(T\) are coincident.
            Given two vertices \form#297, \form#26 which define a line
            segment and width \form#298 perpendicular to the line segent,
            this function computes the coordinates of the vertices
\(V_i\), \(i\in[0-3]\) of the rectangle. The vertices are sorted such that the rectangle may be drawn using the line segmants: \((V_0,V_1)\), \((V_1,V_2)\), \((V_2,V_3)\), \((V_3,V_0)\). Given the line segment \((S,T)\) it can be shown that the vertices which share a line segment passing through \(S\) are:

\[ V_i = S + (-\frac{r (T_x - S_x)}{l}, \frac{r (T_y - S_y)}{l}) \]

and

\[ V_j = S + ( \frac{r (T_x - S_x)}{l}, -\frac{r (T_x - S_y)}{l}) \]

where \(l = ||T - S||\) and \(r = w / 2\). Should \(S\) and \(T\) be coincident then the destination rectangle vertices are left unmodified..
Parameters
sFirst vertex of line segment.
tSecond vertex of line segment.
wline width.
vDestination pointer for the four vertices of the rectangle.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_ADD, WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SUB.

◆ WlzGeomLinePlaneIntersection()

WlzDVertex3 WlzGeomLinePlaneIntersection ( WlzDVertex3  v,
WlzDVertex3  p0,
WlzDVertex3  p1,
WlzDVertex3  p2,
WlzDVertex3  p3,
int *  dstPar 
)

Computes the intersection of a line with a plane.

Returns
Intersection of line with plane.
            Computes the intersection of line (vector \form#261)
            through an off plane vertex ( \form#310) with a plane.
            The plane is defined by three on plane vertices
\((\mathbf{p_0}, \mathbf{p_1}, \mathbf{p_2})\). \(q\) the vertex at the intersection is found by solving the vector equations

\[ (\mathbf{p_1} - \mathbf{p_0}) \times (\mathbf{p_2} - \mathbf{p_0}) . (\mathbf{q} - \mathbf{p_0}) = \mathbf{0}, \mathbf{n} \times (\mathbf{p_3} - \mathbf{q}) = \mathbf{0} \]

These are equivalunt to solving

\[ \left|\begin{array}{ccc} p_{1x} - p_{0x} & p_{1y} - p_{0y} & p_{1z} - p_{0z} \\ p_{2x} - p_{0x} & p_{2y} - p_{0y} & p_{2z} - p_{0z} \\ q_x - p_{0x} & q_y - p_{0y} & q_z - p_{0z} \end{array}\right| = 0, \left|\begin{array}{ccc} \mathbf{i} & \mathbf{j} & \mathbf{k} \\ n_x & n_y & n_z \\ p_{3x} - q_x & p_{3y} - q_y & p_{3z} - q_z \end{array}\right| = 0 \]

Parameters
vUnit vector for ray.
p0First point on the plane.
p1Second point on the plane.
p2Third point on the plane.
p3Off plane point that ray passes through.
dstParDestination value set to 1 if the vector is parrallel to the plane, must not be NULL.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.

◆ WlzGeomLineTriangleIntersect3D()

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.

Returns
Value indicating the position of the vertex with respect to the triangle in 3D:
  • 0 if the vertex is outside the triangle.
  • 1 if the vertex is on an edge of the triangle or line passes through triangle in it's plane.
  • 2 if the vertex is inside the triangle.
    Given a parameterised line

\[ R(t) = O + t D \]

and a triangle specified by it's vertices \((V_0, V_1, V_2)\), the point of intersection may be written in terms of the barycentric coordinates \(u, v\)

\[ O + t D = T(u, v) = (1 - u - v) V_0 + u V_1 + v V_2 \]

Solving for \(t\), \(u\) and \(v\) gives

\[ \left[ \begin{array}{c} t \\ u \\ v \end{array} \right] = \frac{1}{P \cdot E_1} \left[ \begin{array}{c} Q \cdot E_2 \\ P \cdot t \\ Q \cdot D \end{array} \right] \]

where \(E_1 = V_1 - V_0\), \(E_2 = V_2 - V_0\) \(t = O - V_0\), \(P = D \times E_2\) and \(Q = T \times E_1\).
   The \form#325 and \form#6 are only set if the line passes
   through the triangle.
Parameters
orgLine origin, \(O\).
dirLine direction, \(D\) (does not need to be a unit vector).
v0First vertex on triangle.
v1Second vertex on the triangle
v2Third vertex on the triangle
dstParDestination pointer for flag set to 1 if the vector is parrallel to the plane, may be NULL.
dstTDestination pointer for t parameter, may be NULL.
dstUDestination pointer for u parameter, may be NULL.
dstVDestination pointer for v parameter, may be NULL.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, and WLZ_VTX_3_SUB.

◆ WlzGeomLineLineSegmentIntersect3D()

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.

Returns
Integer value which classifies the intersection of the line with the line line segment. Values are:
  • 0 no intersection.
  • 1 intersection along line segment, all points on the line segment are coincident.
  • 2 intersection at end points.
  • 3 intersection at a single point (not end points).
    Given a line parameterised by \form#16:

\[ \mathbf{x} = \mathbf{R_0} + t \mathbf{R_d} \]

and a line segment parameterised by \(s\):

\[ \mathbf{x} = \mathbf{P_0} + s (\mathbf{P_1} - \mathbf{P_0}) \]

their intersection at \(\mathbf{x}\) is given by

\[ \mathbf{P_0} + s (\mathbf{P_1} - \mathbf{P_0}) = \mathbf{R_0} + t \mathbf{R_d} \]

giving

\[ s (\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d} = (\mathbf{R_0} - \mathbf{P_0} \times \mathbf{R_d} t \mathbf{R_d} \times (\mathbf{P_1} - \mathbf{P_0}) = (\mathbf{P_0} - \mathbf{R_0} \times (\mathbf{P_1} - \mathbf{P_0}) \]

\[ s ((\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}) \cdot ((\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}) = ((\mathbf{R_0} - \mathbf{P_0}) \times \mathbf{R_d}) \cdot ((\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}) t (\mathbf{R_d} \times (\mathbf{P_1} - \mathbf{P_0})) \cdot (\mathbf{R_d} \times (\mathbf{P_1} - \mathbf{P_0})) = ((\mathbf{P_0} - \mathbf{R_0} \times (\mathbf{P_1} - \mathbf{P_0})) \cdot (\mathbf{R_d} \times (\mathbf{P_1} - \mathbf{P_0})) \]

\[ s = \frac{\mbox{det}(\mathbf{R_0} - \mathbf{P_0}, \mathbf{R_d}, (\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}} {\|(\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}\|^2} t = \frac{\mbox{det}(\mathbf{R_0} - \mathbf{P_0}, \mathbf{R_d}, (\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}} {\|(\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}\|^2} \]

If the denominator \(\|(\mathbf{P_1} - \mathbf{P_0}) \times \mathbf{R_d}\|^2 = 0\) then the line and the line segment are parrallel, provided that \(\mathbf{P_1} \neq \mathbf{P_0}\) and \(\mathbf{R_d} \neq \mathbf{0}\). The line segment is intersected by the line if \(s\) is in the range [0-1]. The point of intersection \(\mathbf{x}\) is

\[ \mathbf{x} = \mathbf{P_0} + s (\mathbf{P_1} - \mathbf{P_0}) \]

Special cases are considered for \(\mathbf{R_0} = \mathbf{P_0},\mathbf{P_1}\), \(\mathbf{R_0} = 2 \mathbf{P_0} - \mathbf{P_1}\) and \(\mathbf{R_d} = \alpha (\mathbf{P_1} - \mathbf{P_0})\).
Parameters
r0A vertex on the line.
rDDirection of the line.
p0First end vertex of the line segment.
p1Second end vertex of the line segment.
dstNDestination ptr for intersection vertex, may be NULL. The intersection value will be set if the ray and line segment intersect at a single point.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_CROSS, WLZ_VTX_3_SCALE_ADD, WLZ_VTX_3_SQRLEN, WLZ_VTX_3_SUB, WlzGeomVtxOnLine3D(), and WlzGeomVtxOnLineSegment3D().

◆ WlzGeomVtxOnLine3D()

int WlzGeomVtxOnLine3D ( WlzDVertex3  p0,
WlzDVertex3  r0,
WlzDVertex3  rD 
)

Tests whether a vertex is a line.

Returns
If the vertex is on the ray 1, otherwise 0.
            Tests whether the given vertex \form#247 is on the
            given line
\(\mathbf{x} = \mathbf{r_0} + \mathbf{r_D}\). If the point on the line and the given vertex are coincident then obviously the vertex is on the line otherwise the vertex is on the line if the squared length of the cross product of the line direction and the direction of the vertex (with respect to the point on the line) is less than ALG_DBL_TOLLERANCE, ie if

\[ \| (\mathbf{p_0} - \mathbf{r_0}) \times \mathbf{r_D} \| < \epsilon \]

Parameters
p0Given vertex.
r0Point on line.
rDDirection of line.

References ALG_DBL_TOLLERANCE, WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzGeomLineLineSegmentIntersect3D().

◆ WlzGeomBaryCoordsTri2D()

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.

Returns
Non zero if barycentric coordinates valid.
Parameters
p0First vertex of triangle.
p1Second vertex of triangle.
p2Third vertex of triangle.
pXGiven position, which is within (or on) the triangle.
lambdaGiven array for the three lambda values.

References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.

Referenced by WlzGeomInterpolateTri2D().

◆ WlzGeomInterpolateTri2D()

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.

Returns
Interpolated value.
Parameters
p0First vertex of triangle.
p1Second vertex of triangle.
p2Third vertex of triangle.
v0Value at first vertex of triangle.
v1Value at second vertex of triangle.
v2Value at third vertex of triangle.
pXGiven position, which is within (or on) the triangle.

References WlzGeomBaryCoordsTri2D().

Referenced by WlzCMeshCurvToImage(), and WlzCMeshMeshMeshProduct().

◆ WlzGeomInterpolatePoly2D()

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.

Returns
Interpolated value.
Parameters
nNumber of polygon vertices, which is the same as the number of values and generalised barycentric coordinates.
pThe vertices of the convex polygon (ordered counter-clockwise).
vValues at the corresponding polygon vertices.
wUsed to compute and return the generalised barycentric coordinates of the given position. The coordinate values of vertices outside the polygon's convex hull will be zero.
qGiven position, which is must be within or on the convex hull of the polygon.
dstErrDestination error pointer, may be NULL.

References AlcFree(), AlcMalloc(), _WlzGeomPolyListItem2D::v, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WlzConvHullClarkson2D(), and WlzGeomInterpolateConvexPoly2D().

◆ WlzGeomInterpolateConvexPoly2D()

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.

Returns
Interpolated value.
Parameters
nNumber of polygon vertices, which is the same as the number of values and generalised barycentric coordinates.
pThe vertices of the convex polygon (ordered counter-clockwise).
vValues at the corresponding polygon vertices.
wUsed to compute and return the generalised barycentric coordinates of the given position.
qGiven position, which is must be in or on the polygon.

References _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.

Referenced by WlzGeomInterpolatePoly2D().

◆ WlzGeomBaryCoordsTet3D()

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.

Returns
Non zero if barycentric coordinates valid.
Parameters
p0First vertex of tetrahedron.
p1Second vertex of tetrahedron.
p2Third vertex of tetrahedron.
p3Fourth vertex of tetrahedron.
pXGiven position, which is within (or on) the tetrahedron.
lambdaGiven array for the four lambda values.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.

Referenced by WlzGeomInterpolateTet3D().

◆ WlzGeomInterpolateTet3D()

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.

Returns
Interpolated value.
Parameters
p0First vertex of tetrahedron.
p1Second vertex of tetrahedron.
p2Third vertex of tetrahedron.
p3Fourth vertex of tetrahedron.
v0Value at first vertex of tetrahedron.
v1Value at second vertex of tetrahedron.
v2Value at third vertex of tetrahedron.
v3Value at fourth vertex of tetrahedron.
pXGiven position, which is within (or on) the tetrahedron.

References WlzGeomBaryCoordsTet3D().

◆ WlzGeomMap3DTriangleTo2D()

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.

If the 3D coordinates of the vertices of the triangle are \(\mathbf{p_0}\), \(\mathbf{p_1}\) and \(\mathbf{p_2}\) then this function computes the coordinates of vertices in the plane of the triangle, such that: \(\mathbf{q_0} = (0,0)\), \(\mathbf{q_1} = (\|\mathbf{l_1}\|,0)\) and \(\mathbf{q_2} = (\mathbf{l_2}\cdot\mathbf{u}\). Where: \(\mathbf{l_1} = \mathbf{p_1} - \mathbf{p_0}\), \(\mathbf{l_2} = \mathbf{p_2} - \mathbf{p_0}\), \(\mathbf{u} = \frac{\mathbf{1}}{\|\mathbf{l_1}\|} \mathbf{l_1}\), \(\mathbf{n} = \frac{\mathbf{1}} {\|\mathbf{l_1} \times \mathbf{l_2}\|} \mathbf{l_1} \times \mathbf{l_2}\) and \(\mathbf{v} = \mathbf{n} \times \mathbf{l_1}\).

   The first vertex in the plane is not returned because it's
   coordinates are always (0,0).
Parameters
p0First vertex of triangle (origin of the 2D plane).
p1Second vertex of triangle (on the \(\mathbf{u}\) axis in the 2D plane).
p2Third vertex of triangle.
dstQ1Destination pointer for the second vertex in the plane, must not be NULL.
dstQ2Destination pointer for the third vertex in the plane, must not be NULL.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_LENGTH, WLZ_VTX_3_SCALE, and WLZ_VTX_3_SUB.

◆ WlzGeomTriangleAABBIntersect3D()

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).

Returns
The result of the intersection test: 0 - no intersection, 1 - triangle and box touch or intersect.

Given an axis aligned bounding box and a triangle this function tests for an intersection using the Separating Axis Theorem (SAT) which states : Two 3D convex domains do not intersect iff there exists a line, called a separating axis, on which projection intervals of the domains do not intersect. The minimal set of axes that need to be considered is formed by the normals to all faces of the (polyhedral) domains and the cross product of all edge combinations in which one edge is from each polyhedron. For an axis aligned bounding box and a triangle in 3D, this is equivalent to testing for the intersection of the given axis aligned bounding box with the axis aligned bounding box of the triangle, testing for intersections on the axes normal to the face of the triangle and testing for intersection along the cross product of the axis aligned bounding box - triangle edges. The mathematics are simplified by the box being axis aligned.

The algorithm may return false positives when the domains are very close to touching.

Parameters
t0First vertex of the triangle.
t1Second vertex of the triangle.
t2Third vertex of the triangle.
b0Minimum of the bounding box.
b1Maximum of the bounding box.
tstDetermines the actual intersection tests used: 0 - AABB / triangle. 1 - AABB / AABB(triangle) only. 2 - AABB / triangle omitting the AABB / AABB(tetrahedron) test this is probably only useful if the AABB / AABB(triangle) are known to intersect.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshIntersectDom2D5(), and WlzGeoModelGridWSpSet3D().

◆ WlzGeomTriangleTriangleIntersect2D()

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.

Returns
Non zero if the two triangles intersect or touch.
            See WlzGeomTriangleTriangleIntersect2DA().
Parameters
s0First vertex of the first triangle.
s1Second vertex of the first triangle.
s2Third vertex of the first triangle.
t0First vertex of the second triangle.
t1Second vertex of the second triangle.
t2Third vertex of the second triangle.

◆ WlzGeomTriangleTriangleIntersect3D()

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).

Returns
Non zero if the two triangles intersect or touch.
            See WlzGeomTriangleTriangleIntersect3DA().
Parameters
s0First vertex of the first triangle.
s1Second vertex of the first triangle.
s2Third vertex of the first triangle.
t0First vertex of the second triangle.
t1Second vertex of the second triangle.
t2Third vertex of the second triangle.

◆ WlzGeomCurvature()

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.

Returns
Woolz error code.
Parameters
nCNumber of curvatures required: 1 computes the Gaussian curvature, 2 computes both principle curvatures.
dstCDestination pointer for the curvature value(s), must not be NULL. If nC == 2 then the values atr Gaussian followed by mean curvature.
nrmNormal at the first vertex.
nVNumber of vertices, must be >= 3.
vtxArray of vertex positions, the first of which must be the vertex at which the curvature is to be computed. The array contents are modified by this function. Must not be NULL.

References AlcFree(), AlcMalloc(), ALG_DBL_TOLLERANCE, AlgMatrixFree(), AlgMatrixRectNew(), AlgMatrixSVSolve(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrix::rect, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_LENGTH, WLZ_VTX_3_SCALE, WLZ_VTX_3_SET, WLZ_VTX_3_SUB, and WlzErrorFromAlg().

Referenced by WlzCMeshComputeCurvaturesFromNodNorm().

◆ WlzGeometryLSqOPlane()

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.

Returns
Woolz error code.
            The orthogonal distance regression plane is an eigenvector
            problem. This solution is based on one from
            http://mathforum.org which is credited to Doctor George.
            Starting with the distance from a point to a plane we
            wish to find \form#352 such as to minimise

\[ f(a,b,c,d) = \sum{ \frac{a x_i + b y_i + c z_i} {a^2 + b^2 + c^2}} \]

setting \(\frac{\partial f}{\partial d} = 0\) gives

\[ d = -(a x_0 + b y_0 + c z_0) \]

where \((x_0, y_0, z_0)\) is the centroid of the vertices. Using the centroid

\[ f(a,b,c,d) = \sum{ \frac{|a (x_i - x_0) + b (y_i - y_0) + c (z_i - z_0) |} {a^2 + b^2 + c^2}} \]

Define \(v\) \(M\) such that:

\[ v^T = [a \, b \, c] \]

\[ M = \left[ \begin{array}{ccc} x_1 - x_0 & y_1 - y_0 & z_1 - z_0 \\ x_2 - x_0 & y_2 - y_0 & z_2 - z_0 \\ \cdots & \cdots & \cdots \\ x_n - x_0 & y_n - y_0 & z_n - z_0 \end{array} \right] \]

\[ f(v) = \frac{(v^T M^T)(M v)}{v^T v} \]

\[ f(v) = \frac{v^T (M^T M) v}{v^T v} \]

Define \(A\)

\[ A = M^T M \]

The Rayleigh Quotient \(f(v)\) is minimised by the eigenvector of \(A\). However there is no need to compute the eigenvectors of \(A\). The SVD of \(M\) is

\[ M = U S V^T \]

where \(S\) is a diagonal vector containing the singular values of \(M\). The columns of \(V\) are it's singular vectors and \(U\) is an orthogonal matrix.

\[ A = M M^T \]

\[ A = (U S V^T)^T (U S V^T) \]

\[ A = (V S^T U^T) (U S V^T) \]

\[ A = V S^2 V^T \]

The decomposition of \(A\) diagonalises the matrix and gives an eigenvector decomposition. It means that the eigenvectors of \(A\) are the squares of the singular values of \(M\) and the eigenvectors of \(A\) are the singular vectors of \(M\). \(M\) is the covarience matrix.
Parameters
dstNrmDestination pointer for the plane normal.
dstCenDestination pointer for the centroid of the vertices which is on the plane.
nVtxumber of vertices.
vtxVector of vertices.

References AlcFree(), AlcMalloc(), AlgMatrixFree(), AlgMatrixRectNew(), AlgMatrixSVDecomp(), AlgMatrixZero(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrix::rect, _WlzGeomPolyListItem2D::v, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_VTX_3_SUB, WlzCentreOfMassVtx3D(), and WlzErrorFromAlg().

◆ WlzGeomTriangleVtxDistSq3D()

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:

Returns
Square of the minimum distance from the test vertex to the triangle.
                      l1
                      ^ 
           \ R2 |
            \   |
             \  |
              \ |
               \|
                *v2
                |\      
                | \ 
             R3 |  \  R1
                |   \
                |    \ 
                |  R0 \ 
                |      \ v1
         -------*-------*-------> l0
                |v0      \ 
             R4 |    R5   \  R6
                      |          \

\[ T(l_0,l_1) = \vec{v}_0 + l_0 \vec{u}_1 + l_1 \vec{u}_2 \]

where \(\vec{u}_1 = \vec{v}_1 - \vec{v}_0\) and \(\vec{u}_2 = \vec{v}_2 - \vec{v}_0.\) Inside the triangle

\[ l_0,l_1,(l_0 + l_1) \in [0-1] \]

The squared distance from some point \(\vec{u_T}\) with origin at \(\vec{v_0}\) is then found by minimising

\[ \|\vec{T}(l_0,l_1) - \vec{u}_T\|^2 \]

As this is a quadratic in \(L_0\) and \(l_1\) it can easily be minimised.

Parameters
dstPTDestination pointer for the position of vertex in the triangle that is closest to the test vertex .
dstZTDestination pointer, the value of which will be set to a non-zero value if the triangle has zero area. May be NULL.
dstITDestination pointer, the value of which will be set to a non-zero value if the projected test vertex is within the triangle, ie in region zero. May be NULL.
dstL0Destination pointer for the first triangle parameter ( \(l_0\)), may be NULL.
dstL1Destination pointer for the second triangle parameter ( \(l_1\)), may be NULL.
vTThe position of the test vertex.
v0First vertex of the triangle.
v1Second vertex of the triangle.
v2Third vertex of the triangle.

References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_DOT, and WLZ_VTX_3_SUB.

Referenced by WlzCMeshClosePointDom2D5(), and WlzGeomTetrahedronVtxDistSq3D().

◆ WlzGeomTriangleVtxDistSq2D()

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().

Returns
Signed squared distance from test vertex to triangle edge, with the sign such that it's negative if the test vertex is inside the triangle, zero if it's on an edge and positive if the vertex is outside the triangle.

The method used is first to test whether the vertex is outside, on or inside the triangle and then compute the closest point on each edge segment and the squared distance from the test vertex to that point.

Distance to an edge segment is computed using a parametric representation of the triangle edge segments of the form \( Q(t) = t (\vec{v_1} - \vec{v_0}) + \vec{v_0} \) then with \( \vec{v} = \vec{v_1} - \vec{v_0} \) and \( \vec{u} = \vec{v_T} - \vec{v_0} \) \( t = \frac{\vec{u} \cdot \vec{v}}{\|\vec{v}\|^2} \) but with \( t \) clipped to the interval [0-1]. The closest point on the segment is then at \( t \vec{v} + \vec{v_0} \) this is returned if the destination pointer is non-null.

Parameters
dstUDestination pointer for the closest point on a triangle edge, may be NULL.
dstEIDestination pointer for the index of the triangle edge on which the closest point lies, with 0 for v0 - v1, 1 for v1 - v2, 2 for v2 - v0. May be NULL.
vTTest vertex position.
v0First vertex of the triangle.
v1Second vertex of the triangle.
v2Third vertex of the triangle.

References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, WLZ_CLAMP, WLZ_VTX_2_DOT, WLZ_VTX_2_SCALE_ADD, WLZ_VTX_2_SQRLEN, WLZ_VTX_2_SUB, WlzGeomVtxOnLineSegment2D(), and WlzGeomVxInTriangle2D().

◆ WlzGeomTetrahedronVtxDistSq3D()

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().

Returns
Signed squared distance from test vertex to tetrahedron face with the sign such that it's negative if the test vertex is inside the tetrahedron, zero if it's on a face and positive if the vertex is outside the tetrahedron.
Parameters
dstUDestination pointer for the closest point on a tetrahedron face, may be NULL.
dstFIDestination pointer for the index of the tetrahedron face on which the closest point lies. May be NULL.
vTTest vertex position.
v0First vertex of the tetrahedron.
v1Second vertex of the tetrahedron.
v2Third vertex of the tetrahedron.
v3Fourth vertex of the tetrahedron.

References _WlzGeomPolyListItem2D::v, WLZ_MESH_TOLERANCE, WlzGeomTriangleVtxDistSq3D(), WlzGeomVxInTetrahedron(), and WlzGeomVxInTriangle3D().

◆ WlzGeomPolyTriangulate2D()

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.

Returns
Woolz error code.
            The algorithm is based on the paper:
            Kong, X., H. Everett, and G.T. Toussaint. 
            The Graham scan triangulates simple polygons. 
            Pattern Recognition Letters. November 1990, 713-716.

            In this algorithm the input is a simple polygon stored
            as a doubly linked circular list. With next(Pi) and
            prev(Pi) the successor and predecessor of Pi respectively.
            The algorithm produces a set D of diagonals comprising a
            triangulation of P. R is a set containing all the concave
            vertices of P.
1.  pi = p2;
2.  while (pi==Po) do
3.  if(IsAnEar(P,R, prev(pi)) and P is not a triangle // prev(Pi) is an ear
4.  D = D U (prev (prev (Pi)), Pi) // Store a diagonal
5.  P = P - prev(Pi) // Cut the ear
6.  if Pi in R and Pi is a convex vertex // Pi has become convex
7.  R = R - Pi
8.  if prev(Pi) in R and prev(Pi) is a convex vertex // prev(Pi) now convex
9.  R = R - prev(pi)
10. if (prev(Pi)=P0) // next(P0) was cut
11. Pi = next(Pi) // Advance the scan
12. else pi = next(pi) // prev(pi) not an ear or P is triangle. Advance scan.
13. end while
Parameters
nPVtxNumber of polygon vertices.
pVtxThe polygon vertices.
dstNTriDestination pointer for the number of triangles.
dstTriIndices (triples) of the vertices for each triangle. This should be freed using AlcFree().

References AlcMalloc(), _WlzGeomPolyListItem2D::next, _WlzGeomPolyListItem2D::prev, _WlzGeomPolyListItem2D::v, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, and WLZ_ERR_PARAM_NULL.