Scheduled service maintenance on November 22


On Friday, November 22, 2024, between 06:00 CET and 18:00 CET, GIN services will undergo planned maintenance. Extended service interruptions should be expected. We will try to keep downtimes to a minimum, but recommend that users avoid critical tasks, large data uploads, or DOI requests during this time.

We apologize for any inconvenience.

affine.m 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554
  1. % Using 2D or 3D affine matrix to rotate, translate, scale, reflect and
  2. % shear a 2D image or 3D volume. 2D image is represented by a 2D matrix,
  3. % 3D volume is represented by a 3D matrix, and data type can be real
  4. % integer or floating-point.
  5. %
  6. % You may notice that MATLAB has a function called 'imtransform.m' for
  7. % 2D spatial transformation. However, keep in mind that 'imtransform.m'
  8. % assumes y for the 1st dimension, and x for the 2nd dimension. They are
  9. % equivalent otherwise.
  10. %
  11. % In addition, if you adjust the 'new_elem_size' parameter, this 'affine.m'
  12. % is equivalent to 'interp2.m' for 2D image, and equivalent to 'interp3.m'
  13. % for 3D volume.
  14. %
  15. % Usage: [new_img new_M] = ...
  16. % affine(old_img, old_M, [new_elem_size], [verbose], [bg], [method]);
  17. %
  18. % old_img - original 2D image or 3D volume. We assume x for the 1st
  19. % dimension, y for the 2nd dimension, and z for the 3rd
  20. % dimension.
  21. %
  22. % old_M - a 3x3 2D affine matrix for 2D image, or a 4x4 3D affine
  23. % matrix for 3D volume. We assume x for the 1st dimension,
  24. % y for the 2nd dimension, and z for the 3rd dimension.
  25. %
  26. % new_elem_size (optional) - size of voxel along x y z direction for
  27. % a transformed 3D volume, or size of pixel along x y for
  28. % a transformed 2D image. We assume x for the 1st dimension
  29. % y for the 2nd dimension, and z for the 3rd dimension.
  30. % 'new_elem_size' is 1 if it is default or empty.
  31. %
  32. % You can increase its value to decrease the resampling rate,
  33. % and make the 2D image or 3D volume more coarse. It works
  34. % just like 'interp3'.
  35. %
  36. % verbose (optional) - 1, 0
  37. % 1: show transforming progress in percentage
  38. % 2: progress will not be displayed
  39. % 'verbose' is 1 if it is default or empty.
  40. %
  41. % bg (optional) - background voxel intensity in any extra corner that
  42. % is caused by the interpolation. 0 in most cases. If it is
  43. % default or empty, 'bg' will be the average of two corner
  44. % voxel intensities in original data.
  45. %
  46. % method (optional) - 1, 2, or 3
  47. % 1: for Trilinear interpolation
  48. % 2: for Nearest Neighbor interpolation
  49. % 3: for Fischer's Bresenham interpolation
  50. % 'method' is 1 if it is default or empty.
  51. %
  52. % new_img - transformed 2D image or 3D volume
  53. %
  54. % new_M - transformed affine matrix
  55. %
  56. % Example 1 (3D rotation):
  57. % load mri.mat; old_img = double(squeeze(D));
  58. % old_M = [0.88 0.5 3 -90; -0.5 0.88 3 -126; 0 0 2 -72; 0 0 0 1];
  59. % new_img = affine(old_img, old_M, 2);
  60. % [x y z] = meshgrid(1:128,1:128,1:27);
  61. % sz = size(new_img);
  62. % [x1 y1 z1] = meshgrid(1:sz(2),1:sz(1),1:sz(3));
  63. % figure; slice(x, y, z, old_img, 64, 64, 13.5);
  64. % shading flat; colormap(map); view(-66, 66);
  65. % figure; slice(x1, y1, z1, new_img, sz(1)/2, sz(2)/2, sz(3)/2);
  66. % shading flat; colormap(map); view(-66, 66);
  67. %
  68. % Example 2 (2D interpolation):
  69. % load mri.mat; old_img=D(:,:,1,13)';
  70. % old_M = [1 0 0; 0 1 0; 0 0 1];
  71. % new_img = affine(old_img, old_M, [.2 .4]);
  72. % figure; image(old_img); colormap(map);
  73. % figure; image(new_img); colormap(map);
  74. %
  75. % This program is inspired by:
  76. % SPM5 Software from Wellcome Trust Centre for Neuroimaging
  77. % http://www.fil.ion.ucl.ac.uk/spm/software
  78. % Fischer, J., A. del Rio (2004). A Fast Method for Applying Rigid
  79. % Transformations to Volume Data, WSCG2004 Conference.
  80. % http://wscg.zcu.cz/wscg2004/Papers_2004_Short/M19.pdf
  81. %
  82. % - Jimmy Shen (jimmy@rotman-baycrest.on.ca)
  83. %
  84. function [new_img, new_M] = affine(old_img, old_M, new_elem_size, verbose, bg, method)
  85. if ~exist('old_img','var') | ~exist('old_M','var')
  86. error('Usage: [new_img new_M] = affine(old_img, old_M, [new_elem_size], [verbose], [bg], [method]);');
  87. end
  88. if ndims(old_img) == 3
  89. if ~isequal(size(old_M),[4 4])
  90. error('old_M should be a 4x4 affine matrix for 3D volume.');
  91. end
  92. elseif ndims(old_img) == 2
  93. if ~isequal(size(old_M),[3 3])
  94. error('old_M should be a 3x3 affine matrix for 2D image.');
  95. end
  96. else
  97. error('old_img should be either 2D image or 3D volume.');
  98. end
  99. if ~exist('new_elem_size','var') | isempty(new_elem_size)
  100. new_elem_size = [1 1 1];
  101. elseif length(new_elem_size) < 2
  102. new_elem_size = new_elem_size(1)*ones(1,3);
  103. elseif length(new_elem_size) < 3
  104. new_elem_size = [new_elem_size(:); 1]';
  105. end
  106. if ~exist('method','var') | isempty(method)
  107. method = 1;
  108. elseif ~exist('bresenham_line3d.m','file') & method == 3
  109. error([char(10) char(10) 'Please download 3D Bresenham''s line generation program from:' char(10) char(10) 'http://www.mathworks.com/matlabcentral/fileexchange/loadFile.do?objectId=21057' char(10) char(10) 'to test Fischer''s Bresenham interpolation method.' char(10) char(10)]);
  110. end
  111. % Make compatible to MATLAB earlier than version 7 (R14), which
  112. % can only perform arithmetic on double data type
  113. %
  114. old_img = double(old_img);
  115. old_dim = size(old_img);
  116. if ~exist('bg','var') | isempty(bg)
  117. bg = mean([old_img(1) old_img(end)]);
  118. end
  119. if ~exist('verbose','var') | isempty(verbose)
  120. verbose = 1;
  121. end
  122. if ndims(old_img) == 2
  123. old_dim(3) = 1;
  124. old_M = old_M(:, [1 2 3 3]);
  125. old_M = old_M([1 2 3 3], :);
  126. old_M(3,:) = [0 0 1 0];
  127. old_M(:,3) = [0 0 1 0]';
  128. end
  129. % Vertices of img in voxel
  130. %
  131. XYZvox = [ 1 1 1
  132. 1 1 old_dim(3)
  133. 1 old_dim(2) 1
  134. 1 old_dim(2) old_dim(3)
  135. old_dim(1) 1 1
  136. old_dim(1) 1 old_dim(3)
  137. old_dim(1) old_dim(2) 1
  138. old_dim(1) old_dim(2) old_dim(3) ]';
  139. old_R = old_M(1:3,1:3);
  140. old_T = old_M(1:3,4);
  141. % Vertices of img in millimeter
  142. %
  143. XYZmm = old_R*(XYZvox-1) + repmat(old_T, [1, 8]);
  144. % Make scale of new_M according to new_elem_size
  145. %
  146. new_M = diag([new_elem_size 1]);
  147. % Make translation so minimum vertex is moved to [1,1,1]
  148. %
  149. new_M(1:3,4) = round( min(XYZmm,[],2) );
  150. % New dimensions will be the maximum vertices in XYZ direction (dim_vox)
  151. % i.e. compute dim_vox via dim_mm = R*(dim_vox-1)+T
  152. % where, dim_mm = round(max(XYZmm,[],2));
  153. %
  154. new_dim = ceil(new_M(1:3,1:3) \ ( round(max(XYZmm,[],2))-new_M(1:3,4) )+1)';
  155. % Initialize new_img with new_dim
  156. %
  157. new_img = zeros(new_dim(1:3));
  158. % Mask out any changes from Z axis of transformed volume, since we
  159. % will traverse it voxel by voxel below. We will only apply unit
  160. % increment of mask_Z(3,4) to simulate the cursor movement
  161. %
  162. % i.e. we will use mask_Z * new_XYZvox to replace new_XYZvox
  163. %
  164. mask_Z = diag(ones(1,4));
  165. mask_Z(3,3) = 0;
  166. % It will be easier to do the interpolation if we invert the process
  167. % by not traversing the original volume. Instead, we traverse the
  168. % transformed volume, and backproject each voxel in the transformed
  169. % volume back into the original volume. If the backprojected voxel
  170. % in original volume is within its boundary, the intensity of that
  171. % voxel can be used by the cursor location in the transformed volume.
  172. %
  173. % First, we traverse along Z axis of transformed volume voxel by voxel
  174. %
  175. for z = 1:new_dim(3)
  176. if verbose & ~mod(z,10)
  177. fprintf('%.2f percent is done.\n', 100*z/new_dim(3));
  178. end
  179. % We need to find out the mapping from voxel in the transformed
  180. % volume (new_XYZvox) to voxel in the original volume (old_XYZvox)
  181. %
  182. % The following equation works, because they all equal to XYZmm:
  183. % new_R*(new_XYZvox-1) + new_T == old_R*(old_XYZvox-1) + old_T
  184. %
  185. % We can use modified new_M1 & old_M1 to substitute new_M & old_M
  186. % new_M1 * new_XYZvox == old_M1 * old_XYZvox
  187. %
  188. % where: M1 = M; M1(:,4) = M(:,4) - sum(M(:,1:3),2);
  189. % and: M(:,4) == [T; 1] == sum(M1,2)
  190. %
  191. % Therefore: old_XYZvox = old_M1 \ new_M1 * new_XYZvox;
  192. %
  193. % Since we are traverse Z axis, and new_XYZvox is replaced
  194. % by mask_Z * new_XYZvox, the above formula can be rewritten
  195. % as: old_XYZvox = old_M1 \ new_M1 * mask_Z * new_XYZvox;
  196. %
  197. % i.e. we find the mapping from new_XYZvox to old_XYZvox:
  198. % M = old_M1 \ new_M1 * mask_Z;
  199. %
  200. % First, compute modified old_M1 & new_M1
  201. %
  202. old_M1 = old_M; old_M1(:,4) = old_M(:,4) - sum(old_M(:,1:3),2);
  203. new_M1 = new_M; new_M1(:,4) = new_M(:,4) - sum(new_M(:,1:3),2);
  204. % Then, apply unit increment of mask_Z(3,4) to simulate the
  205. % cursor movement
  206. %
  207. mask_Z(3,4) = z;
  208. % Here is the mapping from new_XYZvox to old_XYZvox
  209. %
  210. M = old_M1 \ new_M1 * mask_Z;
  211. switch method
  212. case 1
  213. new_img(:,:,z) = trilinear(old_img, new_dim, old_dim, M, bg);
  214. case 2
  215. new_img(:,:,z) = nearest_neighbor(old_img, new_dim, old_dim, M, bg);
  216. case 3
  217. new_img(:,:,z) = bresenham(old_img, new_dim, old_dim, M, bg);
  218. end
  219. end; % for z
  220. if ndims(old_img) == 2
  221. new_M(3,:) = [];
  222. new_M(:,3) = [];
  223. end
  224. return; % affine
  225. %--------------------------------------------------------------------
  226. function img_slice = trilinear(img, dim1, dim2, M, bg)
  227. img_slice = zeros(dim1(1:2));
  228. TINY = 5e-2; % tolerance
  229. % Dimension of transformed 3D volume
  230. %
  231. xdim1 = dim1(1);
  232. ydim1 = dim1(2);
  233. % Dimension of original 3D volume
  234. %
  235. xdim2 = dim2(1);
  236. ydim2 = dim2(2);
  237. zdim2 = dim2(3);
  238. % initialize new_Y accumulation
  239. %
  240. Y2X = 0;
  241. Y2Y = 0;
  242. Y2Z = 0;
  243. for y = 1:ydim1
  244. % increment of new_Y accumulation
  245. %
  246. Y2X = Y2X + M(1,2); % new_Y to old_X
  247. Y2Y = Y2Y + M(2,2); % new_Y to old_Y
  248. Y2Z = Y2Z + M(3,2); % new_Y to old_Z
  249. % backproject new_Y accumulation and translation to old_XYZ
  250. %
  251. old_X = Y2X + M(1,4);
  252. old_Y = Y2Y + M(2,4);
  253. old_Z = Y2Z + M(3,4);
  254. for x = 1:xdim1
  255. % accumulate the increment of new_X, and apply it
  256. % to the backprojected old_XYZ
  257. %
  258. old_X = M(1,1) + old_X ;
  259. old_Y = M(2,1) + old_Y ;
  260. old_Z = M(3,1) + old_Z ;
  261. % within boundary of original image
  262. %
  263. if ( old_X > 1-TINY & old_X < xdim2+TINY & ...
  264. old_Y > 1-TINY & old_Y < ydim2+TINY & ...
  265. old_Z > 1-TINY & old_Z < zdim2+TINY )
  266. % Calculate distance of old_XYZ to its neighbors for
  267. % weighted intensity average
  268. %
  269. dx = old_X - floor(old_X);
  270. dy = old_Y - floor(old_Y);
  271. dz = old_Z - floor(old_Z);
  272. x000 = floor(old_X);
  273. x100 = x000 + 1;
  274. if floor(old_X) < 1
  275. x000 = 1;
  276. x100 = x000;
  277. elseif floor(old_X) > xdim2-1
  278. x000 = xdim2;
  279. x100 = x000;
  280. end
  281. x010 = x000;
  282. x001 = x000;
  283. x011 = x000;
  284. x110 = x100;
  285. x101 = x100;
  286. x111 = x100;
  287. y000 = floor(old_Y);
  288. y010 = y000 + 1;
  289. if floor(old_Y) < 1
  290. y000 = 1;
  291. y100 = y000;
  292. elseif floor(old_Y) > ydim2-1
  293. y000 = ydim2;
  294. y010 = y000;
  295. end
  296. y100 = y000;
  297. y001 = y000;
  298. y101 = y000;
  299. y110 = y010;
  300. y011 = y010;
  301. y111 = y010;
  302. z000 = floor(old_Z);
  303. z001 = z000 + 1;
  304. if floor(old_Z) < 1
  305. z000 = 1;
  306. z001 = z000;
  307. elseif floor(old_Z) > zdim2-1
  308. z000 = zdim2;
  309. z001 = z000;
  310. end
  311. z100 = z000;
  312. z010 = z000;
  313. z110 = z000;
  314. z101 = z001;
  315. z011 = z001;
  316. z111 = z001;
  317. x010 = x000;
  318. x001 = x000;
  319. x011 = x000;
  320. x110 = x100;
  321. x101 = x100;
  322. x111 = x100;
  323. v000 = double(img(x000, y000, z000));
  324. v010 = double(img(x010, y010, z010));
  325. v001 = double(img(x001, y001, z001));
  326. v011 = double(img(x011, y011, z011));
  327. v100 = double(img(x100, y100, z100));
  328. v110 = double(img(x110, y110, z110));
  329. v101 = double(img(x101, y101, z101));
  330. v111 = double(img(x111, y111, z111));
  331. img_slice(x,y) = v000*(1-dx)*(1-dy)*(1-dz) + ...
  332. v010*(1-dx)*dy*(1-dz) + ...
  333. v001*(1-dx)*(1-dy)*dz + ...
  334. v011*(1-dx)*dy*dz + ...
  335. v100*dx*(1-dy)*(1-dz) + ...
  336. v110*dx*dy*(1-dz) + ...
  337. v101*dx*(1-dy)*dz + ...
  338. v111*dx*dy*dz;
  339. else
  340. img_slice(x,y) = bg;
  341. end % if boundary
  342. end % for x
  343. end % for y
  344. return; % trilinear
  345. %--------------------------------------------------------------------
  346. function img_slice = nearest_neighbor(img, dim1, dim2, M, bg)
  347. img_slice = zeros(dim1(1:2));
  348. % Dimension of transformed 3D volume
  349. %
  350. xdim1 = dim1(1);
  351. ydim1 = dim1(2);
  352. % Dimension of original 3D volume
  353. %
  354. xdim2 = dim2(1);
  355. ydim2 = dim2(2);
  356. zdim2 = dim2(3);
  357. % initialize new_Y accumulation
  358. %
  359. Y2X = 0;
  360. Y2Y = 0;
  361. Y2Z = 0;
  362. for y = 1:ydim1
  363. % increment of new_Y accumulation
  364. %
  365. Y2X = Y2X + M(1,2); % new_Y to old_X
  366. Y2Y = Y2Y + M(2,2); % new_Y to old_Y
  367. Y2Z = Y2Z + M(3,2); % new_Y to old_Z
  368. % backproject new_Y accumulation and translation to old_XYZ
  369. %
  370. old_X = Y2X + M(1,4);
  371. old_Y = Y2Y + M(2,4);
  372. old_Z = Y2Z + M(3,4);
  373. for x = 1:xdim1
  374. % accumulate the increment of new_X and apply it
  375. % to the backprojected old_XYZ
  376. %
  377. old_X = M(1,1) + old_X ;
  378. old_Y = M(2,1) + old_Y ;
  379. old_Z = M(3,1) + old_Z ;
  380. xi = round(old_X);
  381. yi = round(old_Y);
  382. zi = round(old_Z);
  383. % within boundary of original image
  384. %
  385. if ( xi >= 1 & xi <= xdim2 & ...
  386. yi >= 1 & yi <= ydim2 & ...
  387. zi >= 1 & zi <= zdim2 )
  388. img_slice(x,y) = img(xi,yi,zi);
  389. else
  390. img_slice(x,y) = bg;
  391. end % if boundary
  392. end % for x
  393. end % for y
  394. return; % nearest_neighbor
  395. %--------------------------------------------------------------------
  396. function img_slice = bresenham(img, dim1, dim2, M, bg)
  397. img_slice = zeros(dim1(1:2));
  398. % Dimension of transformed 3D volume
  399. %
  400. xdim1 = dim1(1);
  401. ydim1 = dim1(2);
  402. % Dimension of original 3D volume
  403. %
  404. xdim2 = dim2(1);
  405. ydim2 = dim2(2);
  406. zdim2 = dim2(3);
  407. for y = 1:ydim1
  408. start_old_XYZ = round(M*[0 y 0 1]');
  409. end_old_XYZ = round(M*[xdim1 y 0 1]');
  410. [X Y Z] = bresenham_line3d(start_old_XYZ, end_old_XYZ);
  411. % line error correction
  412. %
  413. % del = end_old_XYZ - start_old_XYZ;
  414. % del_dom = max(del);
  415. % idx_dom = find(del==del_dom);
  416. % idx_dom = idx_dom(1);
  417. % idx_other = [1 2 3];
  418. % idx_other(idx_dom) = [];
  419. %del_x1 = del(idx_other(1));
  420. % del_x2 = del(idx_other(2));
  421. % line_slope = sqrt((del_x1/del_dom)^2 + (del_x2/del_dom)^2 + 1);
  422. % line_error = line_slope - 1;
  423. % line error correction removed because it is too slow
  424. for x = 1:xdim1
  425. % rescale ratio
  426. %
  427. i = round(x * length(X) / xdim1);
  428. if i < 1
  429. i = 1;
  430. elseif i > length(X)
  431. i = length(X);
  432. end
  433. xi = X(i);
  434. yi = Y(i);
  435. zi = Z(i);
  436. % within boundary of the old XYZ space
  437. %
  438. if ( xi >= 1 & xi <= xdim2 & ...
  439. yi >= 1 & yi <= ydim2 & ...
  440. zi >= 1 & zi <= zdim2 )
  441. img_slice(x,y) = img(xi,yi,zi);
  442. % if line_error > 1
  443. % x = x + 1;
  444. % if x <= xdim1
  445. % img_slice(x,y) = img(xi,yi,zi);
  446. % line_error = line_slope - 1;
  447. % end
  448. % end % if line_error
  449. % line error correction removed because it is too slow
  450. else
  451. img_slice(x,y) = bg;
  452. end % if boundary
  453. end % for x
  454. end % for y
  455. return; % bresenham