Woolz Image Processing
Version 1.7.5
|
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... | |
WlzDVertex2 WlzGeomTriangleCen2D | ( | WlzDVertex2 | v0, |
WlzDVertex2 | v1, | ||
WlzDVertex2 | v2 | ||
) |
Computes the position of the centroid of the given triangle in 2D.
v0 | First vertex of the triangle. |
v1 | Second vertex of the triangle. |
v2 | Third vertex of the triangle. |
References WLZ_VTX_2_ADD3, and WLZ_VTX_2_SCALE.
Referenced by WlzCMeshElmClosestPosIn2D().
WlzDVertex3 WlzGeomTriangleCen3D | ( | WlzDVertex3 | v0, |
WlzDVertex3 | v1, | ||
WlzDVertex3 | v2 | ||
) |
Computes the position of the centroid of the given triangle in 3D.
v0 | First vertex of the triangle. |
v1 | Second vertex of the triangle. |
v2 | Third vertex of the triangle. |
References WLZ_VTX_3_ADD3, and WLZ_VTX_3_SCALE.
WlzDVertex3 WlzGeomTetrahedronCen3D | ( | WlzDVertex3 | v0, |
WlzDVertex3 | v1, | ||
WlzDVertex3 | v2, | ||
WlzDVertex3 | v3 | ||
) |
Computes the position of the centroid of the given tetrahedron in 3D.
v0 | First vertex of the tetrahedron. |
v1 | Second vertex of the tetrahedron. |
v2 | Third vertex of the tetrahedron. |
v3 | Fourth vertex of the tetrahedron. |
References WLZ_VTX_3_ADD4, and WLZ_VTX_3_SCALE.
Referenced by WlzCMeshElmClosestPosIn3D().
int WlzGeomTriangleCircumcentre | ( | WlzDVertex2 * | ccVx, |
WlzDVertex2 | vx0, | ||
WlzDVertex2 | vx1, | ||
WlzDVertex2 | vx2 | ||
) |
Computes the circumcentre of the given triangle.
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.
ccVx | Destination pointer for the circumcentre. |
vx0 | First vertex of triangle. |
vx1 | Second vertex of triangle. |
vx2 | Third vertex of triangle. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzMeshElemSplit().
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.
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).
p0 | First vertex of triangle. |
p1 | Second vertex of triangle. |
p2 | Third vertex of triangle. |
pP | Given vertex. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.
Referenced by WlzCMeshElmEnclosesPos2D(), WlzGeomPolyTriangulate2D(), WlzGeomTriangleTriangleIntersect2D(), WlzGeomTriangleVtxDistSq2D(), and WlzMeshElemFindVxWalk().
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()).
v0 | First vertex of triangle. |
v1 | Second vertex of triangle. |
v2 | Third vertex of triangle. |
vQ | Given query vertex. |
vPMax | Maximum 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().
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.
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\]
.
v0 | First vertex of tetrahedron. |
v1 | Second vertex of tetrahedron. |
v2 | Third vertex of tetrahedron. |
v3 | Fourth vertex of tetrahedron. |
vP | Given vertex. |
References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
Referenced by WlzCMeshElmEnclosesPos3D(), and WlzGeomTetrahedronVtxDistSq3D().
WlzLong WlzGeomTriangleSnArea2I | ( | WlzIVertex2 | vx0, |
WlzIVertex2 | vx1, | ||
WlzIVertex2 | vx2 | ||
) |
Computes 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.vx0 | First vertex of triangle. |
vx1 | Second vertex of triangle. |
vx2 | Third vertex of triangle. |
References _WlzIVertex2::vtX, and _WlzIVertex2::vtY.
double WlzGeomTriangleSnArea2 | ( | WlzDVertex2 | vx0, |
WlzDVertex2 | vx1, | ||
WlzDVertex2 | vx2 | ||
) |
Computes 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.vx0 | First vertex of triangle. |
vx1 | Second vertex of triangle. |
vx2 | Third vertex of triangle. |
References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzCMeshBoundConform2D(), WlzCMeshCmpElmFeat2D(), WlzCMeshElmSnArea22D(), WlzCMeshFixNegativeElms2D(), WlzCMeshSetElm2D(), WlzCMeshToDomObjValues(), WlzGeometryTrackUpAndDown_s(), WlzGeomPolyTriangulate2D(), WlzGMEdgeTInsertRadial(), WlzMeshClosestNod2D(), WlzMeshElemFindVxWalk(), WlzMeshElemReplace1(), WlzMeshElemReplaceNWithN(), WlzMeshElemSplit(), WlzMeshElemVerify(), WlzMeshFromObjBox(), WlzMeshIDomAdd(), WlzMeshTransformAdapt(), WlzMeshTransformVtx(), WlzMeshVxVecAdd(), and WlzVertexHeapSortIdxFnD3().
double WlzGeomTetraSnVolume6 | ( | WlzDVertex3 | vx0, |
WlzDVertex3 | vx1, | ||
WlzDVertex3 | vx2, | ||
WlzDVertex3 | vx3 | ||
) |
Computes 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.vx0 | First vertex of tetrahedron. |
vx1 | Second vertex of tetrahedron. |
vx2 | Third vertex of tetrahedron. |
vx3 | Forth vertex of tetrahedron. |
References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
Referenced by WlzCMeshBoundConform3D(), WlzCMeshCmpElmFeat3D(), WlzCMeshElmSnVolume63D(), WlzCMeshExpansion(), WlzCMeshFixNegativeElms3D(), WlzCMeshSetElm3D(), and WlzGeomTetraInSphereDiam().
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.
v0 | First vertex of tetrahedron. |
v1 | Second vertex of tetrahedron. |
v2 | Third vertex of tetrahedron. |
v3 | Forth vertex of tetrahedron. |
References WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SUB, and WlzGeomTriangleArea2Sq3I().
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.
v0 | First vertex of tetrahedron. |
v1 | Second vertex of tetrahedron. |
v2 | Third vertex of tetrahedron. |
v3 | Forth vertex of tetrahedron. |
References WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SUB, and WlzGeomTriangleArea2Sq3().
WlzLong WlzGeomTriangleArea2Sq3I | ( | WlzIVertex3 | vx0, |
WlzIVertex3 | vx1, | ||
WlzIVertex3 | vx2 | ||
) |
Computes twice the square of the area of the given 3D 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}\).vx0 | First vertex of triangle \(\mathbf{v_0}\). |
vx1 | Second vertex of triangle \(\mathbf{v_1}\). |
vx2 | Third vertex of triangle \(\mathbf{v_2}\). |
References WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.
Referenced by WlzGeomTetVolZeroI().
double WlzGeomTriangleArea2Sq3 | ( | WlzDVertex3 | vx0, |
WlzDVertex3 | vx1, | ||
WlzDVertex3 | vx2 | ||
) |
Computes twice the square of the area of the given 3D 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}\).vx0 | First vertex of triangle \(\mathbf{v_0}\). |
vx1 | Second vertex of triangle \(\mathbf{v_1}\). |
vx2 | Third 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().
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.
vx0 | First vertex of triangle. |
vx1 | Second vertex of triangle. |
vx2 | Third vertex of triangle. |
gVx | Given vertex to test. |
References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzMeshQueConflictElem().
int WlzGeomLineSegmentsIntersect | ( | WlzDVertex2 | p0, |
WlzDVertex2 | p1, | ||
WlzDVertex2 | q0, | ||
WlzDVertex2 | q1, | ||
WlzDVertex2 * | dstN | ||
) |
Tests to see if the two given line segments intersect.
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.
p0 | 1st vertex of 1st line segment. |
p1 | 2nd vertex 1st line segment. |
q0 | 1st vertex of 2nd line segment. |
q1 | 2nd vertex of 2nd line segment. |
dstN | Destination 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, main(), _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzCMeshFMarNodes3D(), and WlzGeomTriangleTriangleIntersect2D().
int WlzGeomCmpAngle | ( | WlzDVertex2 | p0, |
WlzDVertex2 | p1 | ||
) |
Given two end connected 2D line segments this function compares the CCW angle of the segments.
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)\).
p0 | 1st segment endpoint vertex. |
p1 | 2nd segment endpoint vertex. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzGeomVtxEqual3D(), and WlzGMModelAddVertexToHT().
int WlzGeomVtxEqual2D | ( | WlzDVertex2 | pos0, |
WlzDVertex2 | pos1, | ||
double | tolSq | ||
) |
Checks to see if two verticies are the same within some tollerance.
pos0 | First node position. |
pos1 | Second node position. |
tolSq | Square of tollerance value. |
References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzContourRBFBndObj3D(), WlzGeomVtxOnLineSegment2D(), and WlzGMModelConstructSimplex2N().
int WlzGeomVtxEqual3D | ( | WlzDVertex3 | pos0, |
WlzDVertex3 | pos1, | ||
double | tol | ||
) |
Checks to see if two verticies are the same within some tollerance.
pos0 | First node position. |
pos1 | Second node position. |
tol | Tollerance value. |
References _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WlzGeomCmpAngle().
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.
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.
nV | Number of 3D verticies. |
vP | The 3D verticies. |
idxBuf | Buffer of nV indicies used for sorting the verticies. |
wP | Workspace with nV 2D verticies. |
rV | The radial vector. |
References AlgHeapSortIdx(), WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SCALE.
WlzDVertex3 WlzGeomTriangleNormal | ( | WlzDVertex3 | v0, |
WlzDVertex3 | v1, | ||
WlzDVertex3 | v2 | ||
) |
Computes the unit normal vector perpendicular to the triangle \(v_0, v_1, v_2\).
v0 | First vertex of triangle. |
v1 | Second vertex of triangle. |
v2 | Third 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(), WlzEffWriteObjStl(), WlzGeometryTrackUpAndDown_s(), and WlzGMVertexNormal3D().
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).
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.
t0 | First vertex of triangle. |
t1 | Second vertex of triangle. |
t2 | Third vertex of triangle. |
b0 | Minimum coordinates of axis aligned bounding box. |
b1 | Maximum coordinates of axis aligned bounding box. |
tst | Determines 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.
Referenced by WlzCMeshNewNod3D().
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).
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.
t0 | First vertex of tetrahedron. |
t1 | Second vertex of tetrahedron. |
t2 | Third vertex of tetrahedron. |
t3 | Fourth vertex of tetrahedron. |
b0 | Minimum coordinates of axis aligned bounding box. |
b1 | Maximum coordinates of axis aligned bounding box. |
tst | Determines 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_DOT, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.
Referenced by WlzCMeshNewNod3D().
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.
a | Plane X parameter. |
b | Plane Y parameter. |
c | Plane Z parameter. |
d | Other plane parameter. |
box | Axis 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().
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\).
a | Plane X parameter. |
b | Plane Y parameter. |
c | Plane Z parameter. |
d | Other plane parameter. |
p0 | First end point of the line segment. |
p1 | Second end point of the line segment. |
dstIsn | Destination pointer for point of intersection, may be NULL. |
References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
Referenced by 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.
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.
a | Plane X parameter. |
b | Plane Y parameter. |
c | Plane Z parameter. |
d | Other plane parameter. |
p0 | First triangle vertex. |
p1 | Second triangle vertex. |
p2 | Third triangle vertex. |
dstIsn0 | Destination pointer for first point of intersection, may be NULL. |
dstIsn1 | Destination 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().
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.
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}\).centre | Centre of ellipse. |
sAx | Ellipse semi axes, both \(> 0.0\). |
gPnt | Given point. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.
Referenced by WlzMatchICPWeightMatches().
unsigned int WlzGeomHashVtx3D | ( | WlzDVertex3 | pos, |
double | tol | ||
) |
Computes a hash value from a given 3D double precision position.
pos | Given position. |
tol | Tolerance, \(x = x \pm tol\). |
References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
Referenced by WlzGMModelAddVertexToHT(), WlzGMModelMatchVertexG3D(), and WlzGMModelRemVertex().
unsigned int WlzGeomHashVtx2D | ( | WlzDVertex2 | pos, |
double | tol | ||
) |
Computes a hash value from a given 2D double precision position.
pos | Given position. |
tol | Tolerance, \(x = x \pm tol\). |
References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzGMModelMatchVertexG2D().
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.
pos0 | First vertex. |
pos1 | Second vertex. |
tol | Tolerance, \(x = x \pm tol\). |
References _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.
Referenced by WlzBasisFnIMQ3DFromCPts(), WlzBasisFnMQ3DFromCPts(), and WlzGMVertexCmpSign3D().
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.
pos0 | First vertex. |
pos1 | Second vertex. |
tol | Tolerance, \(x = x \pm tol\). |
References _WlzDVertex2::vtX, _WlzDVertex2::vtY, and WLZ_VTX_2_SUB.
Referenced by WlzBasisFnGauss2DFromCPts(), WlzBasisFnIMQ2DFromCPts(), WlzBasisFnMQ2DFromCPts(), WlzBasisFnTPS2DFromCPts(), and WlzGMVertexCmpSign2D().
WlzDVertex2 WlzGeomUnitVector2D | ( | WlzDVertex2 | vec | ) |
Computes the 2D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\).
vec | Given vector, \(\mathbf{v}\). |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SCALE.
Referenced by WlzContourRBFBndObj3D().
WlzDVertex3 WlzGeomUnitVector3D | ( | WlzDVertex3 | vec | ) |
Computes the 3D unit vector \(\frac{1}{|\mathbf{v}|} \mathbf{v}\).
vec | Given vector, \(\mathbf{v}\). |
References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SCALE.
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.
v1 | Position of vertex, \(\mathbf{p}_1\). |
v0 | Position 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.
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.
v1 | Position of vertex, \(\mathbf{p}_1\). |
v0 | Position 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.
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\).
lPos0 | Vertex at one end of line segment. |
lPos1 | Vertex at other end of line segment. |
pos | Position of point. |
References WLZ_VTX_2_DOT, and WLZ_VTX_2_SUB.
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.
step | Spiral step count. |
Referenced by WlzCMeshClosestNod2D(), WlzCMeshClosestNod2D5(), WlzCMeshClosestNod3D(), WlzCMeshElmClosestPosIn2D(), and WlzCMeshElmClosestPosIn3D().
int WlzGeomItrSpiral2I | ( | int | step, |
int * | pX, | ||
int * | pY | ||
) |
Iterates the given positions coordinates through a 2D expanding integer spiral.
step | Spiral step count, must be zero when this function is called for the for first step. |
pX | Destination pointer for column coordinate. |
pY | Destination pointer for line coordinate. |
Referenced by WlzCMeshClosestNod2D(), and WlzCMeshElmClosestPosIn2D().
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.
step | Spiral step count. |
References cbrt().
int WlzGeomItrSpiral3I | ( | int | step, |
int * | pX, | ||
int * | pY, | ||
int * | pZ | ||
) |
Iterates the given positions coordinates through a 3D expanding integer spiral.
step | Spiral step count, must be zero when this function is called for the for first step. |
pX | Destination pointer for column coordinate. |
pY | Destination pointer for line coordinate. |
pZ | Destination pointer for plane coordinate. |
References cbrt().
Referenced by WlzCMeshClosestNod2D5(), WlzCMeshClosestNod3D(), and WlzCMeshElmClosestPosIn3D().
double WlzGeomDistSq2D | ( | WlzDVertex2 | v0, |
WlzDVertex2 | v1 | ||
) |
Computes square of the Euclidean distance between the given two vertices.
v0 | First of the given vertices. |
v1 | Second of the given vertices. |
References WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.
Referenced by WlzCMeshElmMaxSqEdgLen2D(), and WlzCMeshUpdateMaxSqEdgLen2D().
double WlzGeomDistSq3D | ( | WlzDVertex3 | v0, |
WlzDVertex3 | v1 | ||
) |
Computes square of the Euclidean distance between the given two vertices.
v0 | First of the given vertices. |
v1 | Second of the given vertices. |
References WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.
Referenced by WlzCMeshElmMaxSqEdgLen2D5(), WlzCMeshElmMaxSqEdgLen3D(), WlzCMeshUpdateMaxSqEdgLen2D5(), and WlzCMeshUpdateMaxSqEdgLen3D().
double WlzGeomDist2D | ( | WlzDVertex2 | v0, |
WlzDVertex2 | v1 | ||
) |
Computes the Euclidean distance between the given two vertices.
v0 | First of the given vertices. |
v1 | Second of the given vertices. |
References WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SUB.
Referenced by WlzCMeshFMarNodes3D().
double WlzGeomDist3D | ( | WlzDVertex3 | v0, |
WlzDVertex3 | v1 | ||
) |
Computes the Euclidean distance between the given two vertices.
v0 | First of the given vertices. |
v1 | Second of the given vertices. |
References WLZ_VTX_3_LENGTH, and WLZ_VTX_3_SUB.
Referenced by WlzCMeshFMarNodes3D().
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.
xTr | Transform coordinates for x. |
yTr | Transform coordinates for y. |
dd | Twice the area of the source triangle. |
sVx | Source triangle vertices. |
dVx | Destination triangle vertices. |
thresh | Threshold value for twice the area. |
References _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzCMeshToDomObjValues(), WlzMeshFromObjBox(), and WlzMeshTransformVtx().
int WlzGeomTetraAffineSolveLU | ( | double * | tr, |
WlzDVertex3 * | sVx, | ||
WlzDVertex3 * | dVx | ||
) |
Computes the affine transform coefficients from the source to target tetrahedron.
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.tr | Transform matrix with 4x4 contiguous coefficients which are equivalent to the base storage of the matrix in a WlzAffineTransform. |
sVx | Source tetrahedron vertices. |
dVx | Destination tetrahedron vertices. |
References AlgMatrixLUInvertRaw4(), _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
int WlzGeomTetraAffineSolve | ( | double * | tr, |
WlzDVertex3 * | sVx, | ||
WlzDVertex3 * | dVx, | ||
double | thresh | ||
) |
Computes the affine transform coefficients from the source to target tetrahedron.
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.tr | Transform matrix with 4x4 contiguous coefficients which are equivalent to the base storage of the matrix in a WlzAffineTransform. |
sVx | Source tetrahedron vertices. |
dVx | Destination tetrahedron vertices. |
thresh | Threshold value which is the lower limit of 6 x the source tetrahedron's volume. |
References _WlzDVertex3::vtX, _WlzDVertex3::vtY, and _WlzDVertex3::vtZ.
Referenced by WlzCMeshStrainTensorAtPts(), and WlzCMeshToDomObjValues().
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.
Given the line segment \form#240, \form#241 any position along the segment can be given by a parameter \form#169 (range [0-1]), where \form#277.
obj | Given object. Object must have a valid domain. |
p0 | First vertex. |
p1 | Second vertex. |
tol | Acceptable placement error. |
inside | Non-zero if the returned position should be inside or on the boundary, if zero it will be outside or on the boundary. |
method | Method 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. |
dstStat | Destination 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().
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.
Given the line segment \form#240, \form#241 any position along the segment can be given by a parameter \form#169 (range [0-1]), where \form#277.
obj | Given object. Object must have a valid domain. |
p0 | First vertex. |
p1 | Second vertex. |
tol | Acceptable placement error. |
inside | Non-zero if the returned position should be inside or on the boundary, if zero it will be outside or on the boundary. |
method | Method 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. |
dstStat | Destination 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().
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.
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).
vx0 | First vertex of tetrahedron. |
vx1 | Second vertex of tetrahedron. |
vx2 | Third vertex of tetrahedron. |
vx3 | Forth vertex of tetrahedron. |
References ALG_DBL_TOLLERANCE, WlzGeomTetraSnVolume6(), and WlzGeomTriangleArea2Sq3().
double WlzGeomTetraInSphereRegDiam | ( | double | side | ) |
Given the side length of a regular tetrahedron this function computes the maximum diameter of an inscribed sphere.
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.
side | Length of an edge. |
References ALG_DBL_TOLLERANCE, ALG_M_SQRT2, and ALG_M_SQRT3.
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:
(y) ^ pi/2 | | pi | <------+-------> (x) | 0 | 3pi/2| V
org | Position of the origin. |
dst | Position of the destination. |
dstRad | Destination 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.
Referenced by WlzCentrality().
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.
v0 | First vertex. |
v1 | Second vertex (the common one). |
v2 | Third vertex. |
References ALG_DBL_TOLLERANCE, WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.
Referenced by WlzCMeshFMarNodes3D().
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.
tst | Test vertex. |
seg0 | First vertex of line segment. |
seg1 | Second vertex of line segment. |
tol | Tollerance. |
References _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_SUB, WlzGeomVtxEqual2D(), _WlzDBox2::xMax, _WlzDBox2::xMin, _WlzDBox2::yMax, and _WlzDBox2::yMin.
Referenced by WlzGeomTriangleVtxDistSq2D().
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.
Consider a line segment from a vertex at \form#280 to another vertex at \form#281 , with a third test vertex at \form#282 , 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.pX | Test vertex. |
p0 | First vertex of line segment. |
p1 | Second vertex of line segment. |
dstN | Destination 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().
double WlzGeomArcLength2D | ( | WlzDVertex2 | a, |
WlzDVertex2 | b, | ||
WlzDVertex2 | c | ||
) |
Computes the arc length from a to b traveling CCW on a circle with centre c.
a | Start point. |
b | End point. |
c | Cirecle 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.
int WlzGeomRectFromWideLine | ( | WlzDVertex2 | s, |
WlzDVertex2 | t, | ||
double | w, | ||
WlzDVertex2 * | v | ||
) |
Computes the coordinates of vertices that may be used to draw a wide rectangle.
Given two vertices \form#290, \form#26 which define a line segment and width \form#291 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..s | First vertex of line segment. |
t | Second vertex of line segment. |
w | line width. |
v | Destination pointer for the four vertices of the rectangle. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_ADD, WLZ_VTX_2_LENGTH, and WLZ_VTX_2_SUB.
Referenced by WlzDrawDomainObj().
WlzDVertex3 WlzGeomLinePlaneIntersection | ( | WlzDVertex3 | v, |
WlzDVertex3 | p0, | ||
WlzDVertex3 | p1, | ||
WlzDVertex3 | p2, | ||
WlzDVertex3 | p3, | ||
int * | dstPar | ||
) |
Computes the intersection of a line with a plane.
Computes the intersection of line (vector \form#254) through an off plane vertex ( \form#303) 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 \]
v | Unit vector for ray. |
p0 | First point on the plane. |
p1 | Second point on the plane. |
p2 | Third point on the plane. |
p3 | Off plane point that ray passes through. |
dstPar | Destination value set to 1 if the vector is parrallel to the plane, must not be NULL. |
References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.
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.
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#318 and \form#6 are only set if the line passes through the triangle.
org | Line origin, \(O\). |
dir | Line direction, \(D\) (does not need to be a unit vector). |
v0 | First vertex on triangle. |
v1 | Second vertex on the triangle |
v2 | Third vertex on the triangle |
dstPar | Destination pointer for flag set to 1 if the vector is parrallel to the plane, may be NULL. |
dstT | Destination pointer for t parameter, may be NULL. |
dstU | Destination pointer for u parameter, may be NULL. |
dstV | Destination pointer for v parameter, may be NULL. |
References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_SUB, and WlzGeomLineLineSegmentIntersect3D().
Referenced by WlzCMeshFMarNodes3D().
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.
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})\).r0 | A vertex on the line. |
rD | Direction of the line. |
p0 | First end vertex of the line segment. |
p1 | Second end vertex of the line segment. |
dstN | Destination 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().
Referenced by WlzGeomLineTriangleIntersect3D().
int WlzGeomVtxOnLine3D | ( | WlzDVertex3 | p0, |
WlzDVertex3 | r0, | ||
WlzDVertex3 | rD | ||
) |
Tests whether a vertex is a line.
Tests whether the given vertex \form#240 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 \]
p0 | Given vertex. |
r0 | Point on line. |
rD | Direction of line. |
References ALG_DBL_TOLLERANCE, WLZ_VTX_3_CROSS, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.
Referenced by WlzGeomLineLineSegmentIntersect3D().
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.
p0 | First vertex of triangle. |
p1 | Second vertex of triangle. |
p2 | Third vertex of triangle. |
pX | Given position, which is within (or on) the triangle. |
lambda | Given array for the three lambda values. |
References ALG_DBL_TOLLERANCE, _WlzDVertex2::vtX, and _WlzDVertex2::vtY.
Referenced by WlzCMeshValueTransfer(), and 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.
p0 | First vertex of triangle. |
p1 | Second vertex of triangle. |
p2 | Third vertex of triangle. |
v0 | Value at first vertex of triangle. |
v1 | Value at second vertex of triangle. |
v2 | Value at third vertex of triangle. |
pX | Given position, which is within (or on) the triangle. |
References WlzGeomBaryCoordsTri2D().
Referenced by WlzBasisFnValueMOSPhi(), WlzCMeshCurvToImage(), WlzCMeshMeshMeshProduct(), WlzCMeshProduct(), and WlzCMeshToDomObjValues().
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.
n | Number of polygon vertices, which is the same as the number of values and generalised barycentric coordinates. |
p | The vertices of the convex polygon (ordered counter-clockwise). |
v | Values at the corresponding polygon vertices. |
w | Used 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. |
q | Given position, which is must be within or on the convex hull of the polygon. |
dstErr | Destination error pointer, may be NULL. |
References AlcFree(), AlcMalloc(), WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WlzConvHullClarkson2D(), and 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.
n | Number of polygon vertices, which is the same as the number of values and generalised barycentric coordinates. |
p | The vertices of the convex polygon (ordered counter-clockwise). |
v | Values at the corresponding polygon vertices. |
w | Used to compute and return the generalised barycentric coordinates of the given position. |
q | Given position, which is must be in or on the polygon. |
References _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_VTX_2_DOT, WLZ_VTX_2_LENGTH, WLZ_VTX_2_SQRLEN, and WLZ_VTX_2_SUB.
Referenced by WlzGeomInterpolatePoly2D().
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.
p0 | First vertex of tetrahedron. |
p1 | Second vertex of tetrahedron. |
p2 | Third vertex of tetrahedron. |
p3 | Fourth vertex of tetrahedron. |
pX | Given position, which is within (or on) the tetrahedron. |
lambda | Given array for the four lambda values. |
References ALG_DBL_TOLLERANCE, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, and WLZ_VTX_3_SUB.
Referenced by WlzBasisFnValueMOSPhi(), WlzCMeshValueTransfer(), and 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.
p0 | First vertex of tetrahedron. |
p1 | Second vertex of tetrahedron. |
p2 | Third vertex of tetrahedron. |
p3 | Fourth vertex of tetrahedron. |
v0 | Value at first vertex of tetrahedron. |
v1 | Value at second vertex of tetrahedron. |
v2 | Value at third vertex of tetrahedron. |
v3 | Value at fourth vertex of tetrahedron. |
pX | Given position, which is within (or on) the tetrahedron. |
References WlzGeomBaryCoordsTet3D().
Referenced by WlzBasisFnValueMOSPhi(), WlzCMeshProduct(), and WlzCMeshToDomObjValues().
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).
p0 | First vertex of triangle (origin of the 2D plane). |
p1 | Second vertex of triangle (on the \(\mathbf{u}\) axis in the 2D plane). |
p2 | Third vertex of triangle. |
dstQ1 | Destination pointer for the second vertex in the plane, must not be NULL. |
dstQ2 | Destination 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.
Referenced by WlzCMeshFMarNodes3D().
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).
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.
t0 | First vertex of the triangle. |
t1 | Second vertex of the triangle. |
t2 | Third vertex of the triangle. |
b0 | Minimum of the bounding box. |
b1 | Maximum of the bounding box. |
tst | Determines 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_DOT, WLZ_VTX_3_SQRLEN, and WLZ_VTX_3_SUB.
Referenced by WlzCMeshIntersectDom2D5(), WlzCMeshNewNod3D(), and WlzGeoModelGridWSpSet3D().
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.
See WlzGeomTriangleTriangleIntersect2DA().
s0 | First vertex of the first triangle. |
s1 | Second vertex of the first triangle. |
s2 | Third vertex of the first triangle. |
t0 | First vertex of the second triangle. |
t1 | Second vertex of the second triangle. |
t2 | Third vertex of the second triangle. |
References WlzGeomLineSegmentsIntersect(), and WlzGeomVxInTriangle2D().
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).
See WlzGeomTriangleTriangleIntersect3DA().
s0 | First vertex of the first triangle. |
s1 | Second vertex of the first triangle. |
s2 | Third vertex of the first triangle. |
t0 | First vertex of the second triangle. |
t1 | Second vertex of the second triangle. |
t2 | Third vertex of the second triangle. |
References ALG_DBL_TOLLERANCE, _WlzGeomPolyListItem2D::v, _WlzDVertex2::vtX, _WlzDVertex3::vtX, _WlzDVertex2::vtY, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_VTX_3_CROSS, WLZ_VTX_3_DOT, WLZ_VTX_3_FABS, and WLZ_VTX_3_SUB.
Referenced by WlzGMModelCutDom().
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.
nC | Number of curvatures required: 1 computes the Gaussian curvature, 2 computes both principle curvatures. |
dstC | Destination pointer for the curvature value(s), must not be NULL. If nC == 2 then the values atr Gaussian followed by mean curvature. |
nrm | Normal at the first vertex. |
nV | Number of vertices, must be >= 3. |
vtx | Array 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().
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.
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#345 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.dstNrm | Destination pointer for the plane normal. |
dstCen | Destination pointer for the centroid of the vertices which is on the plane. |
nVtx | umber of vertices. |
vtx | Vector 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().
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:
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.
dstPT | Destination pointer for the position of vertex in the triangle that is closest to the test vertex . |
dstZT | Destination pointer, the value of which will be set to a non-zero value if the triangle has zero area. May be NULL. |
dstIT | Destination 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. |
dstL0 | Destination pointer for the first triangle parameter ( \(l_0\)), may be NULL. |
dstL1 | Destination pointer for the second triangle parameter ( \(l_1\)), may be NULL. |
vT | The position of the test vertex. |
v0 | First vertex of the triangle. |
v1 | Second vertex of the triangle. |
v2 | Third 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().
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().
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.
dstU | Destination pointer for the closest point on a triangle edge, may be NULL. |
dstEI | Destination 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. |
vT | Test vertex position. |
v0 | First vertex of the triangle. |
v1 | Second vertex of the triangle. |
v2 | Third 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().
Referenced by WlzCMeshElmClosestPosIn2D().
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().
dstU | Destination pointer for the closest point on a tetrahedron face, may be NULL. |
dstFI | Destination pointer for the index of the tetrahedron face on which the closest point lies. May be NULL. |
vT | Test vertex position. |
v0 | First vertex of the tetrahedron. |
v1 | Second vertex of the tetrahedron. |
v2 | Third vertex of the tetrahedron. |
v3 | Fourth vertex of the tetrahedron. |
References _WlzGeomPolyListItem2D::v, WLZ_MESH_TOLERANCE, WlzGeomTriangleVtxDistSq3D(), WlzGeomVxInTetrahedron(), and WlzGeomVxInTriangle3D().
Referenced by WlzCMeshElmClosestPosIn3D().
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.
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
nPVtx | Number of polygon vertices. |
pVtx | The polygon vertices. |
dstNTri | Destination pointer for the number of triangles. |
dstTri | Indices (triples) of the vertices for each triangle. This should be freed using AlcFree(). |
References AlcFree(), AlcMalloc(), _WlzGeomPolyListItem2D::cvx, _WlzGeomPolyListItem2D::next, _WlzGeomPolyListItem2D::prev, _WlzGeomPolyListItem2D::v, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WlzGeomTriangleSnArea2(), and WlzGeomVxInTriangle2D().
Referenced by WlzCMeshElmFuse2D().