projectionDD_mex.cpp 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455
  1. /*
  2. %% Author: Rodrigo de Barros Vimieiro
  3. % Date: September, 2018
  4. % rodrigo.vimieiro@gmail.com
  5. % =========================================================================
  6. %{
  7. % -------------------------------------------------------------------------
  8. % projectionDD_mex(data3d,param)
  9. % -------------------------------------------------------------------------
  10. % DESCRIPTION:
  11. % This function calculates the volume projection based on the
  12. % Distance-Driven principle. It works by calculating the overlap in X
  13. % and Y axis of the volume and the detector boundaries.
  14. % The geometry is for DBT with half cone-beam. All parameters are set
  15. % in "ParameterSettings" code.
  16. %
  17. % INPUT:
  18. %
  19. % - data3d = 3D volume for projection
  20. % - param = Parameter of all geometry
  21. %
  22. % OUTPUT:
  23. %
  24. % - proj = projections for each angle.
  25. %
  26. % Reference: Three-Dimensional Digital Tomosynthesis - Yulia
  27. % Levakhina (2014), Cap 3.6 and 3.7.
  28. %
  29. % Original Paper: De Man, Bruno, and Samit Basu. "Distance-driven
  30. % projection and backprojection in three dimensions." Physics in
  31. % Medicine & Biology (2004).
  32. %
  33. % ---------------------------------------------------------------------
  34. % Copyright (C) <2018> <Rodrigo de Barros Vimieiro>
  35. %
  36. % This program is free software: you can redistribute it and/or modify
  37. % it under the terms of the GNU General Public License as published by
  38. % the Free Software Foundation, either version 3 of the License, or
  39. % (at your option) any later version.
  40. %
  41. % This program is distributed in the hope that it will be useful,
  42. % but WITHOUT ANY WARRANTY; without even the implied warranty of
  43. % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  44. % GNU General Public License for more details.
  45. %
  46. % You should have received a copy of the GNU General Public License
  47. % along with this program. If not, see <http://www.gnu.org/licenses/>.
  48. %}
  49. % =========================================================================
  50. %% 3-D Distance Driven Projection Code
  51. */
  52. #include "projectionDD_mex.h"
  53. /****************************************************************************
  54. * function: mexFunction() - Matlab interface function. *
  55. * INPUTS: *
  56. * nlhs - Number of expected output mxArrays, specified as an integer. *
  57. * plhs[] - Array of pointers to the expected mxArray output arguments. *
  58. * nrhs - Number of input mxArrays, specified as an integer. *
  59. * prhs[] - Array of pointers to the mxArray input arguments. *
  60. ****************************************************************************/
  61. void mexFunction(int nlhs, mxArray* plhs[],
  62. int nrhs, const mxArray* prhs[]) {
  63. /* Check for proper number of arguments */
  64. if (nrhs != 2) {
  65. mexErrMsgIdAndTxt("MATLAB:mexcpp:nargin",
  66. "projection_mex requires two input arguments.");
  67. }
  68. if (nlhs != 1) {
  69. mexErrMsgIdAndTxt("MATLAB:mexcpp:nargout",
  70. "projection_mex requires one output argument.");
  71. }
  72. /* Check if the input is of proper type */
  73. if (!mxIsSingle(prhs[0])) {
  74. mexErrMsgIdAndTxt("MATLAB:mexcpp:typeargin",
  75. "First argument has to be single.");
  76. }
  77. if (!mxIsStruct(prhs[1])) {
  78. mexErrMsgIdAndTxt("MATLAB:mexcpp:typeargin",
  79. "Second argument has to be a struct.");
  80. }
  81. float* pVolume = (float*)mxGetPr(prhs[0]);
  82. double* pTubeAngleD = mxGetPr(mxGetField(prhs[1], 0, "tubeDeg"));
  83. double* pDetAngleD = mxGetPr(mxGetField(prhs[1], 0, "detectorDeg"));
  84. const int nProj = (const int)mxGetScalar(mxGetField(prhs[1], 0, "nProj"));
  85. const int nPixX = (const int)mxGetScalar(mxGetField(prhs[1], 0, "nx"));
  86. const int nPixY = (const int)mxGetScalar(mxGetField(prhs[1], 0, "ny"));
  87. const int nSlices = (const int)mxGetScalar(mxGetField(prhs[1], 0, "nz"));
  88. const int nDetX = (const int)mxGetScalar(mxGetField(prhs[1], 0, "nu"));
  89. const int nDetY = (const int)mxGetScalar(mxGetField(prhs[1], 0, "nv"));
  90. const float dx = (const float)mxGetScalar(mxGetField(prhs[1], 0, "dx"));
  91. const float dy = (const float)mxGetScalar(mxGetField(prhs[1], 0, "dy"));
  92. const float dz = (const float)mxGetScalar(mxGetField(prhs[1], 0, "dz"));
  93. const float du = (const float)mxGetScalar(mxGetField(prhs[1], 0, "du"));
  94. const float dv = (const float)mxGetScalar(mxGetField(prhs[1], 0, "dv"));
  95. const float DSD = (const float)mxGetScalar(mxGetField(prhs[1], 0, "DSD"));
  96. const float DDR = (const float)mxGetScalar(mxGetField(prhs[1], 0, "DDR"));
  97. const float DAG = (const float)mxGetScalar(mxGetField(prhs[1], 0, "DAG"));
  98. mwSize NRow = mxGetM(prhs[0]);
  99. mwSize NCol = mxGetN(prhs[0]);
  100. mwSize NElemVol = mxGetNumberOfElements(prhs[0]);
  101. /* Check if the input is in proper size */
  102. if (NRow != nPixY) {
  103. mexErrMsgIdAndTxt("MATLAB:mexcpp:typeargin",
  104. "First argument needs to have the same number of rows as in the configuration file.");
  105. }
  106. if (NCol / nSlices != nPixX) {
  107. mexErrMsgIdAndTxt("MATLAB:mexcpp:typeargin",
  108. "First argument needs to have the same number of columns as in the configuration file.");
  109. }
  110. // Cast double vectors to single
  111. float* const pTubeAngle = (float*)mxMalloc(nProj * sizeof(float));
  112. float* const pDetAngle = (float*)mxMalloc(nProj * sizeof(float));
  113. for (int n = 0; n < nProj; n++) {
  114. pTubeAngle[n] = (float)pTubeAngleD[n];
  115. pDetAngle[n] = (float)pDetAngleD[n];
  116. }
  117. //mexPrintf("Nx:%d Ny:%d Nz:%d \nNu:%d Nv:%d \nDx:%.2f Dy:%.2f Dz:%.2f \nDu:%.2f Dv:%.2f \nDSD:%.2f DDR:%.2f \n", nPixX, nPixY, nSlices, nDetX, nDetY, dx, dy, dz, du, dv, DSD, DDR);
  118. const int nDetXMap = nDetX + 1;
  119. const int nDetYMap = nDetY + 1;
  120. const int nPixXMap = nPixX + 1;
  121. const int nPixYMap = nPixY + 1;
  122. // Memory allocation for projections coordinates
  123. float* const pDetX = (float*)mxMalloc(nDetXMap * sizeof(float));
  124. float* const pDetY = (float*)mxMalloc(nDetYMap * sizeof(float));
  125. float* const pDetZ = (float*)mxMalloc(nDetYMap * sizeof(float));
  126. float* const pObjX = (float*)mxMalloc(nPixXMap * sizeof(float));
  127. float* const pObjY = (float*)mxMalloc(nPixYMap * sizeof(float));
  128. float* const pObjZ = (float*)mxMalloc(nSlices * sizeof(float));
  129. // Memory allocation for mapped coordinates
  130. float* const pDetmY = (float*)mxMalloc(nDetYMap * sizeof(float));
  131. float* const pDetmX = (float*)mxMalloc(nDetYMap * nDetXMap * sizeof(float));
  132. float* const pPixmX = (float*)mxMalloc(nPixXMap * sizeof(float));
  133. float* const pPixmY = (float*)mxMalloc(nPixYMap * sizeof(float));
  134. // Memory allocation for rotated detecto coords
  135. float* const pRdetY = (float*)mxMalloc(nDetYMap * sizeof(float));
  136. float* const pRdetZ = (float*)mxMalloc(nDetYMap * sizeof(float));
  137. // Memory allocation for slice
  138. float* const pSlice = (float*)mxMalloc(nPixX * nPixY * sizeof(float));
  139. // Map detector and object boudaries
  140. for (int k = nDetX, ind = 0; k >= 0; k--, ind++) pDetX[ind] = k * du;
  141. for (int k = -nDetY / 2, ind = 0; k <= nDetY / 2; k++, ind++) pDetY[ind] = k * dv;
  142. for (int k = 0; k <= nDetY; k++) pDetZ[k] = 0;
  143. for (int k = nPixX, ind = 0; k >= 0; k--, ind++) pObjX[ind] = k * dx;
  144. for (int k = -nPixY / 2, ind = 0; k <= nPixY / 2; k++, ind++) pObjY[ind] = k * dy;
  145. for (int k = 0, ind = 0; k <= nSlices; k++, ind++) pObjZ[ind] = k * dz + DAG + dz / 2;
  146. // X - ray tube initial position
  147. float tubeX = 0;
  148. float tubeY = 0;
  149. float tubeZ = DSD;
  150. // Iso - center position
  151. float isoY = 0;
  152. float isoZ = DDR;
  153. // Allocate memory for temp projection variable
  154. float* const pProjt = (float*)mxMalloc(nDetY * nDetX * nProj * sizeof(float));
  155. // Initiate temp proj values with zeros
  156. for (int p = 0; p < nProj; p++)
  157. for (int u = 0; u < nDetX; u++)
  158. for (int v = 0; v < nDetY; v++)
  159. pProjt[(p * nDetY * nDetX) + (u * nDetY) + v] = 0;
  160. // For each projection
  161. for (int p = 0; p < nProj; p++) {
  162. // Get specif tube angle for the projection
  163. float theta = pTubeAngle[p] * M_PI / 180.f;
  164. // Get specif detector angle for the projection
  165. float phi = pDetAngle[p] * M_PI / 180.f;
  166. // Tubre rotation
  167. float rtubeY = ((tubeY - isoY) * (float)cos(theta) - (tubeZ - isoZ) * (float)sin(theta)) + isoY;
  168. float rtubeZ = ((tubeY - isoY) * (float)sin(theta) + (tubeZ - isoZ) * (float)cos(theta)) + isoZ;
  169. // Detector rotation
  170. for (int v = 0; v < nDetYMap; v++) {
  171. pRdetY[v] = ((pDetY[v] - isoY) * (float)cos(phi) - (pDetZ[v] - isoZ) * (float)sin(phi)) + isoY;
  172. pRdetZ[v] = ((pDetY[v] - isoY) * (float)sin(phi) + (pDetZ[v] - isoZ) * (float)cos(phi)) + isoZ;
  173. }
  174. // Map detector onto XY plane(Inside proj loop in case detector rotates)
  175. mapp2xy(pDetmX, pDetmY, tubeX, rtubeY, rtubeZ, pDetX, pRdetY, pRdetZ, nDetXMap, nDetYMap);
  176. // Pixel start index and increment
  177. int detIstart = 0;
  178. int detIinc = 1;
  179. // Mapped detector length
  180. float deltaDetmY = pDetmY[detIstart + detIinc] - pDetmY[detIstart];
  181. // For each slice
  182. for (int z = 0; z < nSlices; z++) {
  183. // Flip X(Img coord is reverse to Global)
  184. for (int x = 0, x_inv = nPixX - 1; x < nPixX; x++, x_inv--)
  185. for (int y = 0; y < nPixY; y++)
  186. pSlice[(x * nPixY) + y] = pVolume[(z * nPixX * nPixY) + (x_inv * nPixY) + y];
  187. // Tmp Z coords value for Y direction
  188. float* const pObjZt = (float*)mxMalloc(nPixYMap * sizeof(float));
  189. // Tmp Pixel X mapped coords
  190. float* const pPixmXt = (float*)mxMalloc(nPixYMap * nPixXMap * sizeof(float));
  191. // Get specif Z coord value for each slice
  192. for (int k = 0; k < nPixYMap; k++) pObjZt[k] = pObjZ[z];
  193. // Map slice onto XY plane
  194. mapp2xy(pPixmXt, pPixmY, tubeX, rtubeY, rtubeZ, pObjX, pObjY, pObjZt, nPixXMap, nPixYMap);
  195. // Flip X(Img coord is reverse to Global)
  196. for (int x = 0, x_inv = nPixXMap - 1; x < nPixXMap; x++, x_inv--)
  197. pPixmX[x] = pPixmXt[x_inv * nPixYMap];
  198. // Free temp variables
  199. mxFree(pObjZt);
  200. mxFree(pPixmXt);
  201. // Pixel start index and increment
  202. int pixIstart = 0;
  203. int pixIinc = 1;
  204. // Mapped pixel length
  205. float deltaPixmX = pPixmX[pixIstart + pixIinc] - pPixmX[pixIstart];
  206. float deltaPixmY = pPixmY[pixIstart + pixIinc] - pPixmY[pixIstart];
  207. // Start pixel and detector indices
  208. int detIndY = detIstart;
  209. int pixIndY = pixIstart;
  210. // Case 1
  211. // Find first detector overlap maped with pixel maped on Y
  212. if (pDetmY[detIndY] - pPixmY[pixIstart] < -deltaDetmY)
  213. while (pDetmY[detIndY] - pPixmY[pixIstart] < -deltaDetmY)
  214. detIndY = detIndY + detIinc;
  215. else
  216. // Case 2
  217. // Find first pixel overlap maped with detector maped on Y
  218. if (pDetmY[detIstart] - pPixmY[pixIndY] > deltaPixmY)
  219. while (pDetmY[detIstart] - pPixmY[pixIndY] > deltaPixmY)
  220. pixIndY = pixIndY + pixIinc;
  221. float moving_left_boundaryY = NULL;
  222. // Get the left coordinate of the first overlap on Y axis
  223. if (pDetmY[detIndY] < pPixmY[pixIndY])
  224. moving_left_boundaryY = pPixmY[pixIndY];
  225. else
  226. moving_left_boundaryY = pDetmY[detIndY];
  227. // Allocate memory for specif row of X map detector coords
  228. float* const pDetmXrow = (float*)mxMalloc(nDetXMap * sizeof(float));
  229. float overLapY = NULL;
  230. // Loop over Y intersections
  231. while ((detIndY < nDetY) && (pixIndY < nPixY)) {
  232. float alpha = (float)atan((pDetmY[detIndY] + (deltaDetmY / 2) - rtubeY) / rtubeZ);
  233. // Case A, when you jump to the next detector boundarie but stay
  234. // in the same pixel
  235. if (pDetmY[detIndY + 1] <= pPixmY[pixIndY + 1])
  236. overLapY = (pDetmY[detIndY + 1] - moving_left_boundaryY) / deltaDetmY; // Normalized overlap Calculation
  237. else
  238. // Case B, when you jump to the next pixel boundarie but stay
  239. // in the same detector
  240. overLapY = (pPixmY[pixIndY + 1] - moving_left_boundaryY) / deltaDetmY; // Normalized overlap Calculation
  241. // ***** X overlap *****
  242. int detIndX = detIstart;
  243. int pixIndX = pixIstart;
  244. // Get row / coll of X flipped, which correspond to that Y overlap det
  245. for (int x = 0, x_inv = nDetXMap - 1; x < nDetXMap; x++, x_inv--)
  246. pDetmXrow[x] = pDetmX[(x_inv * nDetYMap) + detIndY];
  247. // Mapped detecor length on X
  248. float deltaDetmX = pDetmXrow[detIstart + detIinc] - pDetmXrow[detIstart];
  249. // Case 1
  250. // Find first detector overlap maped with pixel maped on X
  251. if (pDetmXrow[detIndX] - pPixmX[pixIstart] < -deltaDetmX)
  252. while (pDetmXrow[detIndX] - pPixmX[pixIstart] < -deltaDetmX)
  253. detIndX = detIndX + detIinc;
  254. else
  255. // Case 2
  256. // Find first pixel overlap maped with detector maped on X
  257. if (pDetmXrow[detIstart] - pPixmX[pixIndX] > deltaPixmX)
  258. while (pDetmXrow[detIstart] - pPixmX[pixIndY] > deltaPixmX)
  259. pixIndX = pixIndX + pixIinc;
  260. float moving_left_boundaryX = NULL;
  261. // Get the left coordinate of the first overlap on X axis
  262. if (pDetmXrow[detIndX] < pPixmX[pixIndX])
  263. moving_left_boundaryX = pPixmX[pixIndX];
  264. else
  265. moving_left_boundaryX = pDetmXrow[detIndX];
  266. // Loop over X intersections
  267. while ((detIndX < nDetX) && (pixIndX < nPixX)) {
  268. float gamma = (float)atan((pDetmXrow[detIndX] + (deltaDetmX / 2) - tubeX) / rtubeZ);
  269. // Case A, when you jump to the next detector boundarie but stay
  270. // in the same pixel
  271. if (pDetmXrow[detIndX + 1] <= pPixmX[pixIndX + 1]) {
  272. float overLapX = (pDetmXrow[detIndX + 1] - moving_left_boundaryX) / deltaDetmX; // Normalized overlap Calculation
  273. pProjt[(p * nDetY * nDetX) + (detIndX * nDetY) + detIndY] += overLapX * overLapY * pSlice[pixIndX * nPixY + pixIndY] * dz / ((float)cos(alpha) * (float)cos(gamma));
  274. detIndX = detIndX + detIinc;
  275. moving_left_boundaryX = pDetmXrow[detIndX];
  276. }
  277. else {
  278. // Case B, when you jump to the next pixel boundarie but stay
  279. // in the same detector
  280. float overLapX = (pPixmX[pixIndX + 1] - moving_left_boundaryX) / deltaDetmX; // Normalized overlap Calculation
  281. pProjt[(p * nDetY * nDetX) + (detIndX * nDetY) + detIndY] += overLapX * overLapY * pSlice[pixIndX * nPixY + pixIndY] * dz / ((float)cos(alpha) * (float)cos(gamma));
  282. pixIndX = pixIndX + pixIinc;
  283. moving_left_boundaryX = pPixmX[pixIndX];
  284. }
  285. }
  286. // ***** Back to Y overlap *****
  287. // Case A, when you jump to the next detector boundarie but stay
  288. // in the same pixel
  289. if (pDetmY[detIndY + 1] <= pPixmY[pixIndY + 1]) {
  290. detIndY = detIndY + detIinc;
  291. moving_left_boundaryY = pDetmY[detIndY];
  292. }
  293. else {
  294. // Case B, when you jump to the next pixel boundarie but stay
  295. // in the same detector
  296. pixIndY = pixIndY + pixIinc;
  297. moving_left_boundaryY = pPixmY[pixIndY];
  298. }
  299. } // Y Overlap loop
  300. // Free memory
  301. mxFree(pDetmXrow);
  302. } // Loop end slices
  303. } // Loop end Projections
  304. // Free memory
  305. mxFree(pSlice);
  306. mxFree(pDetX);
  307. mxFree(pDetY);
  308. mxFree(pDetZ);
  309. mxFree(pObjX);
  310. mxFree(pObjY);
  311. mxFree(pObjZ);
  312. mxFree(pDetmY);
  313. mxFree(pDetmX);
  314. mxFree(pPixmX);
  315. mxFree(pPixmY);
  316. mxFree(pRdetY);
  317. mxFree(pRdetZ);
  318. // Output memory allocation
  319. mwSize dims[3] = { nDetY,nDetX,nProj };
  320. plhs[0] = mxCreateNumericArray(3, dims, mxSINGLE_CLASS, mxREAL);
  321. float* pProj = (float*)mxGetData(plhs[0]);
  322. // Flip projection back X(Img coord is reverse to Global)
  323. for (int p = 0; p < nProj; p++)
  324. for (int x = 0, x_inv = nDetX - 1; x < nDetX; x++, x_inv--)
  325. for (int y = 0; y < nDetY; y++)
  326. pProj[(p * nDetX * nDetY) + (x * nDetY) + y] = pProjt[(p * nDetX * nDetY) + (x_inv * nDetY) + y];
  327. // Free memory
  328. mxFree(pProjt);
  329. return;
  330. }
  331. // Map on XY plane
  332. void mapp2xy(float* const pXmapp,
  333. float* const pYmapp,
  334. float tubeX,
  335. float tubeY,
  336. float tubeZ,
  337. float* const pXcoord,
  338. float* const pYcoord,
  339. float* const pZcoord,
  340. const int nXelem,
  341. const int nYelem) {
  342. for (int x = 0; x < nXelem; x++)
  343. for (int y = 0; y < nYelem; y++) {
  344. int ind = (x * nYelem) + y;
  345. pXmapp[ind] = pXcoord[x] + pZcoord[y] * (pXcoord[x] - tubeX) / (tubeZ - pZcoord[y]);
  346. if (x == 0)
  347. pYmapp[y] = pYcoord[y] + pZcoord[y] * (pYcoord[y] - tubeY) / (tubeZ - pZcoord[y]);
  348. }
  349. return;
  350. }