/*! * \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 New Woolz domain object with gradient grey values or NULL * on error. * \ingroup WlzValuesFilters * \brief Computes the magnitude of the gray values in the * 3 given Woolz 3D domain objects. * \param srcObj0 First object. * \param srcObj1 Second object. * \param srcObj2 Third object. * \param dstErr Destination error pointer, may * be null. */ static WlzObject *WlzGreyMagnitude3D(WlzObject *srcObj0, WlzObject *srcObj1, WlzObject *srcObj2, WlzErrorNum *dstErr) { int idN, idP, nPlanes; WlzDomain dstDom; WlzValues dstVal; WlzObject *dstObj2D, *istObj = NULL, *dstObj = NULL; WlzObject *srcObjA[3], *iObj3DA[3], *iObj2DA[3]; WlzPixelV bgdV; WlzErrorNum errNum = WLZ_ERR_NONE; dstDom.core = NULL; dstVal.core = NULL; iObj3DA[0] = iObj3DA[1] = iObj3DA[2] = NULL; if((srcObj0->domain.core->type != WLZ_PLANEDOMAIN_DOMAIN) || (srcObj1->domain.core->type != WLZ_PLANEDOMAIN_DOMAIN) || (srcObj2->domain.core->type != WLZ_PLANEDOMAIN_DOMAIN)) { errNum = WLZ_ERR_DOMAIN_TYPE; } else if ((srcObj0->values.core->type != WLZ_VOXELVALUETABLE_GREY) || (srcObj1->values.core->type != WLZ_VOXELVALUETABLE_GREY) || (srcObj2->values.core->type != WLZ_VOXELVALUETABLE_GREY)) { errNum = WLZ_ERR_VALUES_TYPE; } else { srcObjA[0] = srcObj0; srcObjA[1] = srcObj1; srcObjA[2] = srcObj2; istObj = WlzIntersectN(3, srcObjA, 0, &errNum); } if(errNum == WLZ_ERR_NONE) { if(istObj->type == WLZ_EMPTY_OBJ) { dstObj = istObj; } else { dstDom = WlzAssignDomain(istObj->domain, NULL); WlzFreeObj(istObj); idN = 0; while((errNum == WLZ_ERR_NONE) && (idN < 3)) { iObj3DA[idN] = WlzMakeMain(WLZ_3D_DOMAINOBJ, dstDom, srcObjA[idN]->values, NULL, NULL, &errNum); ++idN; } if(errNum == WLZ_ERR_NONE) { idP = 0; nPlanes = dstDom.p->lastpl - dstDom.p->plane1 + 1; while((errNum == WLZ_ERR_NONE) && (idP < nPlanes)) { if((dstDom.p->domains + idN)->core == NULL) { (dstVal.vox->values + idP)->core = NULL; } else { dstObj2D = NULL; iObj2DA[0] = iObj2DA[1] = iObj2DA[2] = NULL; idN = 0; /* Make 2D objects. */ while((errNum == WLZ_ERR_NONE) && (idN < 3)) { iObj2DA[idN] = WlzMakeMain(WLZ_2D_DOMAINOBJ, *(iObj3DA[idN]->domain.p->domains + idP), *(iObj3DA[idN]->values.vox->values + idP), NULL, NULL, &errNum); ++idN; } /* Compute magnitude. */ if(errNum == WLZ_ERR_NONE) { dstObj2D = WlzGreyMagnitude2D3(iObj2DA[0], iObj2DA[1], iObj2DA[2], &errNum); } /* If first plane get grey type and background, then create new * 3D object with new voxel values. */ if((idP == 0) && (errNum == WLZ_ERR_NONE)) { bgdV = WlzGetBackground(iObj2DA[0], &errNum); if(errNum == WLZ_ERR_NONE) { dstVal.vox = WlzMakeVoxelValueTb(WLZ_VOXELVALUETABLE_GREY, dstDom.p->plane1, dstDom.p->lastpl, bgdV, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { dstObj = WlzMakeMain(WLZ_3D_DOMAINOBJ, dstDom, dstVal, NULL, NULL, &errNum); } } /* Add 2D values to the 3D objects values. */ if(errNum == WLZ_ERR_NONE) { *(dstVal.vox->values + idP) = WlzAssignValues(dstObj2D->values, NULL); } /* Free temporary 2D object. */ if(dstObj2D) { WlzFreeObj(dstObj2D); } idN = 0; while(idN < 3) { if(iObj2DA[idN]) { WlzFreeObj(iObj2DA[idN]); } ++idN; } } ++idP; } } /* Free 3D intersection gradient objects. */ idN = 0; while(idN < 3) { if(iObj3DA[idN]) { WlzFreeObj(iObj3DA[idN]); } ++idN; } /* Decrement link count to the destination 3D domain. */ (void )WlzFreeDomain(dstDom); } } if((errNum != WLZ_ERR_NONE) && (dstObj != NULL)) { WlzFreeObj(dstObj); dstObj = NULL; } if(dstErr) { *dstErr = errNum; } return(dstObj); }
/*! * \return * \brief Computes the maximal domain and gradient direction of * given Woolz domain object. * * Currently only implemented for 2D domain objects. * The domain is the maximaly suppressed domain and the values * are the encoded gradient direction. The direction is encoding * is from the +ve x-axis counter clockwise in eight steps * with a mask of 0x80, ie directions values are in the * range 128 -> 128 + 7. * \verbatim ^ Y axis (downwards when displayed) | +----------+---------+ | \ | /| | \128 + 2|128 + 1/ | | \ | / | | \ | / | | \ | / | | \ | / | | \ | / | |128 + 3 \ | /128 + 0| | \|/ | +----------O---------+--> X axis | /|\ | |128 + 4 / | \128 + 7| | / | \ | | / | \ | | / | \ | | / | \ | | / | \ | | /128 + 5|128 + 6\ | | / | \| +----------+---------+ \endverbatim * \param grdM Gradient magnitude. * \param grdZ Gradient (partial derivative) * through planes. * \param grdY Gradient (partial derivative) * through lines. * \param grdX Gradient (partial derivative) * through columns. * \param minThrV Minimum gradient value to * consider. * \param dstErr Destination error pointer, may * be null. */ WlzObject *WlzNMSuppress(WlzObject *grdM, WlzObject *grdZ, WlzObject *grdY, WlzObject *grdX, WlzPixelV minThrV, WlzErrorNum *dstErr) { int idN; WlzObject *istObj = NULL, *dstObj = NULL; WlzErrorNum errNum = WLZ_ERR_NONE; WlzObject *tObjs[4]; if(grdM == NULL) { errNum = WLZ_ERR_OBJECT_NULL; } else { switch(grdM->type) { case WLZ_2D_DOMAINOBJ: tObjs[0] = grdM; tObjs[1] = grdY; tObjs[2] = grdX; istObj = WlzIntersectN(3, tObjs, 0, &errNum); if(errNum == WLZ_ERR_NONE) { switch(istObj->type) { case WLZ_2D_DOMAINOBJ: tObjs[1] = tObjs[2] = NULL; tObjs[0] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdM->values, NULL, NULL, &errNum); if(errNum == WLZ_ERR_NONE) { tObjs[1] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdY->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { tObjs[2] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdX->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { dstObj = WlzNMSuppress2D(tObjs[0], tObjs[1], tObjs[2], minThrV, &errNum); } for(idN = 0; idN < 3; ++idN) { if(tObjs[0]) { WlzFreeObj(tObjs[0]); } } break; case WLZ_EMPTY_OBJ: dstObj = WlzMakeEmpty(&errNum); break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } break; case WLZ_3D_DOMAINOBJ: tObjs[0] = grdM; tObjs[1] = grdZ; tObjs[2] = grdY; tObjs[3] = grdX; istObj = WlzIntersectN(4, tObjs, 0, &errNum); if(errNum == WLZ_ERR_NONE) { switch(istObj->type) { case WLZ_2D_DOMAINOBJ: tObjs[1] = tObjs[2] = tObjs[3] = NULL; tObjs[0] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdM->values, NULL, NULL, &errNum); if(errNum == WLZ_ERR_NONE) { tObjs[1] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdZ->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { tObjs[2] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdY->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { tObjs[3] = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, grdX->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { dstObj = WlzNMSuppress3D(tObjs[0], tObjs[1], tObjs[2], tObjs[3], minThrV, &errNum); } for(idN = 0; idN < 4; ++idN) { if(tObjs[0]) { WlzFreeObj(tObjs[0]); } } break; case WLZ_EMPTY_OBJ: dstObj = WlzMakeEmpty(&errNum); break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } break; case WLZ_EMPTY_OBJ: dstObj = WlzMakeEmpty(&errNum); break; default: errNum = WLZ_ERR_OBJECT_TYPE; break; } } if(dstErr) { *dstErr = errNum; } return(dstObj); }
/*! * \return New Woolz domain object with gradient grey values or NULL on * error. * \ingroup WlzValuesFilters * \brief Computes the magnitude of the gray values in the * 3 given Woolz 2D domain objects. * \param srcObj0 First object. * \param srcObj1 Second object. * \param srcObj2 Third object. * \param dstErr Destination error pointer, may be null. */ static WlzObject *WlzGreyMagnitude2D3(WlzObject *srcObj0, WlzObject *srcObj1, WlzObject *srcObj2, WlzErrorNum *dstErr) { int idN, itvLen, bufSz; double **iBufA = NULL; WlzObject *tObj0, *istObj = NULL, *dstObj = NULL; WlzObject *iObjA[3], *tObjA[3]; WlzGreyP tGP0; WlzGreyType dstGType = WLZ_GREY_ERROR; WlzGreyType gTypeA[3]; WlzPixelV dstBgd; WlzPixelV bgdA[3]; WlzValues dstVal; WlzIntervalWSpace dstIWSp; WlzGreyWSpace dstGWSp; WlzIntervalWSpace iIWSpA[3]; WlzGreyWSpace iGWSpA[3]; WlzErrorNum errNum = WLZ_ERR_NONE; *(iObjA + 0) = *(iObjA + 1) = *(iObjA + 2) = NULL; /* Check source objects. */ if((srcObj0 == NULL) || (srcObj1 == NULL) ||(srcObj2 == NULL)) { errNum = WLZ_ERR_OBJECT_NULL;; } else if((srcObj0->type != WLZ_2D_DOMAINOBJ) || (srcObj1->type != WLZ_2D_DOMAINOBJ) || (srcObj2->type != WLZ_2D_DOMAINOBJ)) { errNum = WLZ_ERR_OBJECT_TYPE; } else if((srcObj0->domain.core == NULL) || (srcObj1->domain.core == NULL) || (srcObj2->domain.core == NULL)) { errNum = WLZ_ERR_DOMAIN_NULL; } else if((srcObj0->values.core == NULL) || (srcObj1->values.core == NULL) || (srcObj2->values.core == NULL)) { errNum = WLZ_ERR_VALUES_NULL; } /* Compute the intersection of the source objects. */ if(errNum == WLZ_ERR_NONE) { *(tObjA + 0) = srcObj0; *(tObjA + 1) = srcObj1; *(tObjA + 2) = srcObj2; istObj = WlzIntersectN(3, tObjA, 0, &errNum); } if(errNum == WLZ_ERR_NONE) { *(iObjA + 0) = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, srcObj0->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { *(iObjA + 1) = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, srcObj1->values, NULL, NULL, &errNum); } if(errNum == WLZ_ERR_NONE) { *(iObjA + 2) = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, srcObj2->values, NULL, NULL, &errNum); } /* Get background value and grey types */ idN = 0; while((errNum == WLZ_ERR_NONE) && (idN < 3)) { tObj0 = *(iObjA + idN); *(gTypeA + idN) = WlzGreyTableTypeToGreyType(tObj0->values.core->type, &errNum); if(errNum == WLZ_ERR_NONE) { *(bgdA + idN) = WlzGetBackground(tObj0, &errNum); } ++idN; } /* Promote grey types. */ if(errNum == WLZ_ERR_NONE) { if((*(gTypeA + 0) == WLZ_GREY_DOUBLE) || (*(gTypeA + 1) == WLZ_GREY_DOUBLE) || (*(gTypeA + 2) == WLZ_GREY_DOUBLE)) { dstGType = WLZ_GREY_DOUBLE; } else if((*(gTypeA + 0) == WLZ_GREY_FLOAT) || (*(gTypeA + 1) == WLZ_GREY_FLOAT) || (*(gTypeA + 2) == WLZ_GREY_FLOAT)) { dstGType = WLZ_GREY_FLOAT; } else if((*(gTypeA + 0) == WLZ_GREY_INT) || (*(gTypeA + 1) == WLZ_GREY_INT) || (*(gTypeA + 2) == WLZ_GREY_INT)) { dstGType = WLZ_GREY_INT; } else if((*(gTypeA + 0) == WLZ_GREY_SHORT) || (*(gTypeA + 1) == WLZ_GREY_SHORT) || (*(gTypeA + 2) == WLZ_GREY_SHORT)) { dstGType = WLZ_GREY_SHORT; } else if((*(gTypeA + 0) == WLZ_GREY_UBYTE) || (*(gTypeA + 1) == WLZ_GREY_UBYTE) || (*(gTypeA + 2) == WLZ_GREY_UBYTE)) { dstGType = WLZ_GREY_SHORT; } else { /* RGBA to be done RAB */ errNum = WLZ_ERR_GREY_TYPE; } } /* Make destination object with intersection domain and new values. */ if(errNum == WLZ_ERR_NONE) { (void )WlzValueConvertPixel(&dstBgd, *(bgdA + 0), dstGType); dstVal.v = WlzNewValueTb(*(iObjA + 0), WlzGreyValueTableType(0, WLZ_GREY_TAB_RAGR, dstGType, NULL), dstBgd, &errNum); if(errNum == WLZ_ERR_NONE) { dstObj = WlzMakeMain(WLZ_2D_DOMAINOBJ, istObj->domain, dstVal, NULL, NULL, &errNum); } } if(istObj) { WlzFreeObj(istObj); } /* Make buffers. */ if(errNum == WLZ_ERR_NONE) { bufSz = dstObj->domain.i->lastkl - dstObj->domain.i->kol1 + 1; if(AlcDouble2Malloc(&iBufA, 3, bufSz) != ALC_ER_NONE) { errNum = WLZ_ERR_MEM_ALLOC; } } /* Scan through the objects computing the magnitude. */ if(errNum == WLZ_ERR_NONE) { idN = 0; while((errNum == WLZ_ERR_NONE) && (idN < 3)) { errNum = WlzInitGreyScan(*(iObjA + idN), iIWSpA + idN, iGWSpA + idN); ++idN; } if(errNum == WLZ_ERR_NONE) { errNum = WlzInitGreyScan(dstObj, &dstIWSp, &dstGWSp); } while((errNum == WLZ_ERR_NONE) && ((errNum = WlzNextGreyInterval(iIWSpA + 0)) == WLZ_ERR_NONE) && ((errNum = WlzNextGreyInterval(iIWSpA + 1)) == WLZ_ERR_NONE) && ((errNum = WlzNextGreyInterval(iIWSpA + 2)) == WLZ_ERR_NONE) && ((errNum = WlzNextGreyInterval(&dstIWSp)) == WLZ_ERR_NONE)) { itvLen = dstIWSp.rgtpos - dstIWSp.lftpos + 1; /* Copy intervals to double buffers. */ idN = 0; while(idN < 3) { tGP0.dbp = *(iBufA + idN); WlzValueCopyGreyToGrey(tGP0, 0, WLZ_GREY_DOUBLE, (iGWSpA + idN)->u_grintptr, 0, (iGWSpA + idN)->pixeltype, itvLen); ++idN; } /* Compute magnitude. */ WlzBufMagD3(*(iBufA + 0), *(iBufA + 1), *(iBufA + 2), itvLen); /* Clamp into destination interval. */ tGP0.dbp = *(iBufA + 0); WlzValueClampGreyIntoGrey(dstGWSp.u_grintptr, 0, dstGWSp.pixeltype, tGP0, 0, WLZ_GREY_DOUBLE, itvLen); } if(errNum == WLZ_ERR_EOO) { errNum = WLZ_ERR_NONE; } } /* Free intersection objects. */ idN = 0; while(idN < 3) { if(iObjA[idN]) { WlzFreeObj(iObjA[idN]); } ++idN; } /* Free buffers. */ if(iBufA) { Alc2Free((void **)iBufA); } /* Tidy up on error. */ if(dstObj && (errNum != WLZ_ERR_NONE)) { WlzFreeObj(dstObj); dstObj = NULL; } /* Pass back error status. */ if(dstErr) { *dstErr = errNum; } return(dstObj); }