Woolz Image Processing Version 1.4.0
WlzFeatures

Data Structures

struct  _WlzPoints
 An array of either 2D or 3D points which may have either integral of floating point values. Possible types are: WLZ_POINTS_2I, WLZ_POINTS_2D, WLZ_POINTS_3I and WLZ_POINTS_3D. Typedef: WlzPoints. More...
struct  _WlzRect
 An integer rectangle domain. Side from (l[0],k[0]) to (l[1],k[1]) is a long side. The vertices are cyclic. Typedef: WlzIRect. More...
struct  _WlzFRect
 A single precision floating point rectangle domain. Side from (l[0],k[0]) to (l[1],k[1]) is a long side. The vertices are cyclic. Typedef: WlzFRect. More...

Files

file  WlzArea.c
 

Computes the area of an object.


file  WlzBoundingBox.c
 

Functions for computing the axis aligned bounding box of objects.


file  WlzCCor.c
 

Computes the cross correlation of two objects.


file  WlzCentrality.c
 

Functions for computing the centrality of a feature domain with respect to a boundary domain. See WlzCentrality().


file  WlzCentreOfMass.c
 

Computes the centre of mass of objects.


file  WlzDistMetric.c
 

Functions to compute the Hausdorff distance, mean nearest neighbour and the median nearest neighbour distances between two datasets.


file  WlzGreyStats.c
 

Calculates simple statistics about an object's grey values.


file  WlzGreyValueMixing_s.c
 

Functions to mix the grey values of two woolz object and produce a new object using

\[ o = (1 - x) o_1 + x o_2 \]

.


file  WlzLineArea.c
 

Computes the line area of an object.


file  WlzNMSuppress.c
 

A non-maximal supression filter, which constructs a new domain object using a Canny-like non-maximal suppression algorithm. The domain is the non-maximally suppressed domain and the values are the encoded gradient direction.


file  WlzNObjGreyStats.c
 

Calculates statistics for grey values across the objects of the given compound object.


file  WlzPoints.c
 

Functions for handling point domains.


file  WlzRGBAGreyStats.c
 

Calculates simple quick statistics for a domain object with RGBA values.


file  WlzSampleValuesAndCoords.c
 

Extracts values and coordinates from a Woolz object with grey values using a sampoling function.


file  WlzScalarFeatures.c
 

Functions for extracting scalar features from objects.


file  WlzTensor.c
 

Functions which derive and manipulate tensor quantities.


file  WlzVerticies.c
 

Functions for extracting vertices from objects represented by vertices, eg polylines, boundlists and contours.


Enumerations

enum  _WlzScalarFeatureType {
  WLZ_SCALARFEATURE_VALUE,
  WLZ_SCALARFEATURE_GRADIENT
}
 Scalar features of objects. More...

Functions

WlzLong WlzArea (WlzObject *obj, WlzErrorNum *dstErr)
 Computes the area of an object.
WlzIBox2 WlzBoundingBox2I (WlzObject *inObj, WlzErrorNum *dstErr)
 Computes the 2D integer axis aligned bounding box of the given object.
WlzDBox2 WlzBoundingBox2D (WlzObject *inObj, WlzErrorNum *dstErr)
 Computes the 2D double precision axis aligned bounding box of the given object.
WlzIBox3 WlzBoundingBox3I (WlzObject *inObj, WlzErrorNum *dstErr)
 Computes the integer 3D axis aligned bounding box of the given object.
WlzDBox3 WlzBoundingBox3D (WlzObject *inObj, WlzErrorNum *dstErr)
 Computes the double precision 3D axis aligned bounding box of the given object.
WlzDBox3 WlzBoundingBoxGModel3D (WlzGMModel *model, WlzErrorNum *dstErr)
 Computes the 3D axis aligned bounding box of the 3D geometric model. It is an error if the model is not a 3D model.
WlzDBox2 WlzBoundingBoxGModel2D (WlzGMModel *model, WlzErrorNum *dstErr)
 Computes the 2D axis aligned bounding box of the 2D geometric model. It is an error if the model is not a 2D model.
WlzIBox2 WlzBoundingBox2DTo2I (WlzDBox2 bBox2D)
 Converts a 2D double precision bounding box to an integer bounding box. There is no checking done for underflows or overflows.
WlzDBox2 WlzBoundingBox2ITo2D (WlzIBox2 bBox2I)
 Converts a 2D integer bounding box to an double precision bounding box.
WlzIBox3 WlzBoundingBox3DTo3I (WlzDBox3 bBox3D)
 Converts a 3D double precision bounding box to an integer bounding box. There is no checking done for underflows or overflows.
WlzDBox3 WlzBoundingBox3ITo3D (WlzIBox3 bBox3I)
 Converts a 3D integer bounding box to an double precision bounding box.
WlzFBox3 WlzBoundingBox3DTo3F (WlzDBox3 bBox3D)
 Converts a 3D double bounding box to a single precision bounding box.
WlzDBox3 WlzBoundingBox3FTo3D (WlzFBox3 bBox3F)
 Converts a 3D single bounding box to a double precision bounding box.
WlzIBox2 WlzBoundingBoxUnion2I (WlzIBox2 box0, WlzIBox2 box1)
 Computes the 2D integer bounding box which encloses the given pair of bounding boxes.
WlzFBox2 WlzBoundingBoxUnion2F (WlzFBox2 box0, WlzFBox2 box1)
 Computes the 2D single precision bounding box which encloses the given pair of bounding boxes.
WlzDBox2 WlzBoundingBoxUnion2D (WlzDBox2 box0, WlzDBox2 box1)
 Computes the 2D double precision bounding box which encloses the given pair of bounding boxes.
WlzIBox3 WlzBoundingBoxUnion3I (WlzIBox3 box0, WlzIBox3 box1)
 Computes the 3D integer bounding box which encloses the given pair of bounding boxes.
WlzFBox3 WlzBoundingBoxUnion3F (WlzFBox3 box0, WlzFBox3 box1)
 Computes the 3D single precision bounding box which encloses the given pair of bounding boxes.
WlzDBox3 WlzBoundingBoxUnion3D (WlzDBox3 box0, WlzDBox3 box1)
 Computes the 3D double precision bounding box which encloses the given pair of bounding boxes.
double WlzCCorS2D (WlzObject *obj0, WlzObject *obj1, int unionFlg, int normFlg, WlzErrorNum *dstErr)
 Computes the cross correlation of the two given 2D spatial domain objects in the spatial domain. The algorithm used by this function is not an efficient unless a single cross correlation value is required.
double WlzCentrality (WlzObject *fObj, WlzObject *bObj, int nRay, int binFlg, double *dstMaxR, WlzErrorNum *dstErr)
 Computes the centrality of a feature domain with respect to a boundary domain. The boundary domain must enclose the feature domain.
WlzDVertex2 WlzCentreOfMass2D (WlzObject *srcObj, int binObjFlag, double *dstMass, WlzErrorNum *dstErr)
 Calculates the centre of mass of a Woolz object. If the given object does not have grey values or the binary object flag is set (ie non zero) then every pixel or vertex within the object's domain has the same mass.
WlzDVertex3 WlzCentreOfMass3D (WlzObject *srcObj, int binObjFlag, double *dstMass, WlzErrorNum *dstErr)
 Calculates the centre of mass of a Woolz object. If the given object does not have grey values or the binary object flag is set (ie non zero) then every pixel or vertex within the objects domain has the same mass.
WlzDVertex2 WlzCentreOfMassVtx2D (int nVtx, WlzDVertex2 *vtx)
 Computes the centre of mass of a vector of 2D vertices.
WlzDVertex3 WlzCentreOfMassVtx3D (int nVtx, WlzDVertex3 *vtx)
 Computes the centre of mass of a vector of 3D vertices.
WlzErrorNum WlzDistMetricGM (WlzGMModel *model0, WlzGMModel *model1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the vertices of the given geometric models. See WlzDistMetricVertex2D() for details of the metrics.
WlzErrorNum WlzDistMetricDirGM (WlzGMModel *model0, WlzGMModel *model1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the vertices of the given geometric models. See WlzDistMetricDirVertex2D() for details of the metrics.
WlzErrorNum WlzDistMetricVertex2D (int n0, WlzDVertex2 *vx0, int n1, WlzDVertex2 *vx1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. Each of these distance measures is the maximum of the two possible directed measures:

\[ D = \max{(d(A,B), d(B,A))} \]

Where $D$ is a non-directed distance metric and $d$ is the associated directed distance metric. $A$ and $B$ are the two datasets.

WlzErrorNum WlzDistMetricVertex3D (int n0, WlzDVertex3 *vx0, int n1, WlzDVertex3 *vx1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. See WlzDistMetricVertex2D() for an explaination of the distance metrics.
WlzErrorNum WlzDistMetricDirVertex2D (int n0, WlzDVertex2 *vx0, int n1, WlzDVertex2 *vx1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum mean neighbour distances between the given sets of vertices. The directed Hausdorff distance metric is:

\[ \max_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

The directed mean nearest neighbour distance metric is:

\[ \frac{1}{n_a} \sum_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

Likewise the directed median nearest neighbour distance metric is:

\[ \mathrm{median}_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

Where $A$ and $B$ are the two datasets.

WlzErrorNum WlzDistMetricDirVertex3D (int n0, WlzDVertex3 *vx0, int n1, WlzDVertex3 *vx1, double *dstDistH, double *dstDistM, double *dstDistN, double *dstDistI)
 Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. See WlzDistMetricDirVertex2D() for an explaination of the distance metrics.
WlzDVertex3WlzGeometryTrackUpAndDown_s (int numberOfPixelsZ, int startTrackingFile, int numberOfFilesDownOrUp, double disForInAndOutGuid, double disForInAndOut, unsigned char **TwoDImageFilesNameList, int numOf2DWlzFiles, int downOrUp, int sectionLength_N, int subSubSectionLength_L, int numberOfSampleP_k, char *surfacePointFileName, char *surfaceInPointFileName, char *surfaceOutPointFileName, int startShell, int endShell, int startSection, int endSection, double minDis, WlzErrorNum *dstErr)
 Track a curved path through a set of geometric model shells.
size_t WlzGreySize (WlzGreyType gType)
 Given a grey type, computes the number of bytes required to to store a single grey value of that type. Note the grey type WLZ_GREY_BIT is able to store multiple values in a single byte but one byte is required to store a single bit.
int WlzGreyStats (WlzObject *srcObj, WlzGreyType *dstGType, double *dstMin, double *dstMax, double *dstSum, double *dstSumSq, double *dstMean, double *dstStdDev, WlzErrorNum *dstErr)
 Calculates simple quick statistics for given 2D or 3D domain object with grey values. Pointers provided for results may be NULL without causing an error.
WlzObjectWlzGreyValueMixing_s (WlzObject *sObj, WlzObject *tObj, double xmiddle, WlzErrorNum *dstErr)
 calculate the distance map of a woolz obj
int WlzLineArea (WlzObject *obj, WlzErrorNum *dstErr)
 Calculate the line-area of an object defined as the sum of the line segments bounded by the left hand end of the first interval in a line and the right hand end of the last interval in that line.
WlzObjectWlzMakeMarkers (WlzVertexType vType, int nVtx, WlzVertexP vtx, WlzMarkerType mType, int mSz, WlzErrorNum *dstErr)
 Constructs a domain from the union of marker domains with a marker domain at each of the given vertex positions.
WlzObjectWlzMarkerLattice (WlzObject *gObj, WlzMarkerType mType, int mSz, int mSep, WlzErrorNum *dstErr)
 Creates a new domain object that is a formed from a lattice of markers covering the given domain.
WlzPointsWlzMakePoints (WlzObjectType type, int nVtx, WlzVertexP vtxP, int maxVtx, WlzErrorNum *dstErr)
 Creates a new point domain. A point domain consists of an array of vertices which are treated as seperate points.
WlzErrorNum WlzNObjGreyStats (WlzObject *gObj, int mean, int stddev, int *dstN, WlzObject **dstMinObj, WlzObject **dstMaxObj, WlzObject **dstSumObj, WlzObject **dstSSqObj)
 Computes a collection of objects from the object (which should be a compound array object). The computed objects all have the intersection of the given object domains as their domain and WLZ_GREY_DOUBLE values for the minimum, maximum, sum and sum of squares of the given grey values at each pixel/voxel within the intersection domain.
WlzErrorNum Wlz3DSectionOcc (WlzObject *obj, WlzThreeDViewStruct *vs, double sep, double *dstFirst, double *dstLast, int *dstArraySizeOcc, int **dstArrayOcc)
 Computes an array of integers which correspond to a line through the given object's domain perpendicular to the plane defined by the given view structure. At each point along the line the occupancy area or number of domains intersecting the plane is computed. When the given object is a 3D spatial domain object then areas are computed and when it is a compound array then the occupancy (number of domains) is computed.
WlzObjectWlzDomainOccupancy (WlzObject *gObj, WlzErrorNum *dstErr)
 Computes a new Woolz domain object which has the domain of the given object (union of domains if the object is a compound array) and a value table with integer values representing the number of domains present.
WlzObjectWlzPointsToDomObj (WlzPoints *pnt, double scale, WlzErrorNum *dstErr)
 Creates a domain object which coresponds to the union of the given points.
int WlzRGBAGreyStats (WlzObject *srcObj, WlzRGBAColorSpace colSpc, WlzGreyType *dstGType, double *dstMin, double *dstMax, double *dstSum, double *dstSumSq, double *dstMean, double *dstStdDev, WlzErrorNum *dstErr)
 Calculates simple quick statistics for the domain object with RGBA values. Each component has its statictics computed and entered into the four double[4] arrays.
WlzErrorNum WlzSampleValuesAndCoords (WlzObject *obj, WlzGreyType *dstGType, int *dstNVal, WlzGreyP *dstValP, WlzVertexP *dstCoords, WlzSampleFn samFn, int samFac)
 Allocates buffers for both the grey values and the coordinates of the grey values in the given object. On return these buffers contain the sampled object values and the coordinates of the values.
WlzErrorNum WlzScalarFeatures2D (WlzObject *obj, int *dstNFeat, WlzIVertex2 **dstFeat, WlzScalarFeatureType fType, WlzThresholdType thrHL, WlzPixelV thrV, double filterV, double minDist)
 Finds scalar features within the given object. Where the feature values are all either above or below the given threshold value, depending on the value of thrHL, and the features have the given minimum seperation distance.
WlzObjectWlzCMeshDGTensor (WlzObject *cObj, int invert, WlzErrorNum *dstErr)
 Given a conforming mesh transform this function computes the displacement gradient tensor for each of it's valid elements. Given displacement $\vec{u}(\vec{r})$ with position vector $\vec{r}$ which maps a point from a space $\vec{x}$ to a space $\vec{u}$ the displacement gradient tensor is defined in 3D as

\[ { \newcommand{\pd}[2]{\frac{\partial #1}{\partial #2}} u_{i,j} = \left [ \begin{array}{ccc} \pd{u_0}{x_0} & \pd{u_0}{x_1} & \pd{u_0}{x_2} \\ \pd{u_1}{x_0} & \pd{u_1}{x_1} & \pd{u_1}{x_2} \\ \pd{u_2}{x_0} & \pd{u_2}{x_1} & \pd{u_2}{x_2} \end{array} \right] } \]

with

\[ \Delta r_i = u_{ij} r_j \]

where $\vec{u} = \left[u_0, u_1, u_2\right]^T$ and $\Delta \vec{r} = \left[x_0, x_1, x_2\right]^T$. The displacement gradient tensor matrix is just the rotation and independent scaling part of the affine transform that displaces the element.

WlzObjectWlzCMeshDGTensorAtPts (WlzObject *cObj, int invert, WlzDVertex3 sd, int dither, WlzErrorNum *dstErr)
 Given a conforming mesh transform this function computes the displacement gradient tensor at regular cartesian grid sample points throughout the mesh. The tensor values at the sample points are computed at each point by computing an iverse distance weighted least squares general affine transform for the ring of nodes surrounding the closest node. See WlzCMeshDGTensor() for the description of the tensor.
WlzObjectWlzCMeshStrainTensor (WlzObject *cObj, int invert, WlzErrorNum *dstErr)
 Given a conforming mesh transform this function computes the strain tensor for each of it's valid elements. This function uses WlzCMeshDGTensor() to compute the displacement gradient tensor and the derives the strain tensor $e_{ij}$ from this using:

\[ e_{ij} = \frac{1}{2} (u_{ij} + u_{ji}) \]

.

WlzObjectWlzCMeshStrainTensorAtPts (WlzObject *cObj, int invert, WlzDVertex3 sd, int dither, WlzErrorNum *dstErr)
 Given a conforming mesh transform this function computes the displacement gradient tensor at regular cartesian grid sample points throughout the mesh. The tensor values at the sample points are computed using WlzCMeshDGTensorAtPts(). The strain tensor is then computed from the displacement gradient tensor as for WlzCMeshStrainTensor().
WlzVertexP WlzVerticesFromObj (WlzObject *obj, WlzVertexP *dstNr, int *dstCnt, WlzVertexType *dstType, WlzErrorNum *dstErr)
 Allocates a buffer which it fills with the vertices from the given object. The object must be one of the types that is represented by vertices, eg polylines, boundlists and contours.
WlzVertexP WlzVerticesFromObjBnd (WlzObject *obj, int *dstCnt, WlzVertexType *dstType, WlzErrorNum *dstErr)
 Extracts all vertices that lie on the boundary of the given objects domain.
WlzErrorNum WlzVerticesFromObjBnd2I (WlzObject *obj, int *dstNVtx, WlzIVertex2 **dstVtx)
 Extracts all vertices that lie on the boundary of the given 2D domain object's domain.
WlzErrorNum WlzVerticesFromObjBnd3I (WlzObject *obj, int *dstNVtx, WlzIVertex3 **dstVtx)
 Extracts all vertices that lie on the boundary of the given 3D domain object's domain.
WlzErrorNum WlzVerticesFromObj2I (WlzObject *obj, int *dstNVtx, WlzIVertex2 **dstVtx)
 Extracts all vertices that lie within the given 2D domain object's domain.
WlzErrorNum WlzVerticesFromObj3I (WlzObject *obj, int *dstNVtx, WlzIVertex3 **dstVtx)
 Extracts all vertices that lie within the given 3D domain object's domain.
WlzVertexP WlzVerticesFromGM (WlzGMModel *model, WlzVertexP *dstNr, int **dstVId, int *dstCnt, WlzVertexType *dstType, WlzErrorNum *dstErr)
 Allocates a buffer which it fills with the vertices from a geometric model.
WlzVertexP WlzDVerticesFromGM (WlzGMModel *model, int *dstCnt, WlzVertexType *dstType, WlzErrorNum *dstErr)
 Allocates a buffer which it fills with either 2D or 3D double precission vertices from the geometric model. The indicies of the vertices in the buffer are the same as the indices of the vertices in the model.
WlzVertexP WlzDVerticesFromCMesh (WlzCMeshP mesh, int *dstCnt, WlzVertexType *dstType, int skip, WlzErrorNum *dstErr)
 Allocates a buffer which it fills with either 2D or 3D double precission vertices from the conforming mesh. The indicies of the vertices in the buffer are the same as the indices of the vertices in the model unless the skip non-valid nodes flag is set..
AlcKDTTreeWlzVerticesBuildTree (WlzVertexType vType, int nV, WlzVertexP vtx, int *shfBuf, WlzErrorNum *dstErr)
 Allocates and populates a k-D tree from the given vertices. The vertices are either WlzDVertex2 orWlzDVertex3.

Enumeration Type Documentation

Scalar features of objects.

Enumerator:
WLZ_SCALARFEATURE_VALUE 

Grey value.

WLZ_SCALARFEATURE_GRADIENT 

Gradient of grey values.


Function Documentation

WlzLong WlzArea ( WlzObject obj,
WlzErrorNum dstErr 
)

Computes the area of an object.

Returns:
Area of object or a negative value on error.
Todo:
Implement for objects other than 2D domain objects.
Parameters:
objInput object.
dstErrDestination error pointer, may be NULL.

References _WlzIntervalWSpace::colrmn, _WlzDomain::core, _WlzObject::domain, _WlzObject::type, WLZ_2D_DOMAINOBJ, WLZ_EMPTY_OBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_EOO, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_RASTERDIR_ILIC, WlzInitRasterScan(), and WlzNextInterval().

Referenced by Wlz3DSectionOcc(), Wlz3DViewTransformObj(), WlzCCorS2D(), WlzMakeIntervalValues(), WlzSplitMontageObj(), WlzSplitObj(), WlzVerticesFromObj2I(), and WlzVolume().

WlzIBox2 WlzBoundingBox2I ( WlzObject inObj,
WlzErrorNum dstErr 
)

Computes the 2D integer axis aligned bounding box of the given object.

Returns:
2D integer bounding box.
Parameters:
inObjThe given object.
dstErrDestination error pointer, may be NULL.

References WLZ_ERR_NONE, WlzBoundingBox3I(), _WlzIBox3::xMax, _WlzIBox2::xMax, _WlzIBox3::xMin, _WlzIBox2::xMin, _WlzIBox3::yMax, _WlzIBox2::yMax, _WlzIBox3::yMin, and _WlzIBox2::yMin.

Referenced by WlzMeshFromObjBox(), WlzSplitMontageObj(), and WlzSplitObj().

WlzDBox2 WlzBoundingBox2D ( WlzObject inObj,
WlzErrorNum dstErr 
)

Computes the 2D double precision axis aligned bounding box of the given object.

Returns:
2D double precision bounding box.
Parameters:
inObjThe given object.
dstErrDestination error pointer, may be NULL.

References WLZ_ERR_NONE, WlzBoundingBox3D(), _WlzDBox3::xMax, _WlzDBox2::xMax, _WlzDBox3::xMin, _WlzDBox2::xMin, _WlzDBox3::yMax, _WlzDBox2::yMax, _WlzDBox3::yMin, and _WlzDBox2::yMin.

WlzIBox3 WlzBoundingBox3I ( WlzObject inObj,
WlzErrorNum dstErr 
)

Computes the integer 3D axis aligned bounding box of the given object.

Returns:
3D integer bounding box.
Parameters:
inObjThe given object.
dstErrDestination error pointer, may be NULL.

References _WlzDomain::b, _WlzValues::core, _WlzDomain::core, _WlzDomain::ctr, _WlzObject::domain, _WlzDomain::i, _WlzPlaneDomain::kol1, _WlzIntervalDomain::kol1, _WlzPlaneDomain::lastkl, _WlzIntervalDomain::lastkl, _WlzPlaneDomain::lastln, _WlzIntervalDomain::lastln, _WlzPlaneDomain::lastpl, _WlzPlaneDomain::line1, _WlzIntervalDomain::line1, _WlzValues::obj, _WlzDomain::p, _WlzPlaneDomain::plane1, _WlzDomain::poly, _WlzDomain::t, _WlzCoreDomain::type, _WlzObject::type, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_2D_POLYGON, WLZ_3D_DOMAINOBJ, WLZ_3D_POLYGON, WLZ_3D_WARP_TRANS, WLZ_AFFINE_TRANS, WLZ_BOUNDLIST, WLZ_COMPOUND_ARR_1, WLZ_COMPOUND_ARR_2, WLZ_COMPOUND_LIST_1, WLZ_COMPOUND_LIST_2, WLZ_CONTOUR, WLZ_CONV_HULL, WLZ_CONVOLVE_FLOAT, WLZ_CONVOLVE_INT, WLZ_EMPTY_OBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_VALUES_NULL, WLZ_FMATCHOBJ, WLZ_HISTOGRAM, WLZ_INTERVALDOMAIN_INTVL, WLZ_INTERVALDOMAIN_RECT, WLZ_PLANEDOMAIN_DOMAIN, WLZ_PROPERTY_OBJ, WLZ_RECTANGLE, WLZ_TEXT, WLZ_TRANS_OBJ, WLZ_TRANSFORM_2D_AFFINE, WLZ_WARP_TRANS, WlzBoundingBox3DTo3I(), _WlzIBox3::xMax, _WlzIBox3::xMin, _WlzIBox3::yMax, _WlzIBox3::yMin, _WlzIBox3::zMax, and _WlzIBox3::zMin.

Referenced by WlzBoundingBox2I(), WlzGreyValueMixing_s(), WlzRasterObj(), WlzSplitMontageObj(), and WlzSplitObj().

WlzDBox3 WlzBoundingBox3D ( WlzObject inObj,
WlzErrorNum dstErr 
)

Computes the double precision 3D axis aligned bounding box of the given object.

Returns:
3D double precision bounding box.
Parameters:
inObjThe given object.
dstErrDestination error pointer, may be NULL.

References _WlzDomain::b, _WlzValues::core, _WlzDomain::core, _WlzDomain::ctr, _WlzObject::domain, _WlzDomain::i, _WlzPlaneDomain::kol1, _WlzIntervalDomain::kol1, _WlzPlaneDomain::lastkl, _WlzIntervalDomain::lastkl, _WlzPlaneDomain::lastln, _WlzIntervalDomain::lastln, _WlzPlaneDomain::lastpl, _WlzPlaneDomain::line1, _WlzIntervalDomain::line1, _WlzValues::obj, _WlzDomain::p, _WlzPlaneDomain::plane1, _WlzDomain::poly, _WlzDomain::t, _WlzCompoundArray::type, _WlzCoreDomain::type, _WlzObject::type, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_2D_POLYGON, WLZ_3D_DOMAINOBJ, WLZ_3D_POLYGON, WLZ_3D_WARP_TRANS, WLZ_AFFINE_TRANS, WLZ_BOUNDLIST, WLZ_CMESH_2D, WLZ_CMESH_2D5, WLZ_CMESH_3D, WLZ_COMPOUND_ARR_1, WLZ_COMPOUND_ARR_2, WLZ_COMPOUND_LIST_1, WLZ_COMPOUND_LIST_2, WLZ_CONTOUR, WLZ_CONV_HULL, WLZ_CONVOLVE_FLOAT, WLZ_CONVOLVE_INT, WLZ_EMPTY_OBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_VALUES_NULL, WLZ_FMATCHOBJ, WLZ_HISTOGRAM, WLZ_INTERVALDOMAIN_INTVL, WLZ_INTERVALDOMAIN_RECT, WLZ_PLANEDOMAIN_DOMAIN, WLZ_PROPERTY_OBJ, WLZ_RECTANGLE, WLZ_TEXT, WLZ_TRANS_OBJ, WLZ_TRANSFORM_2D_AFFINE, WLZ_WARP_TRANS, WlzCMeshUpdateBBox2D(), WlzCMeshUpdateBBox2D5(), WlzCMeshUpdateBBox3D(), _WlzDBox3::xMax, _WlzDBox3::xMin, _WlzDBox3::yMax, _WlzDBox3::yMin, _WlzDBox3::zMax, and _WlzDBox3::zMin.

Referenced by WlzBoundingBox2D().

WlzDBox3 WlzBoundingBoxGModel3D ( WlzGMModel model,
WlzErrorNum dstErr 
)

Computes the 3D axis aligned bounding box of the 3D geometric model. It is an error if the model is not a 3D model.

Returns:
3D bounding box.
Parameters:
modelGiven model.
dstErrDestination error pointer, may be NULL.

References _WlzGMModel::child, _WlzGMShell::idx, _WlzGMShell::next, _WlzGMModel::type, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_GMMOD_3D, WLZ_GMMOD_3I, WLZ_GMMOD_3N, WlzGMShellGetGBB3D(), _WlzDBox3::xMax, _WlzDBox3::xMin, _WlzDBox3::yMax, _WlzDBox3::yMin, _WlzDBox3::zMax, and _WlzDBox3::zMin.

Referenced by WlzGeoModelGridWSpSet3D().

WlzDBox2 WlzBoundingBoxGModel2D ( WlzGMModel model,
WlzErrorNum dstErr 
)

Computes the 2D axis aligned bounding box of the 2D geometric model. It is an error if the model is not a 2D model.

Returns:
2D bounding box.
Parameters:
modelGiven model.
dstErrDestination error pointer, may be NULL.

References _WlzGMModel::child, _WlzGMShell::idx, _WlzGMShell::next, _WlzGMModel::type, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_GMMOD_2D, WLZ_GMMOD_2I, WLZ_GMMOD_2N, WlzGMShellGetGBB2D(), _WlzDBox2::xMax, _WlzDBox2::xMin, _WlzDBox2::yMax, and _WlzDBox2::yMin.

WlzIBox2 WlzBoundingBox2DTo2I ( WlzDBox2  bBox2D)

Converts a 2D double precision bounding box to an integer bounding box. There is no checking done for underflows or overflows.

Returns:
Integer bounding box.
Parameters:
bBox2DDouble precision bounding box.

References _WlzDBox2::xMax, _WlzIBox2::xMax, _WlzDBox2::xMin, _WlzIBox2::xMin, _WlzDBox2::yMax, _WlzIBox2::yMax, _WlzDBox2::yMin, and _WlzIBox2::yMin.

WlzDBox2 WlzBoundingBox2ITo2D ( WlzIBox2  bBox2I)

Converts a 2D integer bounding box to an double precision bounding box.

Returns:
Double precision bounding box.
Parameters:
bBox2IDouble precision bounding box.

References _WlzIBox2::xMax, _WlzDBox2::xMax, _WlzIBox2::xMin, _WlzDBox2::xMin, _WlzIBox2::yMax, _WlzDBox2::yMax, _WlzIBox2::yMin, and _WlzDBox2::yMin.

WlzIBox3 WlzBoundingBox3DTo3I ( WlzDBox3  bBox3D)

Converts a 3D double precision bounding box to an integer bounding box. There is no checking done for underflows or overflows.

Returns:
Integer bounding box.
Parameters:
bBox3DDouble precision bounding box.

References _WlzDBox3::xMax, _WlzIBox3::xMax, _WlzDBox3::xMin, _WlzIBox3::xMin, _WlzDBox3::yMax, _WlzIBox3::yMax, _WlzDBox3::yMin, _WlzIBox3::yMin, _WlzDBox3::zMax, _WlzIBox3::zMax, _WlzDBox3::zMin, and _WlzIBox3::zMin.

Referenced by WlzBoundingBox3I().

WlzDBox3 WlzBoundingBox3ITo3D ( WlzIBox3  bBox3I)

Converts a 3D integer bounding box to an double precision bounding box.

Returns:
Double precision bounding box.
Parameters:
bBox3IInteger bounding box.

References _WlzIBox3::xMax, _WlzDBox3::xMax, _WlzIBox3::xMin, _WlzDBox3::xMin, _WlzIBox3::yMax, _WlzDBox3::yMax, _WlzIBox3::yMin, _WlzDBox3::yMin, _WlzIBox3::zMax, _WlzDBox3::zMax, _WlzIBox3::zMin, and _WlzDBox3::zMin.

WlzFBox3 WlzBoundingBox3DTo3F ( WlzDBox3  bBox3D)

Converts a 3D double bounding box to a single precision bounding box.

Returns:
Single precision bounding box.
Parameters:
bBox3DDouble precision bounding box.

References _WlzDBox3::xMax, _WlzFBox3::xMax, _WlzDBox3::xMin, _WlzFBox3::xMin, _WlzDBox3::yMax, _WlzFBox3::yMax, _WlzDBox3::yMin, _WlzFBox3::yMin, _WlzDBox3::zMax, _WlzFBox3::zMax, _WlzDBox3::zMin, and _WlzFBox3::zMin.

WlzDBox3 WlzBoundingBox3FTo3D ( WlzFBox3  bBox3F)

Converts a 3D single bounding box to a double precision bounding box.

Returns:
Double precision bounding box.
Parameters:
bBox3FSingle precision bounding box.

References _WlzFBox3::xMax, _WlzDBox3::xMax, _WlzFBox3::xMin, _WlzDBox3::xMin, _WlzFBox3::yMax, _WlzDBox3::yMax, _WlzFBox3::yMin, _WlzDBox3::yMin, _WlzFBox3::zMax, _WlzDBox3::zMax, _WlzFBox3::zMin, and _WlzDBox3::zMin.

WlzIBox2 WlzBoundingBoxUnion2I ( WlzIBox2  box0,
WlzIBox2  box1 
)

Computes the 2D integer bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 2D integer bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzIBox2::xMax, _WlzIBox2::xMin, _WlzIBox2::yMax, and _WlzIBox2::yMin.

WlzFBox2 WlzBoundingBoxUnion2F ( WlzFBox2  box0,
WlzFBox2  box1 
)

Computes the 2D single precision bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 2D single precision bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzFBox2::xMax, _WlzFBox2::xMin, _WlzFBox2::yMax, and _WlzFBox2::yMin.

WlzDBox2 WlzBoundingBoxUnion2D ( WlzDBox2  box0,
WlzDBox2  box1 
)

Computes the 2D double precision bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 2D double precision bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzDBox2::xMax, _WlzDBox2::xMin, _WlzDBox2::yMax, and _WlzDBox2::yMin.

WlzIBox3 WlzBoundingBoxUnion3I ( WlzIBox3  box0,
WlzIBox3  box1 
)

Computes the 3D integer bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 3D integer bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzIBox3::xMax, _WlzIBox3::xMin, _WlzIBox3::yMax, _WlzIBox3::yMin, _WlzIBox3::zMax, and _WlzIBox3::zMin.

WlzFBox3 WlzBoundingBoxUnion3F ( WlzFBox3  box0,
WlzFBox3  box1 
)

Computes the 3D single precision bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 3D float precision bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzFBox3::xMax, _WlzFBox3::xMin, _WlzFBox3::yMax, _WlzFBox3::yMin, _WlzFBox3::zMax, and _WlzFBox3::zMin.

WlzDBox3 WlzBoundingBoxUnion3D ( WlzDBox3  box0,
WlzDBox3  box1 
)

Computes the 3D double precision bounding box which encloses the given pair of bounding boxes.

Returns:
Axis alligned 3D double precision bounding box.
Parameters:
box0First bounding box.
box1Second bounding box.

References _WlzDBox3::xMax, _WlzDBox3::xMin, _WlzDBox3::yMax, _WlzDBox3::yMin, _WlzDBox3::zMax, and _WlzDBox3::zMin.

double WlzCCorS2D ( WlzObject obj0,
WlzObject obj1,
int  unionFlg,
int  normFlg,
WlzErrorNum dstErr 
)

Computes the cross correlation of the two given 2D spatial domain objects in the spatial domain. The algorithm used by this function is not an efficient unless a single cross correlation value is required.

Returns:
Cross correlation value.
Parameters:
obj0First object. Must have been assigned.
obj1Second object. Must have been assigned.
unionFlgComputes the cross correlation value in the union of the two objects domains if non zero. The default is to use the intersection of the domains.
normFlgNormalise the cross-correlation value by dividing it by the area/volume over which it is computed if this flag is non-zero.
dstErrDestination error pointer, may be NULL.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzGreyValueWSpace::gType, _WlzGreyValueWSpace::gVal, _WlzIntervalWSpace::lftpos, _WlzIntervalWSpace::linpos, _WlzIntervalWSpace::rgtpos, _WlzObject::type, _WlzObject::values, _WlzIVertex2::vtX, _WlzIVertex2::vtY, WLZ_2D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_EOO, WLZ_ERR_GREY_TYPE, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_VALUES_NULL, WLZ_GREY_DOUBLE, WLZ_GREY_FLOAT, WLZ_GREY_INT, WLZ_GREY_SHORT, WLZ_GREY_UBYTE, WLZ_RASTERDIR_ILIC, WlzArea(), WlzFreeObj(), WlzGreyValueFreeWSp(), WlzGreyValueGet(), WlzGreyValueMakeWSp(), WlzInitRasterScan(), WlzIntersect2(), WlzNextInterval(), and WlzUnion2().

double WlzCentrality ( WlzObject fObj,
WlzObject bObj,
int  nRay,
int  binFlg,
double *  dstMaxR,
WlzErrorNum dstErr 
)

Computes the centrality of a feature domain with respect to a boundary domain. The boundary domain must enclose the feature domain.

Returns:
Centrality feature value in range [0.0-1.0].

Rays are projected from the centre of mass of the boundary domain to it's boundary. The maximum value of this distance for a given angle $i$ being $R_i$. At each point in the feature domain intersected by the ray the currrent value of the centrality is updated. The distance from the centre of mass to each feature point is defined to be $r_{i,j}$ and the mass at this point $m_{i,j}$. The centrality of the feature domain $\Omega_f$ with respect to the border domain $Omega_b$ is defined to be

\[ c = \frac{\sum_{i,j}{m_{i,j}(R_i - r_{i,j})}} {\sum_{i,j}{m_{i,j}R}} \]

Parameters:
fObjFeature domain object, $Omega_f$.
bObjBoundary domain object, $Omega_b$.
nRayNumber of equally spaced rays projected from the centre of mass.
binFlgTreat as binary object if non-zero, with all masses having value 1.0.
dstMaxRDestination pointer for maximum boundary radius, may be NULL.
dstErrDestination error pointer, may be NULL.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzObject::type, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_UNIMPLEMENTED, and WLZ_ERR_VALUES_NULL.

WlzDVertex2 WlzCentreOfMass2D ( WlzObject srcObj,
int  binObjFlag,
double *  dstMass,
WlzErrorNum dstErr 
)

Calculates the centre of mass of a Woolz object. If the given object does not have grey values or the binary object flag is set (ie non zero) then every pixel or vertex within the object's domain has the same mass.

Returns:
Coordinates of center of mass.
Parameters:
srcObjGiven object.
binObjFlagBinary object flag.
dstMassDestination pointer for mass, may be NULL if not required.
dstErrDestination pointer for error, may be NULL.

References _WlzObject::domain, _WlzValues::obj, _WlzDomain::t, _WlzObject::type, _WlzObject::values, _WlzDVertex2::vtX, _WlzDVertex2::vtY, WLZ_2D_DOMAINOBJ, WLZ_CONTOUR, WLZ_DBG, WLZ_DBG_LVL_1, WLZ_DBG_LVL_2, WLZ_DBG_LVL_FN, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, and WLZ_TRANS_OBJ.

WlzDVertex3 WlzCentreOfMass3D ( WlzObject srcObj,
int  binObjFlag,
double *  dstMass,
WlzErrorNum dstErr 
)

Calculates the centre of mass of a Woolz object. If the given object does not have grey values or the binary object flag is set (ie non zero) then every pixel or vertex within the objects domain has the same mass.

Returns:
Coordinates of center of mass.
Parameters:
srcObjGiven object.
binObjFlagBinary object flag.
dstMassDestination pointer for mass, may be NULL if not required.
dstErrDestination pointer for error, may be NULL.

References _WlzObject::domain, _WlzValues::obj, _WlzDomain::t, _WlzObject::type, _WlzObject::values, _WlzDVertex2::vtX, _WlzDVertex3::vtX, _WlzDVertex2::vtY, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_CONTOUR, WLZ_DBG, WLZ_DBG_LVL_1, WLZ_DBG_LVL_2, WLZ_DBG_LVL_FN, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, and WLZ_TRANS_OBJ.

WlzDVertex2 WlzCentreOfMassVtx2D ( int  nVtx,
WlzDVertex2 vtx 
)

Computes the centre of mass of a vector of 2D vertices.

Returns:
Coordinates of centre of mass.
Parameters:
nVtxNumber of vertices.
vtxThe vertices.

References WLZ_VTX_2_ADD, WLZ_VTX_2_SCALE, and WLZ_VTX_2_ZERO.

WlzDVertex3 WlzCentreOfMassVtx3D ( int  nVtx,
WlzDVertex3 vtx 
)

Computes the centre of mass of a vector of 3D vertices.

Returns:
Coordinates of centre of mass.
Parameters:
nVtxNumber of vertices.
vtxThe vertices.

References WLZ_VTX_3_ADD, WLZ_VTX_3_SCALE, and WLZ_VTX_3_ZERO.

Referenced by WlzGeometryLSqOPlane().

WlzErrorNum WlzDistMetricGM ( WlzGMModel model0,
WlzGMModel model1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the vertices of the given geometric models. See WlzDistMetricVertex2D() for details of the metrics.

Returns:
Woolz error code.
Parameters:
model0First geometric model.
model1Second geometric model.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References AlcFree(), _WlzGMModel::type, _WlzVertexP::v, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_VERTEX_D2, WLZ_VERTEX_D3, WlzDistMetricVertex2D(), WlzDistMetricVertex3D(), and WlzVerticesFromGM().

WlzErrorNum WlzDistMetricDirGM ( WlzGMModel model0,
WlzGMModel model1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the vertices of the given geometric models. See WlzDistMetricDirVertex2D() for details of the metrics.

Returns:
Woolz error code.
Parameters:
model0First geometric model.
model1Second geometric model.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References AlcFree(), _WlzGMModel::type, _WlzVertexP::v, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_VERTEX_D2, WLZ_VERTEX_D3, WlzDistMetricDirVertex2D(), WlzDistMetricDirVertex3D(), and WlzVerticesFromGM().

WlzErrorNum WlzDistMetricVertex2D ( int  n0,
WlzDVertex2 vx0,
int  n1,
WlzDVertex2 vx1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. Each of these distance measures is the maximum of the two possible directed measures:

\[ D = \max{(d(A,B), d(B,A))} \]

Where $D$ is a non-directed distance metric and $d$ is the associated directed distance metric. $A$ and $B$ are the two datasets.

Returns:
Woolz error code.
Parameters:
n0Number of vertices in the first array.
vx0First array of vertices.
n1Number of vertices in the second array.
vx1Second array of vertices.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References WLZ_ERR_NONE, WLZ_MAX, and WlzDistMetricDirVertex2D().

Referenced by WlzDistMetricGM().

WlzErrorNum WlzDistMetricVertex3D ( int  n0,
WlzDVertex3 vx0,
int  n1,
WlzDVertex3 vx1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. See WlzDistMetricVertex2D() for an explaination of the distance metrics.

Returns:
Woolz error code.
Parameters:
n0Number of vertices in the first array.
vx0First array of vertices.
n1Number of vertices in the second array.
vx1Second array of vertices.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References WLZ_ERR_NONE, WLZ_MAX, and WlzDistMetricDirVertex3D().

Referenced by WlzDistMetricGM().

WlzErrorNum WlzDistMetricDirVertex2D ( int  n0,
WlzDVertex2 vx0,
int  n1,
WlzDVertex2 vx1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum mean neighbour distances between the given sets of vertices. The directed Hausdorff distance metric is:

\[ \max_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

The directed mean nearest neighbour distance metric is:

\[ \frac{1}{n_a} \sum_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

Likewise the directed median nearest neighbour distance metric is:

\[ \mathrm{median}_{a \in A}{\min_{b \in B}{\|a - b\|}} \]

Where $A$ and $B$ are the two datasets.

Returns:
Woolz error code.
Parameters:
n0Number of vertices in the first array.
vx0First array of vertices.
n1Number of vertices in the second array.
vx1Second array of vertices.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References AlcFree(), AlcKDTGetNN(), AlcKDTTreeFree(), AlcMalloc(), AlgRankSelectD(), _WlzVertexP::d2, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WLZ_VERTEX_D2, and WlzVerticesBuildTree().

Referenced by WlzDistMetricDirGM(), and WlzDistMetricVertex2D().

WlzErrorNum WlzDistMetricDirVertex3D ( int  n0,
WlzDVertex3 vx0,
int  n1,
WlzDVertex3 vx1,
double *  dstDistH,
double *  dstDistM,
double *  dstDistN,
double *  dstDistI 
)

Computes any combination of the directed Hausdorff, mean nearest neighbour, median nearest neighbour and minimum nearest neighbour distances between the given sets of vertices. See WlzDistMetricDirVertex2D() for an explaination of the distance metrics.

Returns:
Woolz error code.
Parameters:
n0Number of vertices in the first array.
vx0First array of vertices.
n1Number of vertices in the second array.
vx1Second array of vertices.
dstDistHDestination pointer for the directed Hausdorff distance, may be NULL.
dstDistMDestination pointer for the directed mean nearest neighbour distance, may be NULL.
dstDistNDestination pointer for the directed median nearest neighbour distance, may be NULL.
dstDistIDestination pointer for the minimum nearest neighbour distance, may be NULL.

References AlcFree(), AlcKDTGetNN(), AlcKDTTreeFree(), AlcMalloc(), AlgRankSelectD(), _WlzVertexP::d3, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WLZ_VERTEX_D3, and WlzVerticesBuildTree().

Referenced by WlzDistMetricDirGM(), and WlzDistMetricVertex3D().

WlzDVertex3* WlzGeometryTrackUpAndDown_s ( int  numberOfPixelsZ,
int  startTrackingFile,
int  numberOfFilesDownOrUp,
double  disForInAndOutGuid,
double  disForInAndOut,
unsigned char **  TwoDImageFilesNameList,
int  numOf2DWlzFiles,
int  downOrUp,
int  sectionLength_N,
int  subSubSectionLength_L,
int  numberOfSampleP_k,
char *  surfacePointFileName,
char *  surfaceInPointFileName,
char *  surfaceOutPointFileName,
int  startShell,
int  endShell,
int  startSection,
int  endSection,
double  minDis,
WlzErrorNum dstErr 
)

Track a curved path through a set of geometric model shells.

Returns:
An array of vertices tracking thougha set of geometry model shells
Parameters:
numberOfPixelsZthe input z-derection of pixel number (input)
startTrackingFilenot documented
numberOfFilesDownOrUpnumner of files (?)
disForInAndOutGuiddistance parameter
disForInAndOutdistance parameter
TwoDImageFilesNameList2D image files name list
numOf2DWlzFilesnumber of 2D woolz files
downOrUptracking direction flag
sectionLength_Nlength used to cut a patch (unit in pixel)
subSubSectionLength_Llength, smaller than sectionLength_N used to sample points
numberOfSampleP_knumber of points will be sampled in the above section
surfacePointFileNameFileNameStr to output the surface points
surfaceInPointFileNameFileNameStr to output the in surface points
surfaceOutPointFileNameFileNameStr to output the out surface points
startShellthe n-th Shell to begin with tracking
endShellthe end Shell from where stop tracking
startSectionthe n-th Section to begin with tracking
endSectionthe section to stop tracking.
minDisminimum distance (?)
dstErrerror return.
Source:
WlzGeometryTrackUpAndDown_s.c
Author:
: J. Rao, R. Baldock and B. Hill

References ALC_ER_NONE, AlcFree(), AlcKDTGetNN(), AlcKDTTreeFacts(), AlcMalloc(), _WlzDomain::ctr, _WlzVertexP::d2, _WlzObject::domain, _AlcKDTNode::idx, _AlcPointP::kD, _AlcKDTNode::key, MaxNumOfFiles, _WlzContour::model, NumberToTrack, _WlzDVertex2::vtX, _WlzDVertex3::vtX, _WlzDVertex2::vtY, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_CONTOUR, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WlzFreeObj(), and WlzReadObj().

size_t WlzGreySize ( WlzGreyType  gType)

Given a grey type, computes the number of bytes required to to store a single grey value of that type. Note the grey type WLZ_GREY_BIT is able to store multiple values in a single byte but one byte is required to store a single bit.

Returns:
Number of bytes required to store a single grey value or zero on error (invalid grey type).
Parameters:
gTypeGiven grey type.

References WLZ_GREY_BIT, WLZ_GREY_DOUBLE, WLZ_GREY_FLOAT, WLZ_GREY_INT, WLZ_GREY_LONG, WLZ_GREY_RGBA, WLZ_GREY_SHORT, and WLZ_GREY_UBYTE.

Referenced by WlzBuildObj3B(), WlzCutObjToValBox3D(), WlzFreeTiledValues(), WlzIndexedValueSize(), WlzMakeCuboid(), WlzMakeIndexedValues(), WlzMakeLUTValues(), WlzMakePointValues(), and WlzMakeTiledValuesTiles().

int WlzGreyStats ( WlzObject srcObj,
WlzGreyType dstGType,
double *  dstMin,
double *  dstMax,
double *  dstSum,
double *  dstSumSq,
double *  dstMean,
double *  dstStdDev,
WlzErrorNum dstErr 
)

Calculates simple quick statistics for given 2D or 3D domain object with grey values. Pointers provided for results may be NULL without causing an error.

Returns:
Object area or -1 on error.
Parameters:
srcObjObject from which to calculate the statistics.
dstGTypePointer for grey type.
dstMinPointer for minimum value.
dstMaxPointer for maximum value.
dstSumPointer for sum of values.
dstSumSqPointer for sum of squares of values.
dstMeanMean value.
dstStdDevStandard deviation of values.
dstErrDestination pointer for error number, may be NULL if not required.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzPlaneDomain::domains, _WlzPlaneDomain::lastpl, _WlzDomain::p, _WlzPlaneDomain::plane1, _WlzCoreDomain::type, _WlzObject::type, _WlzVoxelValues::values, _WlzObject::values, _WlzValues::vox, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_DBG, WLZ_DBG_LVL_1, WLZ_DBG_LVL_2, WLZ_DBG_LVL_FN, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_VALUES_DATA, WLZ_ERR_VALUES_NULL, WLZ_PLANEDOMAIN_DOMAIN, WlzFreeObj(), and WlzMakeMain().

Referenced by main(), WlzEffWriteObjAnl(), WlzGetPatchTreeToDepth(), WlzPatchTreeToObject(), and WlzRGBAGreyStats().

WlzObject* WlzGreyValueMixing_s ( WlzObject sObj,
WlzObject tObj,
double  xmiddle,
WlzErrorNum dstErr 
)

calculate the distance map of a woolz obj

Returns:
Woolz obj contain the mixing grey value.
Parameters:
sObjGiven Woolz grey value source object.
tObjGiven Woolz grey value target object.
xmiddlemixing parameters
dstErrDestination error pointer.

References _WlzObject::type, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WlzBoundingBox3I(), WlzCopyObject(), WlzGreyValueFreeWSp(), WlzGreyValueMakeWSp(), _WlzIBox3::xMax, _WlzIBox3::xMin, _WlzIBox3::yMax, _WlzIBox3::yMin, _WlzIBox3::zMax, and _WlzIBox3::zMin.

int WlzLineArea ( WlzObject obj,
WlzErrorNum dstErr 
)

Calculate the line-area of an object defined as the sum of the line segments bounded by the left hand end of the first interval in a line and the right hand end of the last interval in that line.

Returns:
Line area calculated, -1 on error.
Parameters:
objInput object.
dstErrError return.
Source:
WlzLineArea.c

References _WlzDomain::core, _WlzObject::domain, _WlzDomain::i, _WlzIntervalWSpace::intrmn, _WlzIntervalDomain::kol1, _WlzIntervalDomain::lastkl, _WlzIntervalDomain::lastln, _WlzIntervalWSpace::lftpos, _WlzIntervalDomain::line1, _WlzIntervalWSpace::rgtpos, _WlzIntervalDomain::type, _WlzObject::type, WLZ_2D_DOMAINOBJ, WLZ_EMPTY_OBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_EOO, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_INTERVALDOMAIN_INTVL, WLZ_INTERVALDOMAIN_RECT, WLZ_RASTERDIR_ILIC, WlzInitRasterScan(), and WlzNextInterval().

Referenced by WlzNewValueTb().

WlzObject* WlzMakeMarkers ( WlzVertexType  vType,
int  nVtx,
WlzVertexP  vtx,
WlzMarkerType  mType,
int  mSz,
WlzErrorNum dstErr 
)

Constructs a domain from the union of marker domains with a marker domain at each of the given vertex positions.

Returns:
New domain object without values.
Parameters:
nVtxNumber of vertices.
vtxGiven vertices.
mTypeMarker type.
mSzMarker size. This is the radius of a sphere marker. The marker size is ignored for point markers.
dstErrDestination error pointer, may be NULL.

References _WlzVertexP::i2, _WlzVertexP::i3, _WlzVertexP::v, _WlzIVertex3::vtX, _WlzIVertex3::vtY, _WlzIVertex3::vtZ, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_ERR_NONE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WLZ_ERR_PARAM_TYPE, WLZ_MARKER_POINT, WLZ_MARKER_SPHERE, WLZ_VERTEX_I2, WLZ_VERTEX_I3, WlzFreeObj(), WlzMakeEmpty(), WlzMakeSinglePixelObject(), WlzMakeSphereObject(), WlzShiftObject(), and WlzUnion2().

Referenced by WlzMarkerLattice().

WlzObject* WlzMarkerLattice ( WlzObject gObj,
WlzMarkerType  mType,
int  mSz,
int  mSep,
WlzErrorNum dstErr 
)
WlzPoints* WlzMakePoints ( WlzObjectType  type,
int  nVtx,
WlzVertexP  vtxP,
int  maxVtx,
WlzErrorNum dstErr 
)

Creates a new point domain. A point domain consists of an array of vertices which are treated as seperate points.

Returns:
New point domain.
Parameters:
typeType of point vertices, which must be one of WLZ_POINTS_2I, WLZ_POINTS_2D, WLZ_POINTS_3I or WLZ_POINTS_3D.
nVtxNumber of points to copy to the new domain.
vtxPPoints to copy to the new domain. These must be of the correct type. If NULL then no points are copied.
maxVtxThe number of vertices for which space is allocated.
dstErrDestination error pointer, may be NULL.

References AlcCalloc(), _WlzPoints::maxPoints, _WlzPoints::nPoints, _WlzPoints::points, _WlzPoints::type, _WlzVertexP::v, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_POINTS_2D, WLZ_POINTS_2I, WLZ_POINTS_3D, and WLZ_POINTS_3I.

WlzErrorNum WlzNObjGreyStats ( WlzObject gObj,
int  mean,
int  stddev,
int *  dstN,
WlzObject **  dstMinObj,
WlzObject **  dstMaxObj,
WlzObject **  dstSumObj,
WlzObject **  dstSSqObj 
)

Computes a collection of objects from the object (which should be a compound array object). The computed objects all have the intersection of the given object domains as their domain and WLZ_GREY_DOUBLE values for the minimum, maximum, sum and sum of squares of the given grey values at each pixel/voxel within the intersection domain.

Returns:
Woolz error code.
Parameters:
gObjGiven object.
meanCompute mean not sum if non-zero.
stddevCompute standard deviation not sum of squares if non-zero.
dstNDestination pointer for number of objects, may be NULL.
dstMinObjDestination pointer for minimum value object, may be NULL if not required.
dstMaxObjDestination pointer for maximum value object, may be NULL if not required.
dstSumObjDestination pointer for sum of values object, may be NULL if not required.
dstSSqObjDestination pointer for sum of squares object, may be NULL if not required.

References _WlzValues::core, _WlzDomain::core, _WlzGreyV::dbv, _WlzObject::domain, _WlzPixelV::type, _WlzObject::type, _WlzPixelV::v, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_BO_SUBTRACT, WLZ_COMPOUND_ARR_1, WLZ_COMPOUND_ARR_2, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_VALUES_NULL, WLZ_FN_SCALAR_SQR, WLZ_FN_SCALAR_SQRT, WLZ_GREY_DOUBLE, WLZ_GREY_TAB_RAGR, WlzFreeObj(), WlzGreySetValue(), WlzGreyTableType(), WlzImageArithmetic(), WlzIntersectN(), WlzNewObjectValues(), WlzScalarFn(), and WlzScalarMultiply().

WlzErrorNum Wlz3DSectionOcc ( WlzObject obj,
WlzThreeDViewStruct vs,
double  sep,
double *  dstFirst,
double *  dstLast,
int *  dstArraySizeOcc,
int **  dstArrayOcc 
)

Computes an array of integers which correspond to a line through the given object's domain perpendicular to the plane defined by the given view structure. At each point along the line the occupancy area or number of domains intersecting the plane is computed. When the given object is a 3D spatial domain object then areas are computed and when it is a compound array then the occupancy (number of domains) is computed.

Returns:
Woolz error code.
Parameters:
objGiven object with the domain or a compound object with many domains.
vsGiven view structure which defines the cutting plane. The view structure will be initialised within this function.
sepPlane separation distance, must be greater than ALG_DBL_TOLLERANCE.
dstFirstDestination pointer for the first coordinate of the line in a 1D coordinate system perpendicular to the cutting plane. May be NULL.
dstLastDestination pointer for the last coordinate of the line in a 1D coordinate system perpendicular to the cutting plane. May be NULL.
dstArraySizeOccDestination pointer for the size of the occupancy array. May be NULL.
dstArrayOccDestination pointer for the occupancy array. The occupancy array should be freed using AlcFree(). May be NULL.

References AlcCalloc(), AlcFree(), ALG_DBL_TOLLERANCE, _WlzDomain::core, _WlzThreeDViewStruct::dist, _WlzObject::domain, _WlzThreeDViewStruct::maxvals, _WlzThreeDViewStruct::minvals, _WlzCompoundArray::n, _WlzCompoundArray::o, _WlzThreeDViewStruct::ref_obj, _WlzThreeDViewStruct::type, _WlzObject::type, _WlzDVertex3::vtZ, WLZ_3D_DOMAINOBJ, WLZ_3D_VIEW_STRUCT, WLZ_COMPOUND_ARR_1, WLZ_COMPOUND_ARR_2, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOUBLE_DATA, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WLZ_ERR_PARAM_TYPE, WLZ_INTERPOLATION_NEAREST, WlzArea(), WlzAssignObject(), WlzFree3DViewStruct(), WlzFreeObj(), WlzGetSubSectionFromObject(), WlzInit3DViewStruct(), WlzIsEmpty(), WlzMake3DViewStructCopy(), and WlzUnionN().

WlzObject* WlzDomainOccupancy ( WlzObject gObj,
WlzErrorNum dstErr 
)

Computes a new Woolz domain object which has the domain of the given object (union of domains if the object is a compound array) and a value table with integer values representing the number of domains present.

Returns:
New Woolz object or NULL on error.
Parameters:
gObjGiven object.
dstErrDestination error pointer.

References _WlzGreyV::inv, _WlzCompoundArray::n, _WlzCompoundArray::o, _WlzGreyV::shv, _WlzObject::type, _WlzPixelV::type, _WlzGreyV::ubv, _WlzPixelV::v, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_COMPOUND_ARR_1, WLZ_COMPOUND_ARR_2, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_GREY_INT, WLZ_GREY_SHORT, WLZ_GREY_TAB_RAGR, WLZ_GREY_UBYTE, WlzFreeObj(), WlzGreyIncValuesInDomain(), WlzGreySetValue(), WlzGreyTableType(), WlzNewObjectValues(), and WlzUnionN().

WlzObject* WlzPointsToDomObj ( WlzPoints pnt,
double  scale,
WlzErrorNum dstErr 
)

Creates a domain object which coresponds to the union of the given points.

Returns:
New domain object which coresponds to the union of the given points.
Parameters:
pntPoint domain.
scaleScale, which if greater than zero is used as the diameter of a circle or sphere centred on each of the points vertices and a multiplier for the point position.
dstErrDestination error poiter, may be NULL.

References _WlzVertexP::d2, _WlzVertexP::d3, _WlzVertex::d3, _WlzVertexP::i2, _WlzVertexP::i3, _WlzVertex::i3, _WlzPoints::points, _WlzPoints::type, _WlzIVertex3::vtX, _WlzDVertex2::vtX, _WlzIVertex2::vtX, _WlzDVertex3::vtX, _WlzIVertex3::vtY, _WlzDVertex2::vtY, _WlzIVertex2::vtY, _WlzDVertex3::vtY, _WlzIVertex3::vtZ, _WlzDVertex3::vtZ, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_ERR_PARAM_TYPE, WLZ_NINT, WLZ_POINTS_2D, WLZ_POINTS_2I, WLZ_POINTS_3D, WLZ_POINTS_3I, WlzFreeObj(), WlzMakeEmpty(), WlzMakeSinglePixelObject(), WlzMakeSphereObject(), and WlzUnion2().

Referenced by WlzDistanceTransform().

int WlzRGBAGreyStats ( WlzObject srcObj,
WlzRGBAColorSpace  colSpc,
WlzGreyType dstGType,
double *  dstMin,
double *  dstMax,
double *  dstSum,
double *  dstSumSq,
double *  dstMean,
double *  dstStdDev,
WlzErrorNum dstErr 
)

Calculates simple quick statistics for the domain object with RGBA values. Each component has its statictics computed and entered into the four double[4] arrays.

Returns:
Object area or -1 on error.
Parameters:
srcObjGiven object.
colSpcColour space.
dstGTypePointer for grey type.
dstMinArray for the 4 minimum value.
dstMaxArray for the 4 maximum value.
dstSumArray for the 4 sum of values.
dstSumSqArray for the 4 sum of squares of values.
dstMeanArray for the 4 mean values.
dstStdDevArray for the 4 standard deviation values.
dstErrDestination pointer for error, may be NULL.

References _WlzCompoundArray::o, WLZ_ERR_NONE, WLZ_GREY_RGBA, WlzFreeObj(), WlzGreyStats(), and WlzRGBAToCompound().

Referenced by WlzPatchTreeToObject().

WlzErrorNum WlzSampleValuesAndCoords ( WlzObject obj,
WlzGreyType dstGType,
int *  dstNVal,
WlzGreyP dstValP,
WlzVertexP dstCoords,
WlzSampleFn  samFn,
int  samFac 
)

Allocates buffers for both the grey values and the coordinates of the grey values in the given object. On return these buffers contain the sampled object values and the coordinates of the values.

Returns:
Parameters:
objGiven object.
dstGTypeType of grey value.
dstNValNumber of grey values, also the number of coordinates.
dstValPDestination pointer for the values.
dstCoordsDestination pointer for the coordinates, these are always either WlzIVertex2 or WlzIVertex3 depending on the objects dimension.
samFnSampling function.
samFacSampling factor.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzVertexP::i2, _WlzObject::type, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_PARAM_DATA, WLZ_ERR_PARAM_NULL, WLZ_ERR_VALUES_NULL, WLZ_SAMPLEFN_MAX, and WLZ_SAMPLEFN_MIN.

Referenced by WlzScalarFeatures2D().

WlzErrorNum WlzScalarFeatures2D ( WlzObject obj,
int *  dstNFeat,
WlzIVertex2 **  dstFeat,
WlzScalarFeatureType  fType,
WlzThresholdType  thrHL,
WlzPixelV  thrV,
double  filterV,
double  minDist 
)

Finds scalar features within the given object. Where the feature values are all either above or below the given threshold value, depending on the value of thrHL, and the features have the given minimum seperation distance.

Returns:
Woolz error code.

The features are found using the following algorithm:

  • Compute feature values from given object.
  • Threshold and extract a list of features.
  • Sort the list of features by value.
  • Create a kD-tree.
  • While the list of feature values is not empty.
    • Remove item from list of features.
    • Find distance to nearest neighbour in the kD-tree.
    • If the distance is greater than the minimum.
      • Add the feature to the kD-tree.
  • Output features of the kD-tree.

The cost of this function depends on the threshold value and to some extent on the minimum distance.

If feature values are computed then they are normalised using WlzGreyNormalise(), to allow threshold values to be set.

Parameters:
objGiven object.
fTypeType of feature to find.
thrHLHigh or low feature values.
thrVThreshold feature value.
minDistMinimum distance between features.
dstNFeatDestination pointer for the number of scalar features found.
dstFeatDestination pointer for coordinates of the features found.
filterVFilter value for computing features.
dstFeatDestination pointer for the coordinates of the scalar features.

References ALC_ER_NONE, ALC_POINTTYPE_DBL, AlcFree(), AlcKDTGetNN(), AlcKDTInsert(), AlcKDTTreeFree(), AlcKDTTreeNew(), AlcMalloc(), AlgHeapSortCmpIdxDFn(), AlgHeapSortCmpIdxFFn(), AlgHeapSortCmpIdxIFn(), AlgHeapSortCmpIdxLFn(), AlgHeapSortCmpIdxSFn(), AlgHeapSortCmpIdxUFn(), AlgHeapSortIdx(), AlgHeapSortInvCmpIdxDFn(), AlgHeapSortInvCmpIdxFFn(), AlgHeapSortInvCmpIdxIFn(), AlgHeapSortInvCmpIdxLFn(), AlgHeapSortInvCmpIdxSFn(), AlgHeapSortInvCmpIdxUFn(), _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzVertexP::i2, _WlzGreyP::inp, _WlzObject::type, _WlzVertexP::v, _WlzObject::values, WLZ_2D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_GREY_TYPE, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_PARAM_DATA, WLZ_ERR_VALUES_NULL, WLZ_GREY_DOUBLE, WLZ_GREY_FLOAT, WLZ_GREY_INT, WLZ_GREY_LONG, WLZ_GREY_SHORT, WLZ_GREY_UBYTE, WLZ_SAMPLEFN_MAX, WLZ_SAMPLEFN_MIN, WLZ_SCALARFEATURE_GRADIENT, WLZ_SCALARFEATURE_VALUE, WLZ_THRESH_HIGH, WLZ_THRESH_LOW, WlzAssignObject(), WlzFreeObj(), WlzGreyModGradient(), WlzGreyNormalise(), WlzSampleValuesAndCoords(), and WlzThreshold().

WlzObject* WlzCMeshDGTensor ( WlzObject cObj,
int  invert,
WlzErrorNum dstErr 
)

Given a conforming mesh transform this function computes the displacement gradient tensor for each of it's valid elements. Given displacement $\vec{u}(\vec{r})$ with position vector $\vec{r}$ which maps a point from a space $\vec{x}$ to a space $\vec{u}$ the displacement gradient tensor is defined in 3D as

\[ { \newcommand{\pd}[2]{\frac{\partial #1}{\partial #2}} u_{i,j} = \left [ \begin{array}{ccc} \pd{u_0}{x_0} & \pd{u_0}{x_1} & \pd{u_0}{x_2} \\ \pd{u_1}{x_0} & \pd{u_1}{x_1} & \pd{u_1}{x_2} \\ \pd{u_2}{x_0} & \pd{u_2}{x_1} & \pd{u_2}{x_2} \end{array} \right] } \]

with

\[ \Delta r_i = u_{ij} r_j \]

where $\vec{u} = \left[u_0, u_1, u_2\right]^T$ and $\Delta \vec{r} = \left[x_0, x_1, x_2\right]^T$. The displacement gradient tensor matrix is just the rotation and independent scaling part of the affine transform that displaces the element.

Returns:
Conforming mesh object with tensor values attached to elements or NULL on error.
Parameters:
cObjGiven conforming mesh object.
invertInvert if non-zero, by default the tensors are computed for the transform from the mesh to the displaced mesh.
dstErrDestination error pointer, may be NULL.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzObject::type, _WlzObject::values, WLZ_CMESH_2D, WLZ_CMESH_3D, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_UNIMPLEMENTED, and WLZ_ERR_VALUES_NULL.

Referenced by WlzCMeshStrainTensor().

WlzObject* WlzCMeshDGTensorAtPts ( WlzObject cObj,
int  invert,
WlzDVertex3  sd,
int  dither,
WlzErrorNum dstErr 
)

Given a conforming mesh transform this function computes the displacement gradient tensor at regular cartesian grid sample points throughout the mesh. The tensor values at the sample points are computed at each point by computing an iverse distance weighted least squares general affine transform for the ring of nodes surrounding the closest node. See WlzCMeshDGTensor() for the description of the tensor.

Returns:
Points object with tensor values at the point locations or NULL on error.
Parameters:
cObjGiven conforming mesh object.
invertInvert if non-zero, by default the tensors are computed for the transform from the mesh to the displaced mesh.
sdSample distance.
ditherDither the sample locations.
dstErrDestination error pointer, may be NULL.

References _WlzValues::core, _WlzDomain::core, _WlzObject::domain, _WlzObject::type, _WlzObject::values, _WlzDVertex3::vtX, _WlzDVertex3::vtY, _WlzDVertex3::vtZ, WLZ_CMESH_2D, WLZ_CMESH_3D, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_PARAM_DATA, WLZ_ERR_UNIMPLEMENTED, WLZ_ERR_VALUES_NULL, and WLZ_MESH_TOLERANCE.

Referenced by WlzCMeshStrainTensorAtPts().

WlzObject* WlzCMeshStrainTensor ( WlzObject cObj,
int  invert,
WlzErrorNum dstErr 
)

Given a conforming mesh transform this function computes the strain tensor for each of it's valid elements. This function uses WlzCMeshDGTensor() to compute the displacement gradient tensor and the derives the strain tensor $e_{ij}$ from this using:

\[ e_{ij} = \frac{1}{2} (u_{ij} + u_{ji}) \]

.

Returns:
Conforming mesh object with tensor values attached to elements or NULL on error.
Parameters:
cObjGiven conforming mesh object.
invertInvert if non-zero, by default the tensors are computed for the transform from the mesh to the displaced mesh.
dstErrDestination error pointer, may be NULL.

References _WlzObject::type, WLZ_CMESH_2D, WLZ_CMESH_3D, WLZ_ERR_NONE, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_UNIMPLEMENTED, WlzCMeshDGTensor(), and WlzFreeObj().

WlzObject* WlzCMeshStrainTensorAtPts ( WlzObject cObj,
int  invert,
WlzDVertex3  sd,
int  dither,
WlzErrorNum dstErr 
)

Given a conforming mesh transform this function computes the displacement gradient tensor at regular cartesian grid sample points throughout the mesh. The tensor values at the sample points are computed using WlzCMeshDGTensorAtPts(). The strain tensor is then computed from the displacement gradient tensor as for WlzCMeshStrainTensor().

Returns:
Points object with tensor values or NULL on error.
Parameters:
cObjGiven conforming mesh object.
invertInvert if non-zero, by default the tensors are computed for the transform from the mesh to the displaced mesh.
sdSample distance.
ditherDither the sample locations.
dstErrDestination error pointer, may be NULL.

References _WlzObject::type, WLZ_CMESH_2D, WLZ_CMESH_3D, WLZ_ERR_NONE, WLZ_ERR_OBJECT_TYPE, WLZ_ERR_UNIMPLEMENTED, WlzCMeshDGTensorAtPts(), and WlzFreeObj().

WlzVertexP WlzVerticesFromObj ( WlzObject obj,
WlzVertexP dstNr,
int *  dstCnt,
WlzVertexType dstType,
WlzErrorNum dstErr 
)

Allocates a buffer which it fills with the vertices from the given object. The object must be one of the types that is represented by vertices, eg polylines, boundlists and contours.

Returns:
Allocated vertices.
Parameters:
objGiven polygon domain object.
dstNrDestination ptr for normals. The normals will always be either WlzDVertex2 or WlzDVertex3. May be NULL.
dstCntDestination ptr for the number of vertices. Can NOT be NULL.
dstTypeDestination ptr for the type of vertices. Can NOT be NULL.
dstErrDestination error pointer, may be NULL.

References _WlzDomain::b, _WlzDomain::core, _WlzDomain::ctr, _WlzObject::domain, _WlzContour::model, _WlzDomain::poly, _WlzObject::type, _WlzVertexP::v, WLZ_2D_POLYGON, WLZ_BOUNDLIST, WLZ_CONTOUR, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, and WlzVerticesFromGM().

Referenced by WlzRegICPObjs(), WlzRegICPObjsGrd(), WlzRegICPObjWSD2D(), and WlzSnapFit().

WlzVertexP WlzVerticesFromObjBnd ( WlzObject obj,
int *  dstCnt,
WlzVertexType dstType,
WlzErrorNum dstErr 
)

Extracts all vertices that lie on the boundary of the given objects domain.

Returns:
Allocated vertices.
Parameters:
objGiven 2D or 3D domain object.
dstCntDestination ptr for the number of vertices.
dstTypeDestination ptr for the type of vertices, which will always be integer.
dstErrDestination error pointer, may be NULL.

References _WlzDomain::core, _WlzObject::domain, _WlzVertexP::i2, _WlzVertexP::i3, _WlzObject::type, _WlzVertexP::v, WLZ_2D_DOMAINOBJ, WLZ_3D_DOMAINOBJ, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WLZ_VERTEX_I2, WLZ_VERTEX_I3, WlzVerticesFromObjBnd2I(), and WlzVerticesFromObjBnd3I().

WlzErrorNum WlzVerticesFromObjBnd2I ( WlzObject obj,
int *  dstNVtx,
WlzIVertex2 **  dstVtx 
)

Extracts all vertices that lie on the boundary of the given 2D domain object's domain.

Returns:
Error code.
Parameters:
objGiven 2D domain object.
dstNVtxDestination ptr for the number of vertices, MUST NOT be NULL.
dstVtxDestination ptr for the vertices, MUST NOT be NULL.

References AlcMalloc(), _WlzDomain::b, _WlzObject::domain, _WlzVertexP::i2, _WlzVertexP::v, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_VERTEX_I2, WlzFreeObj(), and WlzObjToBoundary().

Referenced by WlzVerticesFromObjBnd().

WlzErrorNum WlzVerticesFromObjBnd3I ( WlzObject obj,
int *  dstNVtx,
WlzIVertex3 **  dstVtx 
)

Extracts all vertices that lie on the boundary of the given 3D domain object's domain.

Returns:
Error code.
Parameters:
objGiven 3D domain object.
dstNVtxDestination ptr for the number of vertices, MUST NOT be NULL.
dstVtxDestination ptr for the vertices, MUST NOT be NULL.

References _WlzDomain::core, _WlzObject::domain, _WlzObject::type, WLZ_3D_DOMAINOBJ, WLZ_6_CONNECTED, WLZ_ERR_DOMAIN_NULL, WLZ_ERR_NONE, WLZ_ERR_OBJECT_NULL, WLZ_ERR_OBJECT_TYPE, WlzDiffDomain(), WlzErosion(), WlzFreeObj(), and WlzVerticesFromObj3I().

Referenced by WlzVerticesFromObjBnd().

WlzErrorNum WlzVerticesFromObj2I ( WlzObject obj,
int *  dstNVtx,
WlzIVertex2 **  dstVtx 
)

Extracts all vertices that lie within the given 2D domain object's domain.

Returns:
Error code.
Parameters:
objGiven 2D domain object.
dstNVtxDestination ptr for the number of vertices, MUST NOT be NULL.
dstVtxDestination ptr for the vertices, MUST NOT be NULL.

References AlcMalloc(), _WlzIntervalWSpace::lftpos, _WlzIntervalWSpace::linpos, _WlzIntervalWSpace::rgtpos, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_EOO, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_RASTERDIR_ILIC, WlzArea(), WlzInitRasterScan(), and WlzNextInterval().

Referenced by WlzMarkerLattice(), and WlzVerticesFromObj3I().

WlzErrorNum WlzVerticesFromObj3I ( WlzObject obj,
int *  dstNVtx,
WlzIVertex3 **  dstVtx 
)

Extracts all vertices that lie within the given 3D domain object's domain.

Returns:
Error code.
Parameters:
objGiven 3D domain object.
dstNVtxDestination ptr for the number of vertices, MUST NOT be NULL.
dstVtxDestination ptr for the vertices, MUST NOT be NULL.

References AlcFree(), AlcMalloc(), _WlzValues::core, _WlzObject::domain, _WlzPlaneDomain::domains, _WlzDomain::p, _WlzPlaneDomain::plane1, _WlzIVertex2::vtX, WLZ_2D_DOMAINOBJ, WLZ_ERR_DOMAIN_DATA, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WlzFreeObj(), WlzMakeMain(), WlzVerticesFromObj2I(), and WlzVolume().

Referenced by WlzMarkerLattice(), and WlzVerticesFromObjBnd3I().

WlzVertexP WlzVerticesFromGM ( WlzGMModel model,
WlzVertexP dstNr,
int **  dstVId,
int *  dstCnt,
WlzVertexType dstType,
WlzErrorNum dstErr 
)

Allocates a buffer which it fills with the vertices from a geometric model.

Returns:
Allocated vertices.
Parameters:
modelGiven geometric model.
dstNrDestination ptr for normals, may be NULL.
dstVIdDestination ptr for the GM vertex indicies, may be NULL.
dstCntDestination ptr for the number of vertices.
dstTypeDestination ptr for the type of vertices.
dstErrDestination error pointer, may be NULL.

References _WlzGMModel::type, _WlzVertexP::v, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_GMMOD_2D, WLZ_GMMOD_2I, WLZ_GMMOD_2N, WLZ_GMMOD_3D, WLZ_GMMOD_3I, and WLZ_GMMOD_3N.

Referenced by WlzDistMetricDirGM(), WlzDistMetricGM(), WlzMatchICPCtr(), and WlzVerticesFromObj().

WlzVertexP WlzDVerticesFromGM ( WlzGMModel model,
int *  dstCnt,
WlzVertexType dstType,
WlzErrorNum dstErr 
)

Allocates a buffer which it fills with either 2D or 3D double precission vertices from the geometric model. The indicies of the vertices in the buffer are the same as the indices of the vertices in the model.

Returns:
Allocated vertices.
Parameters:
modelGiven geometric model.
dstCntDestination ptr for the number of vertices.
dstTypeDestination ptr for the type of vertices.
dstErrDestination error pointer, may be NULL.

References _WlzVertexP::d2, _WlzVertexP::d3, _WlzGMModel::type, _WlzVertexP::v, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_GMMOD_2D, WLZ_GMMOD_2I, WLZ_GMMOD_2N, WLZ_GMMOD_3D, WLZ_GMMOD_3I, WLZ_GMMOD_3N, WLZ_VERTEX_D2, and WLZ_VERTEX_D3.

Referenced by WlzGMFilterGeomLPLM().

WlzVertexP WlzDVerticesFromCMesh ( WlzCMeshP  mesh,
int *  dstCnt,
WlzVertexType dstType,
int  skip,
WlzErrorNum dstErr 
)

Allocates a buffer which it fills with either 2D or 3D double precission vertices from the conforming mesh. The indicies of the vertices in the buffer are the same as the indices of the vertices in the model unless the skip non-valid nodes flag is set..

Returns:
Allocated vertices.
Parameters:
meshGiven mesh.
dstCntDestination ptr for the number of vertices.
dstTypeDestination ptr for the type of vertices.
skipSkip non-valid nodes if non-zero.
dstErrDestination error pointer, may be NULL.

References _WlzVertexP::d2, _WlzVertexP::d3, _WlzCMeshP::m2, _WlzCMeshP::m2d5, _WlzCMeshP::m3, _WlzCMesh2D::type, _WlzCMeshP::v, _WlzVertexP::v, WLZ_CMESH_2D, WLZ_CMESH_2D5, WLZ_CMESH_3D, WLZ_ERR_DOMAIN_TYPE, WLZ_ERR_NONE, WLZ_VERTEX_D2, and WLZ_VERTEX_D3.

Referenced by WlzCMeshLPFilterLM().

AlcKDTTree* WlzVerticesBuildTree ( WlzVertexType  vType,
int  nV,
WlzVertexP  vtx,
int *  shfBuf,
WlzErrorNum dstErr 
)

Allocates and populates a k-D tree from the given vertices. The vertices are either WlzDVertex2 orWlzDVertex3.

Returns:
Woolz error code
Parameters:
vTypeType of vertices.
nVNumber of vertices.
vtxThe vertices.
shfBufWorkspace with at least nV ints used to shuffle vertices for randomized input to the K-D tree.
dstErrDestination error pointer, may be NULL.

References ALC_ER_NONE, ALC_POINTTYPE_DBL, AlcKDTInsert(), AlcKDTTreeNew(), AlgShuffleIdx(), _WlzVertexP::d2, _WlzVertexP::d3, _AlcKDTNode::idx, _WlzDVertex3::vtX, _WlzDVertex2::vtX, _WlzDVertex3::vtY, _WlzDVertex2::vtY, _WlzDVertex3::vtZ, WLZ_ERR_MEM_ALLOC, WLZ_ERR_NONE, WLZ_ERR_PARAM_TYPE, WLZ_VERTEX_D2, and WLZ_VERTEX_D3.

Referenced by WlzDistMetricDirVertex2D(), WlzDistMetricDirVertex3D(), and WlzMatchICPCtr().