WlzErrorNum ChangeFixedPoint(WlzThreeDViewStruct *wlzViewStr, double xa, double ya, double za ) { WlzErrorNum errNum = WLZ_ERR_NONE; double dstA, dstB, dstC, dstD; double de, ci; /* change the fixed point to (0 0 0) which will change the distance parameter */ Wlz3DViewGetPlaneEqn(wlzViewStr, &dstA, &dstB, &dstC, &dstD); de = dstA * xa + dstB * ya + dstC * za + dstD; ci = sqrt( dstA * dstA + dstB * dstB + dstC * dstC ); /* get the distance */ wlzViewStr->dist = -de/ci; /* wlzViewStr->dist = dstD /sqrt( dstA * dstA + dstB * dstB + dstC * dstC ); */ wlzViewStr->fixed.vtX = xa; wlzViewStr->fixed.vtY = ya; wlzViewStr->fixed.vtZ = za; /* to get 3D sections you need the type */ wlzViewStr->type = WLZ_3D_VIEW_STRUCT; return errNum; }
/*! * \return New object with the rojection. * \ingroup WlzTransform * \brief Use the view transform to define a projection from * 3D to 2D and then project the object onto this plane. * The object supplied to this function must be a 3D * spatial domain object (WLZ_3D_DOMAINOBJ) with either * no values or for integration WLZ_GREY_UBYTE values. * Integration will assign each output pixel the sum of * all input voxels mapped via either the domain density * or the voxel density. * The integration is controled by the integrate parameter * with valid values: * WLZ_PROJECT_INT_MODE_NONE - a "shadow domain" without values * is computed, * WLZ_PROJECT_INT_MODE_DOMAIN - the voxels of the domain are * integrated using * \f[ p = \frac{1}{255} n d \f] * WLZ_PROJECT_INT_MODE_VALUES - the voxel values are integrated * using * \f[ p = \frac{1}{255} \sum{l\left[v\right]}. \f] * Where * \f$p\f$ is the projected image value, * \f$n\f$ is the number of voxels projected for \f$p\f$, * \f$d\f$ is the density of domain voxels, * \f$l\f$ is the voxel value density look up table and * \f$v\f$ is a voxel value. * \param obj The given object. * \param vStr Given view structure defining the * projection plane. * \param intMod This may take three values: * WLZ_PROJECT_INT_MODE_NONE, * WLZ_PROJECT_INT_MODE_DOMAIN or * WLZ_PROJECT_INT_MODE_VALUES. * \param denDom Density of domain voxels this value * is not used unless the integration * mode is WLZ_PROJECT_INT_MODE_DOMAIN. * \param denVal Density look up table for object * voxel density values which must be * an array of 256 values. This may be * NULL if the integration mode is not * WLZ_PROJECT_INT_MODE_VALUES. * \param depth If greater than zero, the projection * depth perpendicular to the viewing * plane. * \param dstErr Destination error pointer, may be NULL. */ WlzObject *WlzProjectObjToPlane(WlzObject *obj, WlzThreeDViewStruct *vStr, WlzProjectIntMode intMod, WlzUByte denDom, WlzUByte *denVal, double depth, WlzErrorNum *dstErr) { int nThr = 1, itvVal = 0; WlzIVertex2 prjSz; WlzIBox2 prjBox = {0}; double pln[4]; WlzObject *bufObj = NULL, *prjObj = NULL; WlzThreeDViewStruct *vStr1 = NULL; double **vMat = NULL; WlzValues nullVal; WlzErrorNum errNum = WLZ_ERR_NONE; WlzAffineTransform *rescaleTr = NULL; WlzGreyValueWSpace **gVWSp = NULL; void ***prjAry = NULL; const double eps = 0.000001; #ifdef WLZ_DEBUG_PROJECT3D_TIME struct timeval times[3]; #endif /* WLZ_DEBUG_PROJECT3D_TIME */ nullVal.core = NULL; if(obj == NULL) { errNum = WLZ_ERR_OBJECT_NULL; } else if(obj->type != WLZ_3D_DOMAINOBJ) { errNum = WLZ_ERR_OBJECT_TYPE; } else if(obj->domain.core == NULL) { errNum = WLZ_ERR_DOMAIN_NULL; } else if(obj->domain.core->type != WLZ_PLANEDOMAIN_DOMAIN) { errNum = WLZ_ERR_DOMAIN_TYPE; } else if(vStr == NULL) { errNum = WLZ_ERR_TRANSFORM_NULL; } else if((intMod == WLZ_PROJECT_INT_MODE_VALUES) && (obj->values.core == NULL)) { errNum = WLZ_ERR_VALUES_NULL; } #ifdef WLZ_DEBUG_PROJECT3D_TIME gettimeofday(times + 0, NULL); #endif /* WLZ_DEBUG_PROJECT3D_TIME */ /* Create new view transform without voxel scaling. The voxel scaling * is done after the projection. */ if(errNum == WLZ_ERR_NONE) { if((vStr1 = WlzMake3DViewStruct(WLZ_3D_VIEW_STRUCT, &errNum)) != NULL) { vStr1->fixed = vStr->fixed; vStr1->theta = vStr->theta; vStr1->phi = vStr->phi; vStr1->zeta = vStr->zeta; vStr1->dist = vStr->dist; vStr1->scale = vStr->scale; vStr1->voxelSize[0] = 1.0; vStr1->voxelSize[1] = 1.0; vStr1->voxelSize[2] = 1.0; vStr1->voxelRescaleFlg = 0; vStr1->interp = vStr->interp; vStr1->view_mode = vStr->view_mode; vStr1->up = vStr->up; vStr1->initialised = WLZ_3DVIEWSTRUCT_INIT_NONE; vMat = vStr1->trans->mat; errNum = WlzInit3DViewStructAffineTransform(vStr1); if(errNum == WLZ_ERR_NONE) { errNum = Wlz3DViewStructTransformBB(obj, vStr1); } if(errNum != WLZ_ERR_NONE) { WlzFree3DViewStruct(vStr1); vStr1 = NULL; } } } /* Compute bounding box of the projection. */ if(errNum == WLZ_ERR_NONE) { prjBox.xMin = WLZ_NINT(vStr1->minvals.vtX); prjBox.yMin = WLZ_NINT(vStr1->minvals.vtY); prjBox.xMax = WLZ_NINT(vStr1->maxvals.vtX); prjBox.yMax = WLZ_NINT(vStr1->maxvals.vtY); prjSz.vtX = prjBox.xMax - prjBox.xMin + 1; prjSz.vtY = prjBox.yMax - prjBox.yMin + 1; } /* Compute post projection scaling. */ if((errNum == WLZ_ERR_NONE) && (vStr->voxelRescaleFlg != 0)) { WlzIBox2 sBox; WlzIVertex2 sSz; WlzThreeDViewStruct *vStr2; vStr2 = WlzMake3DViewStruct(WLZ_3D_VIEW_STRUCT, &errNum); if(errNum == WLZ_ERR_NONE) { vStr2->fixed = vStr->fixed; vStr2->theta = vStr->theta; vStr2->phi = vStr->phi; vStr2->zeta = vStr->zeta; vStr2->dist = vStr->dist; vStr2->scale = vStr->scale; vStr2->voxelSize[0] = vStr->voxelSize[0]; vStr2->voxelSize[1] = vStr->voxelSize[1]; vStr2->voxelSize[2] = vStr->voxelSize[2]; vStr2->voxelRescaleFlg = vStr->voxelRescaleFlg; vStr2->interp = vStr->interp; vStr2->view_mode = vStr->view_mode; vStr2->up = vStr->up; vStr2->initialised = WLZ_3DVIEWSTRUCT_INIT_NONE; errNum = WlzInit3DViewStructAffineTransform(vStr2); if(errNum == WLZ_ERR_NONE) { errNum = Wlz3DViewStructTransformBB(obj, vStr2); } if(errNum == WLZ_ERR_NONE) { sBox.xMin = WLZ_NINT(vStr2->minvals.vtX); sBox.yMin = WLZ_NINT(vStr2->minvals.vtY); sBox.xMax = WLZ_NINT(vStr2->maxvals.vtX); sBox.yMax = WLZ_NINT(vStr2->maxvals.vtY); sSz.vtX = sBox.xMax - sBox.xMin + 1; sSz.vtY = sBox.yMax - sBox.yMin + 1; rescaleTr = WlzMakeAffineTransform(WLZ_TRANSFORM_2D_AFFINE, &errNum); } if(errNum == WLZ_ERR_NONE) { double **m; m = rescaleTr->mat; m[0][0] = (sSz.vtX * eps) / (prjSz.vtX * eps); m[1][1] = (sSz.vtY * eps) / (prjSz.vtY * eps); m[0][2] = sBox.xMin - WLZ_NINT(m[0][0] * prjBox.xMin); m[1][2] = sBox.yMin - WLZ_NINT(m[1][1] * prjBox.yMin); } (void )WlzFree3DViewStruct(vStr2); } } /* Compute plane equation, used to clip intervals if depth was given. */ if((errNum == WLZ_ERR_NONE) && (depth > eps)) { Wlz3DViewGetPlaneEqn(vStr1, pln + 0, pln + 1, pln + 2, pln + 3); } /* Create rectangular projection array buffers, one for each thread, * also if integrating values create a grey value workspace per thread. */ if(errNum == WLZ_ERR_NONE) { int idB; #ifdef _OPENMP #pragma omp parallel { #pragma omp master { nThr = omp_get_num_threads(); } } #endif if((prjAry = (void ***)AlcCalloc(nThr, sizeof(void **))) == NULL) { errNum = WLZ_ERR_MEM_ALLOC; } else { if(intMod == WLZ_PROJECT_INT_MODE_NONE) { for(idB = 0; idB < nThr; ++idB) { if(AlcUnchar2Calloc((WlzUByte ***)&(prjAry[idB]), prjSz.vtY, prjSz.vtX) != ALC_ER_NONE) { errNum = WLZ_ERR_MEM_ALLOC; break; } } } else { for(idB = 0; idB < nThr; ++idB) { if(AlcInt2Calloc((int ***)&(prjAry[idB]), prjSz.vtY, prjSz.vtX) != ALC_ER_NONE) { errNum = WLZ_ERR_MEM_ALLOC; break; } } } } if((errNum == WLZ_ERR_NONE) && (intMod == WLZ_PROJECT_INT_MODE_VALUES)) { itvVal = (WlzGreyTableIsTiled(obj->values.core->type) == 0); if(itvVal == 0) { if((gVWSp = AlcCalloc(nThr, sizeof(WlzGreyValueWSpace *))) == NULL) { errNum = WLZ_ERR_MEM_ALLOC; } else { for(idB = 0; idB < nThr; ++idB) { gVWSp[idB] = WlzGreyValueMakeWSp(obj, &errNum); if(gVWSp[idB]->gType != WLZ_GREY_UBYTE) { errNum = WLZ_ERR_GREY_TYPE; break; } } } } } } /* Scan through the 3D domain setting value in the projection array. */ if(errNum == WLZ_ERR_NONE) { int pIdx, pCnt; WlzDomain *doms; WlzValues *vals = NULL; doms = obj->domain.p->domains; if(itvVal) { vals = obj->values.vox->values; } pCnt = obj->domain.p->lastpl - obj->domain.p->plane1 + 1; #ifdef _OPENMP #pragma omp parallel for #endif for(pIdx = 0; pIdx < pCnt; ++pIdx) { int thrId = 0; if((errNum == WLZ_ERR_NONE) && (doms[pIdx].core != NULL)) { WlzObject *obj2; WlzGreyWSpace gWSp; WlzIntervalWSpace iWSp; WlzErrorNum errNum2 = WLZ_ERR_NONE; #ifdef _OPENMP thrId = omp_get_thread_num(); #endif obj2 = WlzMakeMain(WLZ_2D_DOMAINOBJ, doms[pIdx], (vals)? vals[pIdx]: nullVal, NULL, NULL, &errNum2); if(errNum2 == WLZ_ERR_NONE) { if(itvVal) { errNum2 = WlzInitGreyScan(obj2, &iWSp, &gWSp); } else { errNum2 = WlzInitRasterScan(obj2, &iWSp, WLZ_RASTERDIR_ILIC); } } if(errNum2 == WLZ_ERR_NONE) { double plnZ, vMZX, vMZY; WlzIVertex3 p0, p1; p0.vtZ = p1.vtZ = obj->domain.p->plane1 + pIdx; vMZX = (vMat[0][2] * p0.vtZ) + vMat[0][3] - prjBox.xMin; vMZY = (vMat[1][2] * p0.vtZ) + vMat[1][3] - prjBox.yMin; plnZ = (pln[2] * p0.vtZ) + pln[3]; while(((itvVal == 0) && ((errNum2 = WlzNextInterval(&iWSp)) == WLZ_ERR_NONE)) || ((itvVal != 0) && ((errNum2 = WlzNextGreyInterval(&iWSp)) == WLZ_ERR_NONE))) { int skip = 0; WlzDVertex2 q0, q1; p0.vtX = iWSp.lftpos; p1.vtX = iWSp.rgtpos; p0.vtY = p1.vtY = iWSp.linpos; if(depth > eps) { int c; double d0, d1, plnYZ; /* Clip the 3D line segment p0,p1 using the plane equation. */ plnYZ = (pln[1] * p0.vtY) + plnZ; d0 = (pln[0] * p0.vtX) + plnYZ; d1 = (pln[0] * p1.vtX) + plnYZ; c = ((d1 > depth) << 3) | ((d0 > depth) << 2) | ((d1 < -depth) << 1) | (d0 < -depth); if(c) { if((c == 3) || (c == 12)) /* 00-- or ++00 */ { /* Both out of range, so don't render. */ skip = 1; } else { if(fabs(pln[0]) > eps) { double plnX; plnX = -1.0 / pln[0]; if((c & 1) != 0) /* x0x- */ { p0.vtX = plnX * (plnYZ + depth); } else if((c & 4) != 0) /* x+x0 */ { p0.vtX = plnX * (plnYZ - depth); } if((c & 2) != 0) /* 0x-x */ { p1.vtX = plnX * (plnYZ + depth); } else if((c & 8) != 0) /* +x0x */ { p1.vtX = plnX * (plnYZ - depth); } } } } } if(skip == 0) { q0.vtX = (vMat[0][0] * p0.vtX) + (vMat[0][1] * p0.vtY) + vMZX; q0.vtY = (vMat[1][0] * p0.vtX) + (vMat[1][1] * p0.vtY) + vMZY; q1.vtX = (vMat[0][0] * p1.vtX) + (vMat[0][1] * p1.vtY) + vMZX; q1.vtY = (vMat[1][0] * p1.vtX) + (vMat[1][1] * p1.vtY) + vMZY; switch(intMod) { case WLZ_PROJECT_INT_MODE_NONE: { WlzIVertex2 u0, u1; WLZ_VTX_2_NINT(u0, q0); WLZ_VTX_2_NINT(u1, q1); WlzProjectObjLine((WlzUByte **)(prjAry[thrId]), u0, u1); } break; case WLZ_PROJECT_INT_MODE_DOMAIN: { int np, nq; WlzDVertex3 dq; WlzIVertex2 u0, u1; WLZ_VTX_2_NINT(u0, q0); WLZ_VTX_2_NINT(u1, q1); WLZ_VTX_2_SUB(dq, q0, q1); np = denDom * (iWSp.rgtpos - iWSp.lftpos + 1); nq = (int )ceil(WLZ_VTX_2_LENGTH(dq) + eps); WlzProjectObjLineDom((int **)(prjAry[thrId]), np / nq, u0, u1); } break; case WLZ_PROJECT_INT_MODE_VALUES: if(itvVal) { WlzProjectObjLineVal((int **)(prjAry[thrId]), denVal, gWSp.u_grintptr.ubp, NULL, vMat, vMZX, vMZY, p0, p1); } else { WlzProjectObjLineVal((int **)(prjAry[thrId]), denVal, NULL, gVWSp[thrId], vMat, vMZX, vMZY, p0, p1); } break; } } } (void )WlzEndGreyScan(&iWSp, &gWSp); if(errNum2 == WLZ_ERR_EOO) { errNum2 = WLZ_ERR_NONE; } } (void )WlzFreeObj(obj2); if(errNum2 != WLZ_ERR_NONE) { #ifdef _OPENMP #pragma omp critical { #endif if(errNum == WLZ_ERR_NONE) { errNum = errNum2; } #ifdef _OPENMP } #endif } } } } /* Free grey value workspaces if they were created. */ if(gVWSp) { int idB; for(idB = 0; idB < nThr; ++idB) { WlzGreyValueFreeWSp(gVWSp[idB]); } AlcFree(gVWSp); } if(errNum == WLZ_ERR_NONE) { int idB; size_t idC, bufSz; WlzGreyP buf0, buf1; WlzIVertex2 prjOrg; prjOrg.vtX = prjBox.xMin; prjOrg.vtY = prjBox.yMin; bufSz = prjSz.vtX * prjSz.vtY; for(idB = 1; idB < nThr; ++idB) { if(intMod == WLZ_PROJECT_INT_MODE_NONE) { buf0.ubp = ((WlzUByte ***)(prjAry))[0][0], buf1.ubp = ((WlzUByte ***)(prjAry))[idB][0]; for(idC = 0; idC < bufSz; ++idC) { buf0.ubp[idC] += buf1.ubp[idC]; } } else { buf0.inp = ((int ***)(prjAry))[0][0], buf1.inp = ((int ***)(prjAry))[idB][0]; for(idC = 0; idC < bufSz; ++idC) { buf0.inp[idC] += buf1.inp[idC]; } } } switch(intMod != WLZ_PROJECT_INT_MODE_NONE) { buf0.inp = ((int ***)(prjAry))[0][0]; for(idC = 0; idC < bufSz; ++idC) { buf0.inp[idC] /= 256; } } if(intMod == WLZ_PROJECT_INT_MODE_NONE) { bufObj = WlzAssignObject( WlzFromArray2D((void **)(prjAry[0]), prjSz, prjOrg, WLZ_GREY_UBYTE, WLZ_GREY_UBYTE, 0.0, 1.0, 1, 0, &errNum), NULL); } else { bufObj = WlzAssignObject( WlzFromArray2D((void **)(prjAry[0]), prjSz, prjOrg, WLZ_GREY_INT, WLZ_GREY_INT, 0.0, 1.0, 1, 0, &errNum), NULL); } } /* Free the projection array(s). */ if(prjAry) { int idB; for(idB = 0; idB < nThr; ++idB) { (void )Alc2Free((prjAry[idB])); } AlcFree(prjAry); } /* Make return object using threshold. */ if(errNum == WLZ_ERR_NONE) { WlzPixelV tV; WlzObject *tObj = NULL; tV.type = WLZ_GREY_UBYTE; tV.v.ubv = 1; tObj = WlzAssignObject( WlzThreshold(bufObj, tV, WLZ_THRESH_HIGH, &errNum), NULL); if(tObj) { if(intMod == WLZ_PROJECT_INT_MODE_NONE) { prjObj = WlzMakeMain(tObj->type, tObj->domain, nullVal, NULL, NULL, &errNum); } else { prjObj = WlzMakeMain(tObj->type, tObj->domain, tObj->values, NULL, NULL, &errNum); } } (void )WlzFreeObj(tObj); } (void )WlzFreeObj(bufObj); (void )WlzFree3DViewStruct(vStr1); /* Scale image. */ if(rescaleTr != NULL) { if(errNum == WLZ_ERR_NONE) { WlzObject *tObj = NULL; tObj = WlzAffineTransformObj(prjObj, rescaleTr, WLZ_INTERPOLATION_NEAREST, &errNum); (void )WlzFreeObj(prjObj); prjObj = tObj; } (void )WlzFreeAffineTransform(rescaleTr); } #ifdef WLZ_DEBUG_PROJECT3D_TIME gettimeofday(times + 1, NULL); ALC_TIMERSUB(times + 1, times + 0, times + 2); (void )fprintf(stderr, "WlzGetProjectionFromObject: Elapsed time = %g\n", times[2].tv_sec + (0.000001 * times[2].tv_usec)); #endif /* WLZ_DEBUG_PROJECT3D_TIME */ if(dstErr) { *dstErr = errNum; } return(prjObj); }