int main(int argc, char **argv) { int tI, idN, option, con = WLZ_0_CONNECTED, nLo = 0, nHi = 0, maxSep = 1024, nObj = 0, ok = 1, usage = 0; char tC; double tD, mrkMass = 1.0, rad = 0.0; int tR[4]; WlzPixelV gV, bV; WlzBlobMark mrk = WLZ_BLOBMARK_CIRCLE; WlzObject *inObj = NULL, *outObj = NULL, *mrkObj = NULL; WlzObject **lObj = NULL; FILE *fP = NULL; char *inObjFileStr, *outObjFileStr; WlzErrorNum errNum = WLZ_ERR_NONE; const char *errMsg; static char optList[] = "c:g:G:hm:n:N:o:r:x:", fileStrDef[] = "-"; opterr = 0; memset(&gV, 0, sizeof(WlzPixelV)); bV.type = WLZ_GREY_UBYTE; bV.v.ubv = 0; gV.type = WLZ_GREY_ERROR; inObjFileStr = fileStrDef; outObjFileStr = fileStrDef; while((usage == 0) && ((option = getopt(argc, argv, optList)) != -1)) { switch(option) { case 'c': if(sscanf(optarg, "%d", &tI) != 1) { usage = 1; } else { switch(tI) { case 4: con = WLZ_4_CONNECTED; break; case 6: con = WLZ_6_CONNECTED; break; case 8: con = WLZ_8_CONNECTED; break; case 18: con = WLZ_18_CONNECTED; break; case 26: con = WLZ_26_CONNECTED; break; default: usage = 1; break; } } break; case 'g': switch(gV.type) { case WLZ_GREY_UBYTE: if((sscanf(optarg, "%d", &tI) != 1) || (tI < 0) || (tI > 255)) { usage = 1; } else { gV.v.ubv = tI; } break; case WLZ_GREY_SHORT: if((sscanf(optarg, "%d", &tI) != 1) || (tI < SHRT_MIN) || (tI > SHRT_MAX)) { usage = 1; } else { gV.v.shv = tI; } break; case WLZ_GREY_INT: if(sscanf(optarg, "%d", &tI) != 1) { usage = 1; } else { gV.v.inv = tI; } break; case WLZ_GREY_FLOAT: if((sscanf(optarg, "%lg", &tD) != 1) || (tD < -(FLT_MAX)) || (tD > FLT_MAX)) { usage = 1; } else { gV.v.flv = tD; } break; case WLZ_GREY_DOUBLE: if(sscanf(optarg, "%lg", &tD) != 1) { usage = 1; } else { gV.v.dbv = tD; } break; case WLZ_GREY_RGBA: tR[3] = 255; tR[0] = tR[1] = tR[2] = 0; if((sscanf(optarg, "%d,%d,%d,%d", &(tR[0]), &(tR[1]), &(tR[2]), &(tR[3])) == 0) || (tR[0] < 0) || (tR[0] > 255) || (tR[1] < 0) || (tR[1] > 255) || (tR[2] < 0) || (tR[2] > 255) || (tR[3] < 0) || (tR[3] > 255)) { usage = 1; } else { WLZ_RGBA_RGBA_SET(gV.v.rgbv, tR[0], tR[1], tR[2], tR[3]); } break; default: usage = 1; break; } break; case 'G': if(sscanf(optarg, "%c", &tC) != 1) { usage = 1; } switch(tC) { case 'v': gV.type = WLZ_GREY_ERROR; break; case 'u': gV.type = WLZ_GREY_UBYTE; break; case 's': gV.type = WLZ_GREY_SHORT; break; case 'i': gV.type = WLZ_GREY_INT; break; case 'f': gV.type = WLZ_GREY_FLOAT; break; case 'd': gV.type = WLZ_GREY_DOUBLE; break; case 'r': gV.type = WLZ_GREY_RGBA; break; default: usage = 1; break; } break; case 'm': if((sscanf(optarg, "%d", &tI) != 1) || ((tI != WLZ_BLOBMARK_CIRCLE) && (tI != WLZ_BLOBMARK_SQUARE))) { usage = 1; } else { mrk = (WlzBlobMark )tI; } break; case 'n': if((sscanf(optarg, "%d", &nLo) != 1) || (nLo < 0)) { usage = 1; } break; case 'N': if((sscanf(optarg, "%d", &nHi) != 1) || (nHi < 0)) { usage = 1; } break; case 'o': outObjFileStr = optarg; break; case 'r': if((sscanf(optarg, "%lg", &rad) != 1) || (rad < 0.0)) { usage = 1; } break; case 'x': if((sscanf(optarg, "%d", &maxSep) != 1) || (maxSep < 1)) { usage = 1; } case 'h': /* FALLTHROUGH */ default: usage = 1; break; } } if((usage == 0) && (nLo > nHi) && (nHi != 0)) { usage = 1; } if((usage == 0) && (optind < argc)) { if((optind + 1) != argc) { usage = 1; } else { inObjFileStr = *(argv + optind); } } ok = (usage == 0); /* Read input domain object. */ if(ok) { if((inObjFileStr == NULL) || (*inObjFileStr == '\0') || ((fP = (strcmp(inObjFileStr, "-")? fopen(inObjFileStr, "r"): stdin)) == NULL) || ((inObj = WlzAssignObject(WlzReadObj(fP, &errNum), NULL)) == NULL) || (errNum != WLZ_ERR_NONE)) { ok = 0; } if(fP) { if(strcmp(inObjFileStr, "-")) { (void )fclose(fP); } fP = NULL; } } /* Check object type and connectivity. */ if(ok) { switch(inObj->type) { case WLZ_2D_DOMAINOBJ: switch(con) { case WLZ_0_CONNECTED: con = WLZ_8_CONNECTED; break; case WLZ_4_CONNECTED: /* FALLTHROUGH */ case WLZ_8_CONNECTED: break; default: ok = 0; errNum = WLZ_ERR_PARAM_DATA; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Connectivity for 2D must be 4 or 8 (%s).\n", *argv, errMsg); break; } break; case WLZ_3D_DOMAINOBJ: switch(con) { case WLZ_0_CONNECTED: con = WLZ_26_CONNECTED; break; case WLZ_6_CONNECTED: /* FALLTHROUGH */ case WLZ_18_CONNECTED: /* FALLTHROUGH */ case WLZ_26_CONNECTED: break; default: ok = 0; errNum = WLZ_ERR_PARAM_DATA; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Connectivity for 3D must be 6, 18 or 26 (%s).\n", *argv, errMsg); break; } break; default: ok = 0; errNum = WLZ_ERR_OBJECT_TYPE; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Input object must either a 2 or 3D domain object (%s).\n", *argv, errMsg); break; } } /* Make basic marker with centre at the origin. */ if(ok) { double mrkRad; if(rad > 0.5) { mrkRad = rad; } else { mrkRad = 127; } if(mrk == WLZ_BLOBMARK_SQUARE) { mrkObj = WlzMakeCuboidObject(inObj->type, mrkRad, mrkRad, mrkRad, 0, 0, 0, &errNum); } else /* mrk = WLZ_BLOBMARK_CIRCLE */ { mrkObj = WlzMakeSphereObject(inObj->type, mrkRad, 0, 0, 0, &errNum); } if(mrkObj == NULL) { ok = 0; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Failed to create basic marker object (%s).\n", *argv, errMsg); } else { mrkMass = WlzVolume(mrkObj, NULL); } } /* Label the given domain. */ if(ok) { errNum = WlzLabel(inObj, &nObj, &lObj, maxSep, 1, con); if((errNum != WLZ_ERR_NONE) || (nObj == 0)) { ok = 0; if(errNum == WLZ_ERR_NONE) { errNum = WLZ_ERR_DOMAIN_DATA; } (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Failed to split the given object into separate regions (%s)\n", *argv, errMsg); } } /* Work through the separate object list removing small/large objects * according to the low and high thresholds. */ if(ok) { int idM; for(idN = 0, idM = 0; idN < nObj; ++idN) { int vol; vol = WlzVolume(lObj[idN], &errNum); if(errNum == WLZ_ERR_NONE) { if(((nLo > 0) && (vol < nLo)) || ((nHi > 0) && (vol > nHi))) { (void )WlzFreeObj(lObj[idN]); } else { lObj[idM] = lObj[idN]; ++idM; } } } nObj = idM; if(nObj == 0) { ok = 0; errNum = WLZ_ERR_DOMAIN_DATA; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Failed to find and separate regions (%s)\n", *argv, errMsg); } } /* Build a marker object by adding a mark at the centre of mass of each * separate fragment. */ if(ok) { WlzObject *obj0 = NULL; idN = 0; obj0 = WlzMakeEmpty(&errNum); while((errNum == WLZ_ERR_NONE) && (idN < nObj)) { double mass; WlzDVertex3 com; WlzObject *obj1 = NULL, *obj2 = NULL; WlzAffineTransform *tr = NULL; com = WlzCentreOfMass3D(lObj[idN], 1, &mass, &errNum); if(errNum == WLZ_ERR_NONE) { double s; if(rad < 0.5) { double t; t = mass / mrkMass; if(inObj->type == WLZ_2D_DOMAINOBJ) { s = sqrt(t); } else /* inObj->type == WLZ_3D_DOMAINOBJ */ { s = cbrt(t); } } else { s = 1.0; } tr = (inObj->type == WLZ_2D_DOMAINOBJ)? WlzAffineTransformFromPrimVal( WLZ_TRANSFORM_2D_AFFINE, com.vtX, com.vtY, 0.0, s, 0.0, 0.0, 0.0, 0.0, 0.0, 0, &errNum): WlzAffineTransformFromPrimVal( WLZ_TRANSFORM_3D_AFFINE, com.vtX, com.vtY, com.vtZ, s, 0.0, 0.0, 0.0, 0.0, 0.0, 0, &errNum); } if(errNum == WLZ_ERR_NONE) { obj1 = WlzAffineTransformObj(mrkObj, tr, WLZ_INTERPOLATION_NEAREST, &errNum); } if(errNum == WLZ_ERR_NONE) { obj2 = WlzUnion2(obj0, obj1, &errNum); } if(errNum == WLZ_ERR_NONE) { (void )WlzFreeObj(obj0); obj0 = obj2; obj2 = NULL; } (void )WlzFreeObj(obj1); (void )WlzFreeObj(obj2); (void )WlzFreeAffineTransform(tr); ++idN; } if(errNum == WLZ_ERR_NONE) { WlzValues val; WlzObjectType vTT; val.core = NULL; if(gV.type != WLZ_GREY_ERROR) { vTT = WlzGreyTableType(WLZ_GREY_TAB_RAGR, gV.type, NULL); if(inObj->type == WLZ_2D_DOMAINOBJ) { val.v = WlzNewValueTb(obj0, vTT, bV, &errNum); } else /* inObj->type == WLZ_3D_DOMAINOBJ */ { val.vox = WlzNewValuesVox(obj0, vTT, bV, &errNum); } } if(errNum == WLZ_ERR_NONE) { outObj = WlzMakeMain(inObj->type, obj0->domain, val, NULL, NULL, &errNum); } if((errNum == WLZ_ERR_NONE) && (gV.type != WLZ_GREY_ERROR)) { errNum = WlzGreySetValue(outObj, gV); } } } if(ok) { errNum = WLZ_ERR_WRITE_EOF; if(((fP = (strcmp(outObjFileStr, "-")? fopen(outObjFileStr, "w"): stdout)) == NULL) || ((errNum = WlzWriteObj(fP, outObj)) != WLZ_ERR_NONE)) { ok = 0; (void )WlzStringFromErrorNum(errNum, &errMsg); (void )fprintf(stderr, "%s: Failed to write output object (%s).\n", *argv, errMsg); } if(fP && strcmp(outObjFileStr, "-")) { (void )fclose(fP); } } (void )WlzFreeObj(inObj); if(lObj != NULL) { for(idN = 0; idN < nObj; ++idN) { (void )WlzFreeObj(lObj[idN]); } AlcFree(lObj); } (void )WlzFreeObj(outObj); if(usage) { (void )fprintf(stderr, "Usage: %s%sExample: %s%s", *argv, " [-c#] [-g#] [-G#] [-h] [-m#] [-n#] [-N#]\n" " [-o<output object>] [-r#]] [-x#] [<input object>]\n" "Options:\n" " -c Connectivity: 4, 6, 8, 18 or 26 connected (default 8 for 2D\n" " domains and 26 for 3D domains).\n" " -g Grey value for marker. This is a single number for all except\n" " RGBA (colour) grey values. RGBA components must be separated by\n" " by a comma.\n" " -G Grey value type for marker specified by letter:\n" " v no grey values (default).\n" " u unsigned byte grey values.\n" " s short grey values.\n" " i int grey values.\n" " f int grey values.\n" " d int grey values.\n" " r red, green, blue, alpha grey values.\n" " -h Help, prints usage message.\n" " -m Marker type specified by a number:\n" " 1 circle/sphere (default)\n" " 2 square/cube\n" " -n Threshold minimum area/volume of blob for a marker (default\n" " >= 1).\n" " -N Threshold maximum area/volume of blob for a marker. If zero\n" " there is no upper limit. (default 0).\n" " -o Output object file.\n" " -r Marker radius. Attempts to keep the same area/volume if zero.\n" " (default 0).\n" " -x Maximum number of separate regions in the object (default 1024).\n" "Reads a spatial domain object and replaces each spatialy separate\n" "region with a marker placed at the centre of mass of the region.\n" "All files are read from the standard input and written to the standard\n" "output unless filenames are given.\n" "If grey values are required then the grey value type must be set before\n" "the actual grey value.\n", *argv, " -o out.wlz -n 4 -r 10 -G r -g 200,100,0,255 in.wlz\n" "A spatial domain object is read from the file in.wlz and each\n" "spatialy separate region of the domain is replaced by a circle or\n" "sphere of radius 10 (pixels). All small regions with less than four\n" "(pixels voxels) is ignored. The output object (with grey values set\n" "to orange) is written to the file out.wlz.\n"); } return(!ok); }
/*! * \return Distance object which shares the given foreground object's * domain and has integer distance values, null on error. * \ingroup WlzMorphologyOps * \brief Computes the distance of every pixel/voxel in the foreground * object from the reference object. * * A distance transform maps all position within a forground * domain to their distances from a reference domain. * The distance transforms implemented within this function * use efficient morphological primitives. * * Given two domains, * \f$\Omega_r\f$ the reference domain and \f$\Omega_f\f$ * the domain specifying the region of interest, * a domain with a thin shell \f$\Omega_i\f$ * is iteratively expanded from it's initial domain * corresponding to the reference domain \f$\Omega_r\f$. * At each iteration * \f$\Omega_i\f$ is dilated and clipped * by it's intersection with \f$\Omega_f\f$ until \f$\Omega_i\f$ * becomes the null domain \f$\emptyset\f$. * At each iteration the current distance is recorded in a value * table which * covers the domain \f$\Omega_f\f$. * * An octagonal distance scheme may be used in which * the distance metric is alternated between 4 and 8 * connected for 2D and 6 and 26 connectivities in 3D. * See: G. Borgefors. "Distance Transformations in Arbitrary * Dimensions" CVGIP 27:321-345, 1984. * * An approximate Euclidean distance transform may be computed * by: Scaling the given foreground and reference objects using * the given approximation scale parameter, dilating the * reference domain using a sphere with a radius having the same * value as the scale parameter and then finaly sampling the * scaled distances. * \param forObj Foreground object. * \param refObj Reference object. * \param dFn Distance function which must be * appropriate to the dimension of * the foreground and reference objects. * \param dParam Parameter required for distance * function. Currently only * WLZ_APX_EUCLIDEAN_DISTANCE requires a * parameter. In this case the parameter * is the approximation scale. * \param dstErr Destination error pointer, may be NULL. */ WlzObject *WlzDistanceTransform(WlzObject *forObj, WlzObject *refObj, WlzDistanceType dFn, double dParam, WlzErrorNum *dstErr) { int idP, lastP, dim, notDone = 1; double scale; WlzObject *tmpObj, *sObj = NULL, *sForObj = NULL, *sRefObj = NULL, *dilObj = NULL, *dstObj = NULL, *difObj = NULL, *curItrObj = NULL; WlzObject *bothObj[2]; WlzDomain *difDoms; WlzPixelV dstV, bgdV; WlzValues *difVals; WlzAffineTransform *tr = NULL; WlzConnectType con; WlzObjectType dstGType; WlzErrorNum errNum = WLZ_ERR_NONE; WlzValues difVal, dstVal, nullVal; /* By defining WLZ_DIST_TRANSFORM_ENV these normalization parameters * are read from the environment. This is useful for optimization. */ #ifndef WLZ_DIST_TRANSFORM_ENV const #endif /* ! WLZ_DIST_TRANSFORM_ENV */ /* These normalizarion factors have been choosen to minimize the sum of * squares of the deviation of the distance values from Euclidean values * over a radius 100 circle or sphere, where the distances are computed * from the circumference of the sphere towards it's centre. The values * were established by experiment. */ double nrmDist4 = 0.97, nrmDist6 = 0.91, nrmDist8 = 1.36, nrmDist18 = 1.34, nrmDist26 = 1.60; #ifdef WLZ_DIST_TRANSFORM_ENV double val; char *envStr; if(((envStr = getenv("WLZ_DIST_TRANSFORM_NRMDIST4")) != NULL) && (sscanf(envStr, "%lg", &val) == 1)) { nrmDist4 = val; } if(((envStr = getenv("WLZ_DIST_TRANSFORM_NRMDIST6")) != NULL) && (sscanf(envStr, "%lg", &val) == 1)) { nrmDist6 = val; } if(((envStr = getenv("WLZ_DIST_TRANSFORM_NRMDIST8")) != NULL) && (sscanf(envStr, "%lg", &val) == 1)) { nrmDist8 = val; } if(((envStr = getenv("WLZ_DIST_TRANSFORM_NRMDIST18")) != NULL) && (sscanf(envStr, "%lg", &val) == 1)) { nrmDist18 = val; } if(((envStr = getenv("WLZ_DIST_TRANSFORM_NRMDIST26")) != NULL) && (sscanf(envStr, "%lg", &val) == 1)) { nrmDist26 = val; } #endif /* WLZ_DIST_TRANSFORM_ENV */ scale = dParam; nullVal.core = NULL; /* Check parameters. */ if((forObj == NULL) || (refObj == NULL)) { errNum = WLZ_ERR_OBJECT_NULL; } else if(((forObj->type != WLZ_2D_DOMAINOBJ) && (forObj->type != WLZ_3D_DOMAINOBJ)) || ((refObj->type != WLZ_POINTS) && (refObj->type != forObj->type))) { errNum = WLZ_ERR_OBJECT_TYPE; } else if((forObj->domain.core == NULL) || (refObj->domain.core == NULL)) { errNum = WLZ_ERR_DOMAIN_NULL; } if(errNum == WLZ_ERR_NONE) { bgdV.type = WLZ_GREY_INT; bgdV.v.inv = 0; dstV.type = WLZ_GREY_DOUBLE; dstV.v.dbv = 0.0; switch(forObj->type) { case WLZ_2D_DOMAINOBJ: switch(dFn) { case WLZ_4_DISTANCE: /* FALLTHROUGH */ case WLZ_8_DISTANCE: /* FALLTHROUGH */ case WLZ_OCTAGONAL_DISTANCE: /* FALLTHROUGH */ case WLZ_APX_EUCLIDEAN_DISTANCE: dim = 2; break; default: errNum = WLZ_ERR_PARAM_DATA; break; } break; case WLZ_3D_DOMAINOBJ: switch(dFn) { case WLZ_6_DISTANCE: /* FALLTHROUGH */ case WLZ_18_DISTANCE: /* FALLTHROUGH */ case WLZ_26_DISTANCE: /* FALLTHROUGH */ case WLZ_OCTAGONAL_DISTANCE: /* FALLTHROUGH */ case WLZ_APX_EUCLIDEAN_DISTANCE: dim = 3; break; default: errNum = WLZ_ERR_PARAM_DATA; break; } break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } if(errNum == WLZ_ERR_NONE) { switch(dFn) { case WLZ_4_DISTANCE: con = WLZ_4_CONNECTED; break; case WLZ_6_DISTANCE: con = WLZ_6_CONNECTED; break; case WLZ_8_DISTANCE: con = WLZ_8_CONNECTED; break; case WLZ_18_DISTANCE: con = WLZ_18_CONNECTED; break; case WLZ_26_DISTANCE: con = WLZ_26_CONNECTED; break; case WLZ_OCTAGONAL_DISTANCE: con = (dim == 2)? WLZ_8_CONNECTED: WLZ_26_CONNECTED; break; case WLZ_APX_EUCLIDEAN_DISTANCE: con = (dim == 2)? WLZ_8_CONNECTED: WLZ_26_CONNECTED; if(scale < 1.0) { errNum = WLZ_ERR_PARAM_DATA; } break; case WLZ_EUCLIDEAN_DISTANCE: errNum = WLZ_ERR_UNIMPLEMENTED; break; default: errNum = WLZ_ERR_PARAM_DATA; break; } } /* Create scaled domains and a sphere domain for structual erosion if the * distance function is approximate Euclidean. */ if(errNum == WLZ_ERR_NONE) { if(dFn == WLZ_APX_EUCLIDEAN_DISTANCE) { tr = (dim == 2)? WlzAffineTransformFromScale(WLZ_TRANSFORM_2D_AFFINE, scale, scale, 0.0, &errNum): WlzAffineTransformFromScale(WLZ_TRANSFORM_3D_AFFINE, scale, scale, scale, &errNum); if(errNum == WLZ_ERR_NONE) { tmpObj = WlzMakeMain(forObj->type, forObj->domain, nullVal, NULL, NULL, &errNum); if(tmpObj) { sForObj = WlzAssignObject( WlzAffineTransformObj(tmpObj, tr, WLZ_INTERPOLATION_NEAREST, &errNum), NULL); (void )WlzFreeObj(tmpObj); } } if(errNum == WLZ_ERR_NONE) { if(refObj->type == WLZ_POINTS) { sRefObj = WlzPointsToDomObj(refObj->domain.pts, scale, &errNum); } else /* type == WLZ_2D_DOMAINOBJ || type == WLZ_3D_DOMAINOBJ */ { tmpObj = WlzMakeMain(refObj->type, refObj->domain, nullVal, NULL, NULL, &errNum); if(errNum == WLZ_ERR_NONE) { sRefObj = WlzAssignObject( WlzAffineTransformObj(tmpObj, tr, WLZ_INTERPOLATION_NEAREST, &errNum), NULL); } } (void )WlzFreeObj(tmpObj); } if(errNum == WLZ_ERR_NONE) { sObj = WlzAssignObject( WlzMakeSphereObject(forObj->type, scale, 0.0, 0.0, 0.0, &errNum), NULL); } (void )WlzFreeAffineTransform(tr); } else { sForObj = WlzAssignObject( WlzMakeMain(forObj->type, forObj->domain, nullVal, NULL, NULL, &errNum), NULL); if(errNum == WLZ_ERR_NONE) { if(refObj->type == WLZ_POINTS) { sRefObj = WlzPointsToDomObj(refObj->domain.pts, 1.0, &errNum); } else { sRefObj = WlzAssignObject( WlzMakeMain(refObj->type, refObj->domain, nullVal, NULL, NULL, &errNum), NULL); } } } } /* Create new values for the computed distances. */ if(errNum == WLZ_ERR_NONE) { dstGType = WlzGreyTableType(WLZ_GREY_TAB_RAGR, WLZ_GREY_INT, NULL); if(dim == 2) { dstVal.v = WlzNewValueTb(sForObj, dstGType, bgdV, &errNum); } else { dstVal.vox = WlzNewValuesVox(sForObj, dstGType, bgdV, &errNum); } } /* Create a distance object using the foreground object's domain and * the new distance values. */ if(errNum == WLZ_ERR_NONE) { dstObj = WlzMakeMain(sForObj->type, sForObj->domain, dstVal, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { bothObj[0] = sForObj; errNum = WlzGreySetValue(dstObj, dstV); } /* Dilate the reference object while setting the distances in each * dilated shell. */ while((errNum == WLZ_ERR_NONE) && notDone) { if(dFn == WLZ_APX_EUCLIDEAN_DISTANCE) { dstV.v.dbv += 1.0; } else { switch(con) { case WLZ_4_CONNECTED: dstV.v.dbv += nrmDist4; break; case WLZ_6_CONNECTED: dstV.v.dbv += nrmDist6; break; case WLZ_8_CONNECTED: dstV.v.dbv += nrmDist8; break; case WLZ_18_CONNECTED: dstV.v.dbv += nrmDist18; break; case WLZ_26_CONNECTED: dstV.v.dbv += nrmDist26; break; default: errNum = WLZ_ERR_CONNECTIVITY_TYPE; break; } } if(dFn == WLZ_APX_EUCLIDEAN_DISTANCE) { dilObj = WlzStructDilation(sRefObj, sObj, &errNum); } else { dilObj = WlzDilation(sRefObj, con, &errNum); } if(errNum == WLZ_ERR_NONE) { switch(sForObj->type) { case WLZ_2D_DOMAINOBJ: curItrObj = WlzAssignObject( WlzIntersect2(dilObj, sForObj, &errNum), NULL); break; case WLZ_3D_DOMAINOBJ: bothObj[1] = dilObj; curItrObj = WlzAssignObject( WlzIntersectN(2, bothObj, 1, &errNum), NULL); break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } (void)WlzFreeObj(dilObj); /* Create difference object for the expanding shell. */ if(errNum == WLZ_ERR_NONE) { difObj = WlzDiffDomain(curItrObj, sRefObj, &errNum); } if((difObj == NULL) || WlzIsEmpty(difObj, &errNum)) { notDone = 0; } else { /* Assign the distance object's values to the difference object * and set all it's values to the current distance. */ if(errNum == WLZ_ERR_NONE) { switch(sForObj->type) { case WLZ_2D_DOMAINOBJ: difObj->values = WlzAssignValues(dstObj->values, NULL); errNum = WlzGreySetValue(difObj, dstV); break; case WLZ_3D_DOMAINOBJ: /* 3D is more complex than 2D: Need to create a temporary * voxel valuetable and assign the individual 2D values. */ difVal.vox = WlzMakeVoxelValueTb(WLZ_VOXELVALUETABLE_GREY, difObj->domain.p->plane1, difObj->domain.p->lastpl, bgdV, NULL, &errNum); if(errNum == WLZ_ERR_NONE) { difObj->values = WlzAssignValues(difVal, NULL); difDoms = difObj->domain.p->domains; difVals = difObj->values.vox->values; idP = difObj->domain.p->plane1; lastP = difObj->domain.p->lastpl; while(idP <= lastP) { if((*difDoms).core) { dstVal = dstObj->values.vox->values[idP - dstObj->domain.p->plane1]; *difVals = WlzAssignValues(dstVal, NULL); } ++idP; ++difDoms; ++difVals; } if(difObj->domain.p->lastpl > difObj->domain.p->plane1) { errNum = WlzGreySetValue(difObj, dstV); } } break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } (void )WlzFreeObj(sRefObj); sRefObj = WlzAssignObject(curItrObj, NULL); (void )WlzFreeObj(curItrObj); } (void )WlzFreeObj(difObj); difObj = NULL; if(dFn == WLZ_OCTAGONAL_DISTANCE) { /* Alternate connectivities for octagonal distance. */ if(dim == 2) { con = (con == WLZ_4_CONNECTED)? WLZ_8_CONNECTED: WLZ_4_CONNECTED; } else /* dim == 3 */ { con = (con == WLZ_6_CONNECTED)? WLZ_26_CONNECTED: WLZ_6_CONNECTED; } } } (void )WlzFreeObj(sObj); (void )WlzFreeObj(sForObj); (void )WlzFreeObj(sRefObj); (void )WlzFreeObj(curItrObj); if((errNum == WLZ_ERR_NONE) && (dFn == WLZ_APX_EUCLIDEAN_DISTANCE)) { tmpObj = WlzDistSample(dstObj, dim, scale, &errNum); (void )WlzFreeObj(dstObj); dstObj = tmpObj; } if(errNum != WLZ_ERR_NONE) { (void )WlzFreeObj(dstObj); dstObj = NULL; } if(dstErr) { *dstErr = errNum; } return(dstObj); }
/*! * \return Woolz error code. * \ingroup WlzBinaryOps * \brief Computes a metric which quantifies the extent to which * the domain of one of the given objects is offset from * the domain of the second. This is a symetric metric, ie * \f$OST(\Omega_0, \Omega_1) \equiv OST(\Omega_0, \Omega_1)\f$. * A domain is computed which is equidistant from the domains * of the two given objects, is within maxDist of each object's * domain and is within the convex hull of the union of the * domains of the two given objects. Within this domain the * 1st, 2nd and 3rd quantiles of the distance * (\f$q_0\f$, \f$q_1\f$ and \f$q_2\f$) are found. * The object's domains are classified as offset if * \f[ frac{q_1}{q_1 + (q_1 - q_0) + (q_2 - q_1)} \geq 0.5 \f] * ie * \f[ frac{q_1}{q_1 + q_2 - q_0} \geq 0.5 \f] * Small equi-distant domains with a volume less than half * the maximum distance do not classify the relationship as * an overlap. * \param o Array with the two given spatial * domain objects, must not be NULL * and nor must the objects. * \param t Array of temporary objects as in * WlzRCCTOIdx. * \param maxDist Maximum distance for offset. This * is used to compute a distance object, * large distances will significantly * increase the processing time. * \param dQ0 Destination pointer for 1st quantile * offset distance, must not be NULL. * \param dQ1 Destination pointer for 2nd quantile * (ie median) offset distance, * must not be NULL. * \param dQ2 Destination pointer for 3rd quantile * offset distance, must not be NULL. */ static WlzErrorNum WlzRCCOffset( WlzObject **o, WlzObject **t, int maxDist, int *dQ0, int *dQ1, int *dQ2) { int empty = 0; int q[3] = {0}; int *dHist = NULL; WlzObject *eObj = NULL; WlzErrorNum errNum = WLZ_ERR_NONE; if(((o[0]->type == WLZ_2D_DOMAINOBJ) && (o[1]->type == WLZ_2D_DOMAINOBJ)) || ((o[0]->type == WLZ_3D_DOMAINOBJ) && (o[1]->type == WLZ_3D_DOMAINOBJ))) { if((o[0]->domain.core == NULL) || (o[1]->domain.core == NULL)) { errNum = WLZ_ERR_DOMAIN_NULL; } } else { errNum = WLZ_ERR_OBJECT_TYPE; } /* Compute distance transforms of the two given objects out to a given * maximum distance and then using these distances the equi-distant * domain between these two objects. The values of the eqi-distant object * are those of the distance between the objects.*/ if(errNum == WLZ_ERR_NONE) { int i; WlzObject *sObj = NULL, *cObj = NULL; WlzObject *dObj[2]; dObj[0] = dObj[1] = NULL; /* Create structuring element with which to dilate the given object * domains(by maxDist). */ sObj = WlzAssignObject( WlzMakeSphereObject(o[0]->type, maxDist, 0, 0, 0, &errNum), NULL); /* Create domain or convex hull of the union of the two given object * domains. */ if(errNum == WLZ_ERR_NONE) { WlzObject *uObj = NULL, *xObj = NULL; errNum = WlzRCCMakeT(o, t, WLZ_RCCTOIDX_O0O1U); if(errNum == WLZ_ERR_NONE) { uObj = WlzAssignObject(t[WLZ_RCCTOIDX_O0O1U], NULL); } if(errNum == WLZ_ERR_NONE) { xObj = WlzAssignObject( WlzObjToConvexHull(uObj, &errNum), NULL); } if(errNum == WLZ_ERR_NONE) { cObj = WlzAssignObject( WlzConvexHullToObj(xObj, o[0]->type, &errNum), NULL); } (void )WlzFreeObj(xObj); (void )WlzFreeObj(uObj); } /* Dilate the two given objects and find the ntersection of the * dilated domains with each other and the convex hull computed * above. Within his domain compute the distances. */ if(errNum == WLZ_ERR_NONE) { for(i = 0; i < 2; ++i) { WlzObject *tObj = NULL, *rObj = NULL; tObj = WlzAssignObject( WlzStructDilation(o[i], sObj, &errNum), NULL); if(errNum == WLZ_ERR_NONE) { rObj = WlzAssignObject( WlzIntersect2(tObj, cObj, &errNum), NULL); } (void )WlzFreeObj(tObj); if(errNum == WLZ_ERR_NONE) { dObj[i] = WlzAssignObject( WlzDistanceTransform(rObj, o[!i], WLZ_OCTAGONAL_DISTANCE, 0.0, maxDist, &errNum), NULL); } (void )WlzFreeObj(rObj); if(errNum == WLZ_ERR_NONE) { WlzPixelV bgdV; bgdV.type = WLZ_GREY_INT; bgdV.v.inv = maxDist; errNum = WlzSetBackground(dObj[i], bgdV); } if(errNum != WLZ_ERR_NONE) { break; } } } /* Find the domain which is equi-distant from the two given objects, * within the xDist range and within the convex hull of the union of * the two given object's domains. */ (void )WlzFreeObj(sObj); sObj = NULL; if(errNum == WLZ_ERR_NONE) { WlzLong vol = 0; WlzObject *qObj = NULL, *tObj = NULL; qObj = WlzAssignObject( WlzImageArithmetic(dObj[0], dObj[1], WLZ_BO_EQ, 0, &errNum), NULL); if(errNum == WLZ_ERR_NONE) { WlzPixelV thrV; thrV.type = WLZ_GREY_INT; thrV.v.inv = 1; tObj = WlzAssignObject( WlzThreshold(qObj, thrV, WLZ_THRESH_HIGH, &errNum), NULL); } /* Check that the eqi-distant domain is of a reasonable size ie has * a area or volume greater than half the maximum distance. */ if(errNum == WLZ_ERR_NONE) { vol = WlzVolume(tObj, &errNum); if((maxDist / 2) >= vol) { empty = 1; } } if((errNum == WLZ_ERR_NONE) && !empty) { WlzObject *mObj; WlzPixelV tmpV; tmpV.type = WLZ_GREY_INT; tmpV.v.inv = 0; mObj = WlzAssignObject( WlzGreyTemplate(dObj[0], tObj, tmpV, &errNum), NULL); if(errNum == WLZ_ERR_NONE) { tmpV.v.inv = 1; eObj = WlzAssignObject( WlzThreshold(mObj, tmpV, WLZ_THRESH_HIGH, &errNum), NULL); } (void )WlzFreeObj(mObj); } (void )WlzFreeObj(tObj); (void )WlzFreeObj(qObj); if((errNum == WLZ_ERR_NONE) && !empty) { WlzLong vol; vol = WlzVolume(eObj, &errNum); if((maxDist / 2) >= vol) { empty = 1; } } } (void )WlzFreeObj(cObj); (void )WlzFreeObj(sObj); (void )WlzFreeObj(dObj[0]); (void )WlzFreeObj(dObj[1]); } /* Compute a quantised distance histogram in which equi-distant distances * are quantized to integer values. */ if((errNum == WLZ_ERR_NONE) && !empty) { if((dHist = (int *)AlcCalloc(maxDist + 1, sizeof(int))) == NULL) { errNum = WLZ_ERR_MEM_ALLOC; } } if((errNum == WLZ_ERR_NONE) && !empty) { if(eObj->type == WLZ_2D_DOMAINOBJ) { errNum = WlzRCCCompDistHist2D(maxDist, dHist, eObj); } else { errNum = WlzRCCCompDistHist3D(maxDist, dHist, eObj); } #ifdef WLZ_RCC_DEBUG_OST { FILE *fP; fP = fopen("WLZ_RCC_DEBUG_OST.wlz", "w"); (void )WlzWriteObj(fP, eObj); (void )fclose(fP); } #endif } WlzFreeObj(eObj); if((errNum == WLZ_ERR_NONE) && !empty) { int i, n, s, nq; /* Compute the median, first and third quantile offset distances, * the ratio of median to the median plus inner inter-quantile range. */ n = 0; for(i = 0; i < maxDist; ++i) { n += dHist[i]; } i = 0; s = 0; nq = n / 4; while(s < nq) { s += dHist[i++]; } q[0] = i; nq = n / 2; while(s < nq) { s += dHist[i++]; } q[1] = i; nq = (3 * n) / 4; while(s < nq) { s += dHist[i++]; } q[2] = i; } AlcFree(dHist); *dQ0 = q[0]; *dQ1 = q[1]; *dQ2 = q[2]; return(errNum); }
/*! * \return New domain object which coresponds to the union of * the given points. * \ingroup WlzFeatures * \brief Creates a domain object which coresponds to the union of * the given points. * \param pnt Point domain. * \param scale Scale, 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. * \param dstErr Destination error poiter, may be NULL. */ WlzObject *WlzPointsToDomObj(WlzPoints *pnt, double scale, WlzErrorNum *dstErr) { int idP; WlzObjectType dType; WlzObject *tObj0, *tObj1, *tObj2, *dObj = NULL; WlzVertex pos; WlzErrorNum errNum = WLZ_ERR_NONE; idP = 0; if(scale > DBL_EPSILON) { pos.d3.vtZ = 0.0; } else { pos.i3.vtZ = 0; } if(pnt == NULL) { errNum = WLZ_ERR_DOMAIN_NULL; } else { switch(pnt->type) { case WLZ_POINTS_2I: /* FALLTHROUGH */ case WLZ_POINTS_2D: dType = WLZ_2D_DOMAINOBJ; break; case WLZ_POINTS_3I: /* FALLTHROUGH */ case WLZ_POINTS_3D: dType = WLZ_3D_DOMAINOBJ; break; default: errNum = WLZ_ERR_DOMAIN_TYPE; break; } } if(errNum == WLZ_ERR_NONE) { tObj0 = WlzMakeEmpty(&errNum); } while((errNum == WLZ_ERR_NONE) && (idP < pnt->nPoints)) { if(scale > DBL_EPSILON) { switch(pnt->type) { case WLZ_POINTS_2I: pos.d3.vtX = pnt->points.i2[idP].vtX; pos.d3.vtY = pnt->points.i2[idP].vtY; break; case WLZ_POINTS_2D: pos.d3.vtX = pnt->points.d2[idP].vtX; pos.d3.vtY = pnt->points.d2[idP].vtY; pos.d3.vtZ = 0.0; break; case WLZ_POINTS_3I: pos.d3.vtX = pnt->points.i3[idP].vtX; pos.d3.vtY = pnt->points.i3[idP].vtY; pos.d3.vtZ = pnt->points.i3[idP].vtZ; break; case WLZ_POINTS_3D: pos.d3 = pnt->points.d3[idP]; break; default: errNum = WLZ_ERR_PARAM_TYPE; break; } if(errNum == WLZ_ERR_NONE) { tObj1 = WlzMakeSphereObject(dType, scale / 2.0, scale * pos.d3.vtX, scale * pos.d3.vtY, scale * pos.d3.vtZ, &errNum); } } else { switch(pnt->type) { case WLZ_POINTS_2I: pos.i3.vtX = pnt->points.i2[idP].vtX; pos.i3.vtY = pnt->points.i2[idP].vtY; break; case WLZ_POINTS_2D: pos.i3.vtX = WLZ_NINT(pnt->points.d2[idP].vtX); pos.i3.vtY = WLZ_NINT(pnt->points.d2[idP].vtY); break; case WLZ_POINTS_3I: pos.i3 = pnt->points.i3[idP]; break; case WLZ_POINTS_3D: pos.i3.vtX = WLZ_NINT(pnt->points.d3[idP].vtX); pos.i3.vtY = WLZ_NINT(pnt->points.d3[idP].vtY); pos.i3.vtZ = WLZ_NINT(pnt->points.d3[idP].vtZ); break; default: errNum = WLZ_ERR_PARAM_TYPE; break; } if(errNum == WLZ_ERR_NONE) { tObj1 = WlzMakeSinglePixelObject(dType, pos.i3.vtX, pos.i3.vtY, pos.i3.vtZ, &errNum); } } if(errNum == WLZ_ERR_NONE) { tObj2 = WlzUnion2(tObj0, tObj1, &errNum); } (void )WlzFreeObj(tObj0); (void )WlzFreeObj(tObj1); tObj0 = tObj2; ++idP; } if(errNum == WLZ_ERR_NONE) { dObj = tObj0; } else { (void )WlzFreeObj(tObj0); } if(dstErr) { *dstErr = errNum; } return(dObj); }