spm_diffeo.m 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455
  1. function varargout = spm_diffeo(varargin)
  2. % Mex function called for image registration stuff
  3. %
  4. %_______________________________________________________________________
  5. %
  6. % FORMAT u = spm_diffeo('vel2mom', v, param)
  7. % v - velocity (flow) field n1*n2*n3*3.
  8. % param - 8 parameters (settings)
  9. % - [1][2][3] Voxel sizes
  10. % - [4][5][6][7][8] Regularisation parameters
  11. % - [4] Absolute displacements need to be penalised by a tiny
  12. % amount. The first element encodes the amount of
  13. % penalty on these. Ideally, absolute displacements
  14. % should not be penalised, but it is usually necessary
  15. % for technical reasons.
  16. % - [5] The `membrane energy' of the deformation is penalised,
  17. % usually by a relatively small amount. This penalises
  18. % the sum of squares of the derivatives of the velocity
  19. % field (ie the sum of squares of the elements of the
  20. % Jacobian tensors).
  21. % - [6] The `bending energy' is penalised (3rd element). This
  22. % penalises the sum of squares of the 2nd derivatives of
  23. % the velocity.
  24. % - [7][8] Linear elasticity regularisation is also included.
  25. % The first parameter (mu) is similar to that for
  26. % linear elasticity, except it penalises the sum of
  27. % squares of the Jacobian tensors after they have been
  28. % made symmetric (by averaging with the transpose).
  29. % This term essentially penalises length changes,
  30. % without penalising rotations.
  31. % The final term also relates to linear elasticity,
  32. % and is the weight that denotes how much to penalise
  33. % changes to the divergence of the velocities (lambda).
  34. % This divergence is a measure of the rate of volumetric
  35. % expansion or contraction.
  36. % u - `momentum' field n1*n2*n3*3.
  37. %
  38. % Convert a velocity field to a momentum field by u = A*v, where
  39. % A is the large sparse matrix encoding some form of regularisation.
  40. % v and m are single precision floating point.
  41. %
  42. %_______________________________________________________________________
  43. %
  44. % FORMAT v = spm_diffeo('mom2vel',g, param)
  45. % v - the solution n1*n2*n3*3
  46. % g - parameterisation of first derivatives
  47. % param - 10 parameters (settings)
  48. % - [1][2][3] Voxel sizes
  49. % - [4][5][6][7][8] Regularisation settings (see vel2mom).
  50. % - [9] Number of Full Multigrid cycles.
  51. % - [10] Number of relaxation iterations per cycle.
  52. %
  53. % Solve equations using a Full Multigrid method. See Press et al
  54. % for more information.
  55. % v = inv(A)*g
  56. % g and v are both single precision floating point.
  57. %
  58. %_______________________________________________________________________
  59. %
  60. % FORMAT v = spm_diffeo('fmg',H, g, param)
  61. % v - the solution n1*n2*n3*3
  62. % H - parameterisation of 2nd derivatives
  63. % g - parameterisation of first derivatives
  64. % param - 10 parameters (settings)
  65. % - [1][2][3] Voxel sizes
  66. % - [4][5][6][7][8] Regularisation settings (see vel2mom).
  67. % - [9] Number of Full Multigrid cycles.
  68. % - [10] Number of relaxation iterations per cycle.
  69. %
  70. % Solve equations using a Full Multigrid method, but using Hessian of
  71. % the matching term. See Press et al for more information.
  72. % v = inv(A+H)*g
  73. % H, g and v are all single precision floating point.
  74. %
  75. %_______________________________________________________________________
  76. %
  77. % FORMAT v = spm_diffeo('cgs',H, g, param)
  78. % v - the solution
  79. % H - parameterisation of 2nd derivatives
  80. % g - parameterisation of first derivatives
  81. % param - 10 parameters (settings)
  82. % - [1][2][3] Voxel sizes
  83. % - [4][5][6][7][8] Regularisation settings (see vel2mom).
  84. % - [9] Tolerance. Indicates required degree of accuracy.
  85. % - [10] Maximum number of iterations.
  86. %
  87. % This is for solving a set of equations using a conjugate gradient
  88. % solver. This method is less efficient than the Full Multigrid, and
  89. % is included for illustrative purposes.
  90. % v = inv(A+H)*g
  91. % H, g and v are all single precision floating point.
  92. %
  93. %_______________________________________________________________________
  94. %
  95. % FORMAT F = spm_diffeo('kernel',d,prm)
  96. % d - image dimensions
  97. % prm - 8 parameters (settings).
  98. % These are described above (for 'vel2mom').
  99. % F - The differential operator encoded as an image (or images).
  100. % Convolving a velocity field by this will give the momentum.
  101. %
  102. % Note that a smaller (3D) kernel is obtained when the linear
  103. % elasticity settings are all zero. If any of the linear
  104. % elasticity settings are non-zero, the resulting kernel is
  105. % represented by a 5D array. For the 3D form, the voxel sizes
  106. % need to be incorporated as an additional scaling of the kernel.
  107. % See the code in spm_shoot_greens.m for an illustration.
  108. %
  109. %_______________________________________________________________________
  110. %
  111. % FORMAT y3 = spm_diffeo('comp',y1,y2)
  112. % y1, y2 - deformation fields n1*n2*n3*3.
  113. % y3 - deformation field field n1*n2*n3*3.
  114. %
  115. % Composition of two deformations y3 = y1(y2)
  116. % y1, y2 and y3 are single precision floating point.
  117. %
  118. %
  119. % FORMAT [y3,J3] = spm_diffeo('comp', y1, y2, J1, J2)
  120. % y1, y2 - deformation fields n1*n2*n3*3.
  121. % y3 - deformation field n1*n2*n3*3.
  122. % J1, J2 - Jacobian tensor fields n1*n2*n3*3*3.
  123. % J3 - Jacobian tensor field n1*n2*n3*3*3.
  124. %
  125. % Composition of two deformations, with their Jacobian fields.
  126. % All fields are single precision floating point.
  127. %
  128. %_______________________________________________________________________
  129. %
  130. % FORMAT iy = spm_diffeo('invdef',y,d,M1,M2);
  131. %
  132. % iy - inverted deformation field of size d(1)*d(2)*d(3)*3.
  133. % y - original deformation field.
  134. % M1 - An affine mapping from mm to voxels in the co-ordinate
  135. % system of the inverse deformation field.
  136. % M2 - An affine mapping from voxels to mm in the co-ordinate
  137. % system of the forward deformation field.
  138. %
  139. % Inversion of a deformation field.
  140. %
  141. % The field is assumed to consist of a piecewise affine transformations,
  142. % whereby each cube jointing 8 neighbouring voxels contains eight
  143. % tetrahedra. The mapping within each tetrahedron is assumed to be
  144. % affine.
  145. %
  146. % Reference:
  147. % J. Ashburner, J. Andersson and K. J. Friston (2000).
  148. % "Image Registration using a Symmetric Prior - in Three-Dimensions".
  149. % Human Brain Mapping 9(4):212-225 (appendix).
  150. %_______________________________________________________________________
  151. %
  152. % FORMAT [f,dfx,dfy,dfz] = spm_diffeo('bsplins', c, y,d)
  153. % c - input image(s) of B-spline coefficients n1*n2*n3*n4
  154. % - see 'bsplinc'
  155. % y - points to sample n1*n2*n3*3
  156. % d(1:3) - degree of B-spline (from 0 to 7) along different dimensions
  157. % - these must be same as used by 'bsplinc'
  158. % d(4:6) - 1/0 to indicate wrapping along the dimensions
  159. %
  160. % f - output image n1*n2*n3*n4
  161. % dfx,dfy,dfz - sampled first derivatives
  162. %
  163. % c, f and y are single precision floating point.
  164. %
  165. % This function takes B-spline basis coefficients from spm_bsplinc,
  166. % and re-convolves them with B-splines centred at the new sample points.
  167. %
  168. % Note that nearest neighbour interpolation is used instead of 0th
  169. % degree B-splines, and the derivatives of trilinear interpolation are
  170. % returned insted of those of 1st degree B-splines. The difference is
  171. % extremely subtle.
  172. %
  173. % c, f and y are single precision floating point.
  174. %
  175. % References:
  176. % M. Unser, A. Aldroubi and M. Eden.
  177. % "B-Spline Signal Processing: Part I-Theory,"
  178. % IEEE Transactions on Signal Processing 41(2):821-832 (1993).
  179. %
  180. % M. Unser, A. Aldroubi and M. Eden.
  181. % "B-Spline Signal Processing: Part II-Efficient Design and
  182. % Applications,"
  183. % IEEE Transactions on Signal Processing 41(2):834-848 (1993).
  184. %
  185. % M. Unser.
  186. % "Splines: A Perfect Fit for Signal and Image Processing,"
  187. % IEEE Signal Processing Magazine, 16(6):22-38 (1999)
  188. %
  189. % P. Thevenaz and T. Blu and M. Unser.
  190. % "Interpolation Revisited"
  191. % IEEE Transactions on Medical Imaging 19(7):739-758 (2000).
  192. %
  193. %_______________________________________________________________________
  194. %
  195. % FORMAT c = spm_diffeo('bsplinc',f,d)
  196. % f - an image
  197. % d(1:3) - degree of B-spline (from 0 to 7) along different dimensions
  198. % d(4:6) - 1/0 to indicate wrapping along the dimensions
  199. % c - returned volume of B-spline coefficients
  200. %
  201. % This function deconvolves B-splines from f, returning
  202. % coefficients, c. These coefficients are then passed to 'bsplins'
  203. % in order to sample the data using B-spline interpolation.
  204. %
  205. %_______________________________________________________________________
  206. %
  207. % FORMAT f2 = spm_diffeo('samp', f1, y)
  208. % f1 - input image(s) n1*n2*n3*n4
  209. % y - points to sample n1*n2*n3*3
  210. % f2 - output image n1*n2*n3*n4
  211. %
  212. % Sample a function according to a deformation using trilinear interp.
  213. % f2 = f1(y)
  214. % f1, f2 and y are single precision floating point.
  215. % Uses boundary condiditions that wrap around (circulant - identical to
  216. % the 'pullc' option - but retained for backward compatibility).
  217. %
  218. %_______________________________________________________________________
  219. %
  220. % FORMAT f2 = spm_diffeo('pull', f1, y)
  221. % f1 - input image(s) n1*n2*n3*n4
  222. % y - points to sample n1*n2*n3*3
  223. % f2 - output image n1*n2*n3*n4
  224. %
  225. % Sample a function according to a deformation using trilinear interp.
  226. % f2 = f1(y)
  227. % f1, f2 and y are single precision floating point.
  228. % Values sampled outside the field of view of f1 are assigned a value
  229. % of NaN.
  230. %
  231. %_______________________________________________________________________
  232. %
  233. % FORMAT f2 = spm_diffeo('pullc', f1, y)
  234. % f1 - input image(s) n1*n2*n3*n4
  235. % y - points to sample n1*n2*n3*3
  236. % f2 - output image n1*n2*n3*n4
  237. %
  238. % Sample a function according to a deformation using trilinear interp.
  239. % f2 = f1(y)
  240. % f1, f2 and y are single precision floating point.
  241. % Uses boundary condiditions that wrap around (circulant - identical to
  242. % the 'samp' option).
  243. %
  244. %_______________________________________________________________________
  245. %
  246. % FORMAT f2 = spm_diffeo('push', f1, y)
  247. % f1 - input image(s) n1*n2*n3*n4
  248. % y - points to sample n1*n2*n3*3
  249. % f2 - output image n1*n2*n3*n4
  250. %
  251. % Push values of a function according to a deformation. Note that the
  252. % deformation should be the inverse of the one used with 'samp' or
  253. % 'bsplins'. f1, f2 and y are single precision floating point.
  254. % Voxels in f1 that would be pushed outside the field of view of f2
  255. % are ignored.
  256. %
  257. %_______________________________________________________________________
  258. %
  259. % FORMAT f2 = spm_diffeo('pushc', f1, y)
  260. % f1 - input image(s) n1*n2*n3*n4
  261. % y - points to sample n1*n2*n3*3
  262. % f2 - output image n1*n2*n3*n4
  263. %
  264. % Push values of a function according to a deformation, but using
  265. % circulant boundary conditions. Data wraps around (circulant).
  266. % f1, f2 and y are single precision floating point.
  267. %
  268. %_______________________________________________________________________
  269. %
  270. % FORMAT ut = spm_diffeo('pushg', u0, y)
  271. % u0 - input momentum n1*n2*n3*3
  272. % y - points to sample n1*n2*n3*3
  273. % ut - output momentum n1*n2*n3*3
  274. %
  275. % FORMAT ut = spm_diffeo('pushg', u0, y)
  276. % u0 - input momentum n1*n2*n3*3
  277. % y - points to sample n1*n2*n3*3
  278. % J - Jacobian tensor field of y n1*n2*n3*3*3
  279. % ut - output momentum n1*n2*n3*3
  280. %
  281. % Push values of a momentum field according to a deformation using
  282. % circulant boundary conditions. This essentially computes
  283. % (Ad_y)^* u = |det dy| (dy)^T u(y), which is a key to the
  284. % EPdiff equations used for geodesic shooting.
  285. % u0, ut and y are single precision floating point.
  286. %
  287. %_______________________________________________________________________
  288. %
  289. % FORMAT f2 = spm_diffeo('resize', f1, dim)
  290. % f1 - input fields n1*n2*n3*n4
  291. % f2 - output field dim1*dim2*dim3*n4
  292. % dim - output dimensions
  293. %
  294. % Resize a field according to dimensions dim. This is a component of
  295. % the multigrid approach, and is used for prolongation.
  296. %
  297. %_______________________________________________________________________
  298. %
  299. % FORMAT v2 = spm_diffeo('restrict', v1)
  300. % v1 - input fields n1*n2*n3*n4
  301. % v2 - output field dim1*dim2*dim3*n4
  302. %
  303. % Restricts a field such that its dimensions are approximately half
  304. % their original. This is a component of the multigrid approach.
  305. %
  306. %_______________________________________________________________________
  307. %
  308. % FORMAT J = spm_diffeo('def2jac',y)
  309. % y - Deformation field
  310. % J - Jacobian tensor field of y
  311. %
  312. % Compute Jacobian tensors from a deformation.
  313. %
  314. %_______________________________________________________________________
  315. %
  316. % FORMAT J = spm_diffeo('def2det',y)
  317. % y - Deformation field
  318. % j - Jacobian determinant field of y
  319. %
  320. % Compute Jacobian determinants from a deformation.
  321. %
  322. %_______________________________________________________________________
  323. %
  324. % FORMAT j = spm_diffeo('det',J)
  325. % J - Jacobian tensor field
  326. % j - Jacobian determinant field
  327. %
  328. % Compute determinants of Jacobian tensors.
  329. %
  330. %_______________________________________________________________________
  331. %
  332. % FORMAT dv = spm_diffeo('div',v)
  333. % v - velocity field
  334. % dv - divergences of velocity field
  335. %
  336. % Computes divergence from velocity field. This is indicative of rates
  337. % of volumetric expansion/contraction.
  338. %
  339. %_______________________________________________________________________
  340. %
  341. % FORMAT [y,J] = spm_diffeo('smalldef',v,s)
  342. % v - velocity field
  343. % s - scaling factor
  344. % y - small deformation
  345. % J - approximate Jacobian tensors of small deformation (computed via
  346. % a matrix exponsntial of the Jacobians of the velocity field).
  347. %
  348. % This function is used for each time step of geodesic shooting. It may
  349. % change in future to use some form of Pade approximation of the
  350. % small deformation.
  351. %
  352. %_______________________________________________________________________
  353. %
  354. % FORMAT v3 = spm_diffeo('brc', v1, v2)
  355. % v1, v2, v3 - flow fields n1*n2*n3*3
  356. %
  357. % Lie Bracket. Useful for many things
  358. % e.g. Baker-Campbell-Haussdorf series expansion.
  359. % The Lie bracket is denoted by
  360. % v3 = [v1,v2]
  361. % and on scalar fields, is computed by
  362. % v3 = J1*v2 - J2*v1, where J1 and J2 are the Jacobian
  363. % tensor fields. For matrices, the Lie bracket is simply
  364. % [A,B] = A*B-B*A
  365. %
  366. %_______________________________________________________________________
  367. %
  368. % FORMAT t = spm_diffeo('trapprox',H, param)
  369. % v - the solution n1*n2*n3*3
  370. % H - parameterisation of 2nd derivatives
  371. % param - 10 parameters (settings)
  372. % - [1][2][3] Voxel sizes
  373. % - [4][5][6][7][8] Regularisation settings (see vel2mom).
  374. % t - approximation of [trace((L+H)\L) trace((L+H)\H)];
  375. %
  376. % Generate an approximation of Trace((L+H)\L) and Trace((L+H)\H) for
  377. % to give a ball-park figure for the "degrees of freedom" in Laplace
  378. % approximations. L is the regulariser in sparse matrix form. The
  379. % approximation is a poor one, which assumes all the off-diagonals of L
  380. % are 0.
  381. % H is single precision floating point.
  382. %
  383. %_______________________________________________________________________
  384. %
  385. % FORMAT v = spm_diffeo('dartel',v,g,f,param)
  386. % v - flow field n1*n2*n3*3 (single precision float)
  387. % g - first image n1*n2*n3*n4 (single precision float)
  388. % f - second image n1*n2*n3*n4 (single precision float)
  389. % param - 9 parameters (settings)
  390. % - [1] Regularisation type, can take values of
  391. % - 0 Linear elasticity
  392. % - 1 Membrane energy
  393. % - 2 Bending energy
  394. % - [2][3][4] Regularisation parameters
  395. % - For "membrane energy", the parameters are
  396. % lambda, unused and id.
  397. % - For "linear elasticity", the parameters are
  398. % mu, lambda, and id
  399. % - For "bending energy", the parameters are
  400. % lambda, id1 and id2, such that regularisation is by
  401. % (-lambda*\grad^2 + id1)^2 + id2
  402. % - [5] Levenberg-Marquardt regularisation
  403. % - [6] Number of Full Multigrid cycles
  404. % - [7] Number of relaxation iterations per cycle
  405. % - [8] K, such that 2^K time points are used to
  406. % generate the deformations. A value of zero
  407. % indicates a small deformation model.
  408. % - [9] code of 0, 1 or 2.
  409. % 0 - asymmetric sums of squares objective function.
  410. % 1 - symmetric sums of squares objective function.
  411. % 2 - assumes multinomial distribution, where template
  412. % encodes the means and interpolation of template
  413. % done using logs and softmax function.
  414. %
  415. % This is for performing a single iteration of the Dartel optimisation.
  416. % All velocity fields and images are represented by single precision floating
  417. % point values. Images can be scalar fields, in which case the objective
  418. % function is the sum of squares difference. Alternatively, images can be
  419. % vector fields, in which case the objective function is the sum of squares
  420. % difference between each scalar field + the sum of squares difference
  421. % between one minus the sum of the scalar fields.
  422. %
  423. %_______________________________________________________________________
  424. %
  425. % FORMAT [y,J] = spm_diffeo('Exp', v, param)
  426. % v - flow field
  427. % J - Jacobian. Usually a tensor field of Jacobian matrices, but can
  428. % be a field of Jacobian determinants.
  429. % param - 2 (or 3) parameters.
  430. % [1] K, the number of recursions (squaring steps), such
  431. % that exponentiation is done using an Euler-like
  432. % integration with 2^K time steps.
  433. % [2] a scaling parameter.
  434. % If there is a third parameter, and it is set to 1, then
  435. % the J will be the Jacobian determinants.
  436. %
  437. % A flow field is "exponentiated" to generate a deformation field
  438. % using a scaling and squaring approach. See the work of Arsigny
  439. % et al, or Cleve Moler's "19 Dubious Ways" papers.
  440. %
  441. %_______________________________________________________________________
  442. %
  443. % Note that the boundary conditions are circulant throughout.
  444. % Interpolation is trilinear, except for the resize function
  445. % which uses a 2nd degree B-spline (without first deconvolving).
  446. %
  447. %_______________________________________________________________________
  448. % Copyright (C) 2012 Wellcome Trust Centre for Neuroimaging
  449. % John Ashburner
  450. % $Id: spm_diffeo.m 7460 2018-10-29 15:55:12Z john $
  451. %-This is merely the help file for the compiled routine
  452. error('spm_diffeo.c not compiled - see Makefile')