calc_dose.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362
  1. /* calc_dose.cpp */
  2. /* The dose at all voxels in the grid dose_mask is calculated using a convolution
  3. method that uses a polyenergetic kernel. Inhomogeneities are accounted for by kernel
  4. scaling, and an inverse square correction is applied after the convolution of the terma
  5. grid with the kernel, rather than being applied directly to the terma grid before
  6. the convolution. */
  7. #include <stdio.h>
  8. #include <stdlib.h>
  9. #include <malloc.h>
  10. #include <math.h>
  11. #include <string.h>
  12. #include "defs.h"
  13. extern char errstr[200]; // error string that all routines have access to
  14. int calc_dose(FLOAT_GRID *density, FLOAT_GRID *terma,FLOAT_GRID *dose, KERNEL *kern,
  15. BEAM *bm, FLOAT_GRID *dose_mask)
  16. {
  17. int M, N, Q; // dimensions of CT array
  18. int baseindex;
  19. float SAD;
  20. float length;
  21. float *ip, *jp, *kp;
  22. float dx, dy, dz;
  23. float delr; // convolution step size
  24. float *x, *y, *z; // vectors of CT coordinates
  25. float *phi, *theta, *sinphi, *cosphi, *sintheta, *costheta;
  26. float one;
  27. float rho;
  28. one = (float)1.0;
  29. // copy CT dimensions and voxel sizes for shorter references later
  30. M = density->x_count;
  31. N = density->y_count;
  32. Q = density->z_count;
  33. dx = density->inc.x;
  34. dy = density->inc.y;
  35. dz = density->inc.z;
  36. // copy vectors describing the beam's eye view coordinate system as well
  37. ip = fvector(0,2); // ip and jp span the beam's eye view
  38. jp = fvector(0,2);
  39. kp = fvector(0,2); // beam direction
  40. // create the unit vector describing the beam direction
  41. for (int j=0;j<3;j++) ip[j] = bm->ip[j];
  42. for (int j=0;j<3;j++) jp[j] = bm->jp[j];
  43. for (int j=0;j<3;j++) kp[j] = bm->kp[j];
  44. // vectors describing the location of each voxel
  45. x = fvector(0,M-1);
  46. y = fvector(0,N-1);
  47. z = fvector(0,Q-1);
  48. // lookup table for vectors in polar coordinates
  49. phi = fvector(0,NPHI-1);
  50. theta = fvector(0,NTHETA-1);
  51. sinphi = fvector(0,NPHI-1);
  52. cosphi = fvector(0,NPHI-1);
  53. sintheta = fvector(0,NTHETA-1);
  54. costheta = fvector(0,NTHETA-1);
  55. //kernel with fewer elements for faster calculation
  56. //see defs.h
  57. KERNEL smallkern;
  58. smallkern.radial_boundary = (float *)malloc(kern->nradii*sizeof(float));
  59. smallkern.angular_boundary = (float *)malloc(kern->ntheta*sizeof(float));
  60. //small kernel dimensions
  61. smallkern.nradii = N_KERNEL_RADII; //same as original kernel
  62. smallkern.ntheta = NTHETA; //must divide evenly into N_KERNEL_ANGLES
  63. SAD = bm->SAD;
  64. for (int i=0;i<N_KERNEL_CATEGORIES;i++)
  65. if ( (smallkern.matrix[i] =
  66. (float *) calloc(smallkern.ntheta*smallkern.nradii,sizeof(float))) == NULL)
  67. {
  68. sprintf(errstr,"Cannot allocate space for matrix %d\n",i);
  69. return(FAILURE);
  70. }
  71. if ( (smallkern.total_matrix =
  72. (float *) calloc(smallkern.ntheta*smallkern.nradii,sizeof(float))) == NULL)
  73. {
  74. sprintf(errstr,"Cannot allocate space for total matrix\n");
  75. return(FAILURE);
  76. }
  77. //set up boundaries
  78. for (int i=0;i<smallkern.ntheta;i++)
  79. smallkern.angular_boundary[i] = ( (float) i + 1) * (float)180.0/(float) smallkern.ntheta;
  80. for (int i=0;i<smallkern.nradii;i++)
  81. smallkern.radial_boundary[i] = kern->radial_boundary[i];
  82. //initialise
  83. for (int i=0;i<smallkern.nradii;i++)
  84. for (int j=0;j<smallkern.ntheta;j++)
  85. {
  86. KERNEL_TOTAL_VALUE(&smallkern,i,j) = (float)0.0;
  87. for (int l=0;l<N_KERNEL_CATEGORIES;l++)
  88. KERNEL_VALUE(&smallkern,l,i,j) = (float)0.0;
  89. }
  90. //create kernel values
  91. for (int i=0;i<smallkern.nradii;i++)
  92. for (int j=0;j<smallkern.ntheta;j++)
  93. {
  94. //first angular index in original kernel for this element
  95. baseindex = j*N_KERNEL_ANGLES/NTHETA;
  96. //for each category, sum values from original kernel
  97. for (int l=0;l<N_KERNEL_CATEGORIES;l++)
  98. for (int k=0;k<N_KERNEL_ANGLES/NTHETA;k++)
  99. KERNEL_VALUE(&smallkern,l,i,j) += KERNEL_VALUE(kern,l,i,baseindex+k);
  100. //and for total kernel
  101. for (int k=0;k<N_KERNEL_ANGLES/NTHETA;k++)
  102. KERNEL_TOTAL_VALUE(&smallkern,i,j) += KERNEL_TOTAL_VALUE(kern,i,baseindex+k);
  103. }
  104. //Make cumulative kernel (with radius)
  105. //this is what is used for the dose calculation
  106. for (int p=0;p<smallkern.ntheta;p++)
  107. for (int r=0;r<smallkern.nradii;r++)
  108. {
  109. for (int i=0;i<N_KERNEL_CATEGORIES;i++)
  110. if (r > 0)
  111. KERNEL_VALUE(&smallkern,i,r,p)
  112. = KERNEL_VALUE(&smallkern,i,r-1,p) + KERNEL_VALUE(&smallkern,i,r,p);
  113. if (r > 0)
  114. KERNEL_TOTAL_VALUE(&smallkern,r,p)
  115. = KERNEL_TOTAL_VALUE(&smallkern,r-1,p) + KERNEL_TOTAL_VALUE(&smallkern,r,p);
  116. }
  117. // fill the coordinate vectors
  118. for (int i=0;i<M;i++) x[i] = density->start.x + (float)i*dx;
  119. for (int j=0;j<N;j++) y[j] = density->start.y + (float)j*dy;
  120. for (int k=0;k<Q;k++) z[k] = density->start.z + (float)k*dz;
  121. // fill in the polar coordinates vectors
  122. for (int q=0;q<NPHI;q++)
  123. {
  124. phi[q] = ((float)q + (float)0.5)*(float)2.0*(float)PI/(float)NPHI;
  125. sinphi[q] = (float)sin(phi[q]);
  126. cosphi[q] = (float)cos(phi[q]);
  127. }
  128. // Theta is subtracted from PI is because direction of travel along kernel ray is actually opposite of
  129. // direction along which energy is radiating, so need to use a source voxel direction that
  130. // is reflected about horizontal. This can be thought of as the kernel inversion line.
  131. for (int p=0;p<smallkern.ntheta;p++)
  132. if (p == 0)
  133. {
  134. theta[p] = (float)0.5*smallkern.angular_boundary[0]*(float)PI/(float)180.0;
  135. sintheta[p] = (float)sin((float)PI - theta[p]);
  136. costheta[p] = (float)cos((float)PI - theta[p]);
  137. }
  138. else
  139. {
  140. theta[p] = (float)0.5*(smallkern.angular_boundary[p-1] + smallkern.angular_boundary[p])*(float)PI/(float)180.0;
  141. sintheta[p] = (float)sin((float)PI - theta[p]);
  142. costheta[p] = (float)cos((float)PI - theta[p]);
  143. }
  144. // store the sines and cosines in a lookup table
  145. // the step size for the convolution integration is the smallest voxel side length
  146. if (dx <= dy && dx <= dz)
  147. delr = (float)2.0*dx;
  148. else if (dy <= dx && dy <= dz)
  149. delr = (float)2.0*dy;
  150. else
  151. delr = (float)2.0*dz;
  152. //calculate dose at each point
  153. //done from deposition (catcher's) point of view
  154. #pragma omp parallel for
  155. for (int k=0;k<Q; k++)
  156. for (int j=0;j<N; j++)
  157. for (int i=0;i<M; i++)
  158. if (GRID_VALUE(dose_mask,i,j,k) > 0) // only calculate dose inside dose mask
  159. {
  160. int i, j, k, l;
  161. int p, q, r;
  162. float del_x, del_y, del_z;
  163. float current_x, current_y, current_z;
  164. float cumval, last_cumval;
  165. int inter_i, inter_j, inter_k;
  166. float inter_x, inter_y, inter_z;
  167. float r_eff, delr_eff, inter_r_eff;
  168. int r_index;
  169. float t, u, v, f;
  170. float kval;
  171. // do the integral for the point in the ROI
  172. for (p=0;p<smallkern.ntheta;p++) //polar
  173. for (q=0;q<NPHI;q++) //azimuthal
  174. {
  175. //initialise position of current voxel
  176. current_x = x[i];
  177. current_y = y[j];
  178. current_z = z[k];
  179. //initialise effective radius along kernel direction
  180. r_eff = 0.0;
  181. //initialise cumulative kernel value for this direction
  182. last_cumval = 0.0;
  183. //Using reciprocity technique where dose at point A due to point B
  184. //is dose at point B due to point A
  185. //x ,y, z increments along ray
  186. del_x = delr*(ip[0]*cosphi[q]*sintheta[p] + jp[0]*sinphi[q]*sintheta[p]
  187. + kp[0]*costheta[p]);
  188. del_y = delr*(ip[1]*cosphi[q]*sintheta[p] + jp[1]*sinphi[q]*sintheta[p]
  189. + kp[1]*costheta[p]);
  190. del_z = delr*(ip[2]*cosphi[q]*sintheta[p] + jp[2]*sinphi[q]*sintheta[p]
  191. + kp[2]*costheta[p]);
  192. //initialise physical radius
  193. r = 0;
  194. do
  195. {
  196. //interaction point is at mid-point of curent increment
  197. inter_x = current_x + (float)0.5*del_x;
  198. inter_y = current_y + (float)0.5*del_y;
  199. inter_z = current_z + (float)0.5*del_z;
  200. //voxel containing interaction point
  201. inter_i = (int) ((inter_x - density->start.x)/dx);
  202. inter_j = (int) ((inter_y - density->start.y)/dy);
  203. inter_k = (int) ((inter_z - density->start.z)/dz);
  204. // stop the integral if interaction point is outside the dose calculation limits
  205. if ( (inter_i < 0) || (inter_i + 1 >= M)
  206. || (inter_j < 0) || (inter_j + 1 >= N)
  207. || (inter_k < 0) || (inter_k + 1 >= Q)
  208. || (GRID_VALUE(dose_mask,inter_i,inter_j,inter_k) <= 0.0))
  209. break;
  210. // Position of the end of the increment. Interaction point is at the
  211. // midpoint.
  212. current_x += del_x;
  213. current_y += del_y;
  214. current_z += del_z;
  215. //effective distance increment
  216. delr_eff = delr*GRID_VALUE(density,inter_i,inter_j,inter_k);
  217. //effective radius of interaction point
  218. inter_r_eff = r_eff + (float)0.5*delr_eff;
  219. r_eff += delr_eff;
  220. // trilinear interpolation method of the terma contribution, f
  221. // relative differences between the interaction point and the lower voxel bound
  222. t = (inter_x - x[inter_i])/dx;
  223. u = (inter_y - y[inter_j])/dy;
  224. v = (inter_z - z[inter_k])/dz;
  225. f = GRID_VALUE(terma,inter_i,inter_j,inter_k)*(one-t)*(one-u)*(one-v)
  226. + GRID_VALUE(terma,inter_i,inter_j,inter_k+1)*(one-t)*(one-u)*v
  227. + GRID_VALUE(terma,inter_i,inter_j+1,inter_k+1)*(one-t)*u*v
  228. + GRID_VALUE(terma,inter_i+1,inter_j+1,inter_k+1)*t*u*v
  229. + GRID_VALUE(terma,inter_i,inter_j+1,inter_k)*(one-t)*u*(one-v)
  230. + GRID_VALUE(terma,inter_i+1,inter_j+1,inter_k)*t*u*(one-v)
  231. + GRID_VALUE(terma,inter_i+1,inter_j,inter_k+1)*t*(one-u)*v
  232. + GRID_VALUE(terma,inter_i+1,inter_j,inter_k)*t*(one-u)*(one-v);
  233. /*
  234. // interpolate density at the interaction point
  235. rho = GRID_VALUE(density,inter_i,inter_j,inter_k)*(one-t)*(one-u)*(one-v)
  236. + GRID_VALUE(density,inter_i,inter_j,inter_k+1)*(one-t)*(one-u)*v
  237. + GRID_VALUE(density,inter_i,inter_j+1,inter_k+1)*(one-t)*u*v
  238. + GRID_VALUE(density,inter_i+1,inter_j+1,inter_k+1)*t*u*v
  239. + GRID_VALUE(density,inter_i,inter_j+1,inter_k)*(one-t)*u*(one-v)
  240. + GRID_VALUE(density,inter_i+1,inter_j+1,inter_k)*t*u*(one-v)
  241. + GRID_VALUE(density,inter_i+1,inter_j,inter_k+1)*t*(one-u)*v
  242. + GRID_VALUE(density,inter_i+1,inter_j,inter_k)*t*(one-u)*(one-v); */
  243. // perform kernel lookup for r_eff, r_index is the kernel index of the first
  244. // bin boundary below the effective radius of the voxel
  245. r_index = binSearch(smallkern.radial_boundary,inter_r_eff,smallkern.nradii);
  246. // interpolate to obtain the effective cumulative kernel value
  247. if (r_index == -1) // radius is between zero and the first bin boundary
  248. {
  249. // fractional difference between inter_r_eff and the first bin boundary
  250. t = inter_r_eff/smallkern.radial_boundary[0];
  251. cumval = (1-t)*KERNEL_TOTAL_VALUE(&smallkern,0,p);
  252. }
  253. else if (r_index >= smallkern.nradii-1) // overshot the kernel bin boundaries
  254. {
  255. cumval = KERNEL_TOTAL_VALUE(&smallkern,smallkern.nradii-1,p);
  256. }
  257. else // inter_r_eff is between the first upper bin boundary and the last
  258. {
  259. t = (inter_r_eff - smallkern.radial_boundary[r_index])
  260. /(smallkern.radial_boundary[r_index + 1] - smallkern.radial_boundary[r_index]);
  261. cumval = (1-t)*KERNEL_TOTAL_VALUE(&smallkern,r_index,p)
  262. + t*KERNEL_TOTAL_VALUE(&smallkern,r_index+1,p);
  263. }
  264. kval = cumval - last_cumval;
  265. last_cumval = cumval;
  266. // Kernel value to use is current increment in cumulative value
  267. // Note that this is the fractional dose deposited at i,j,k due to
  268. // terma in an effective increment (increment*density) along the kernel ray at the
  269. // interaction point. The value comes from the fractional dose deposited in an
  270. // effective increment along the ray at i,j,k due to terma at the current
  271. // interaction point.
  272. //Increment dose value.
  273. //Note that to include the effect of density at i,j,k on ENERGY deposited at i,j,k,
  274. //there should be a density term, but this cancels on conversion of energy to dose as indicated below
  275. // if (GRID_VALUE(terma,inter_i,inter_j,inter_k) > 0.001)
  276. GRID_VALUE(dose,i,j,k) += f*kval;
  277. r++;
  278. }
  279. while (r<10000);
  280. } //p,q
  281. GRID_VALUE(dose,i,j,k)/= NPHI;
  282. }
  283. //Inverse square correction to dose
  284. //This works better than applying the inverse square correction to terma
  285. //See Papanikolaou and Mackie 1993
  286. for (int k=0;k<Q;k++)
  287. for (int j=0;j<N;j++)
  288. for (int i=0;i<M;i++)
  289. if (GRID_VALUE(dose,i,j,k) > 0.0)
  290. {
  291. // squared difference between the source and the voxel
  292. length = (float)pow(x[i] - bm->y_vec[0],2.0f) + (float)pow(y[j] - bm->y_vec[1],2.0f)
  293. + (float)pow(z[k] - bm->y_vec[2],2.0f);
  294. if (length > 0.0)
  295. GRID_VALUE(dose,i,j,k) *= SAD*SAD/length;
  296. }
  297. // free vectors
  298. free_fvector(ip,0,2);
  299. free_fvector(jp,0,2);
  300. free_fvector(kp,0,2);
  301. free_fvector(x,0,M-1);
  302. free_fvector(y,0,N-1);
  303. free_fvector(z,0,Q-1);
  304. free(smallkern.angular_boundary);
  305. free(smallkern.radial_boundary);
  306. // free(poly_kernel.total_matrix);
  307. for (int j=0;j<N_KERNEL_CATEGORIES;j++)
  308. free(smallkern.matrix[j]);
  309. return(SUCCESS);
  310. }