radon_transform.py 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536
  1. import numpy as np
  2. from scipy.interpolate import interp1d
  3. from scipy.constants import golden_ratio
  4. from scipy.fft import fft, ifft, fftfreq, fftshift
  5. from ._warps import warp
  6. from ._radon_transform import sart_projection_update
  7. from .._shared.utils import convert_to_float
  8. from warnings import warn
  9. from functools import partial
  10. __all__ = ['radon', 'order_angles_golden_ratio', 'iradon', 'iradon_sart']
  11. def radon(image, theta=None, circle=True, *, preserve_range=False):
  12. """
  13. Calculates the radon transform of an image given specified
  14. projection angles.
  15. Parameters
  16. ----------
  17. image : ndarray
  18. Input image. The rotation axis will be located in the pixel with
  19. indices ``(image.shape[0] // 2, image.shape[1] // 2)``.
  20. theta : array, optional
  21. Projection angles (in degrees). If `None`, the value is set to
  22. np.arange(180).
  23. circle : bool, optional
  24. Assume image is zero outside the inscribed circle, making the
  25. width of each projection (the first dimension of the sinogram)
  26. equal to ``min(image.shape)``.
  27. preserve_range : bool, optional
  28. Whether to keep the original range of values. Otherwise, the input
  29. image is converted according to the conventions of `img_as_float`.
  30. Also see https://scikit-image.org/docs/dev/user_guide/data_types.html
  31. Returns
  32. -------
  33. radon_image : ndarray
  34. Radon transform (sinogram). The tomography rotation axis will lie
  35. at the pixel index ``radon_image.shape[0] // 2`` along the 0th
  36. dimension of ``radon_image``.
  37. References
  38. ----------
  39. .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
  40. Imaging", IEEE Press 1988.
  41. .. [2] B.R. Ramesh, N. Srinivasa, K. Rajgopal, "An Algorithm for Computing
  42. the Discrete Radon Transform With Some Applications", Proceedings of
  43. the Fourth IEEE Region 10 International Conference, TENCON '89, 1989
  44. Notes
  45. -----
  46. Based on code of Justin K. Romberg
  47. (https://www.clear.rice.edu/elec431/projects96/DSP/bpanalysis.html)
  48. """
  49. if image.ndim != 2:
  50. raise ValueError('The input image must be 2-D')
  51. if theta is None:
  52. theta = np.arange(180)
  53. image = convert_to_float(image, preserve_range)
  54. if circle:
  55. shape_min = min(image.shape)
  56. radius = shape_min // 2
  57. img_shape = np.array(image.shape)
  58. coords = np.array(np.ogrid[: image.shape[0], : image.shape[1]], dtype=object)
  59. dist = ((coords - img_shape // 2) ** 2).sum(0)
  60. outside_reconstruction_circle = dist > radius**2
  61. if np.any(image[outside_reconstruction_circle]):
  62. warn(
  63. 'Radon transform: image must be zero outside the '
  64. 'reconstruction circle'
  65. )
  66. # Crop image to make it square
  67. slices = tuple(
  68. (
  69. slice(int(np.ceil(excess / 2)), int(np.ceil(excess / 2) + shape_min))
  70. if excess > 0
  71. else slice(None)
  72. )
  73. for excess in (img_shape - shape_min)
  74. )
  75. padded_image = image[slices]
  76. else:
  77. diagonal = np.sqrt(2) * max(image.shape)
  78. pad = [int(np.ceil(diagonal - s)) for s in image.shape]
  79. new_center = [(s + p) // 2 for s, p in zip(image.shape, pad)]
  80. old_center = [s // 2 for s in image.shape]
  81. pad_before = [nc - oc for oc, nc in zip(old_center, new_center)]
  82. pad_width = [(pb, p - pb) for pb, p in zip(pad_before, pad)]
  83. padded_image = np.pad(image, pad_width, mode='constant', constant_values=0)
  84. # padded_image is always square
  85. if padded_image.shape[0] != padded_image.shape[1]:
  86. raise ValueError('padded_image must be a square')
  87. center = padded_image.shape[0] // 2
  88. radon_image = np.zeros((padded_image.shape[0], len(theta)), dtype=image.dtype)
  89. for i, angle in enumerate(np.deg2rad(theta)):
  90. cos_a, sin_a = np.cos(angle), np.sin(angle)
  91. R = np.array(
  92. [
  93. [cos_a, sin_a, -center * (cos_a + sin_a - 1)],
  94. [-sin_a, cos_a, -center * (cos_a - sin_a - 1)],
  95. [0, 0, 1],
  96. ]
  97. )
  98. rotated = warp(padded_image, R, clip=False)
  99. radon_image[:, i] = rotated.sum(0)
  100. return radon_image
  101. def _sinogram_circle_to_square(sinogram):
  102. diagonal = int(np.ceil(np.sqrt(2) * sinogram.shape[0]))
  103. pad = diagonal - sinogram.shape[0]
  104. old_center = sinogram.shape[0] // 2
  105. new_center = diagonal // 2
  106. pad_before = new_center - old_center
  107. pad_width = ((pad_before, pad - pad_before), (0, 0))
  108. return np.pad(sinogram, pad_width, mode='constant', constant_values=0)
  109. def _get_fourier_filter(size, filter_name):
  110. """Construct the Fourier filter.
  111. This computation lessens artifacts and removes a small bias as
  112. explained in [1], Chap 3. Equation 61.
  113. Parameters
  114. ----------
  115. size : int
  116. filter size. Must be even.
  117. filter_name : str
  118. Filter used in frequency domain filtering. Filters available:
  119. ramp, shepp-logan, cosine, hamming, hann. Assign None to use
  120. no filter.
  121. Returns
  122. -------
  123. fourier_filter: ndarray
  124. The computed Fourier filter.
  125. References
  126. ----------
  127. .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
  128. Imaging", IEEE Press 1988.
  129. """
  130. n = np.concatenate(
  131. (
  132. np.arange(1, size / 2 + 1, 2, dtype=int),
  133. np.arange(size / 2 - 1, 0, -2, dtype=int),
  134. )
  135. )
  136. f = np.zeros(size)
  137. f[0] = 0.25
  138. f[1::2] = -1 / (np.pi * n) ** 2
  139. # Computing the ramp filter from the fourier transform of its
  140. # frequency domain representation lessens artifacts and removes a
  141. # small bias as explained in [1], Chap 3. Equation 61
  142. fourier_filter = 2 * np.real(fft(f)) # ramp filter
  143. if filter_name == "ramp":
  144. pass
  145. elif filter_name == "shepp-logan":
  146. # Start from first element to avoid divide by zero
  147. omega = np.pi * fftfreq(size)[1:]
  148. fourier_filter[1:] *= np.sin(omega) / omega
  149. elif filter_name == "cosine":
  150. freq = np.linspace(0, np.pi, size, endpoint=False)
  151. cosine_filter = fftshift(np.sin(freq))
  152. fourier_filter *= cosine_filter
  153. elif filter_name == "hamming":
  154. fourier_filter *= fftshift(np.hamming(size))
  155. elif filter_name == "hann":
  156. fourier_filter *= fftshift(np.hanning(size))
  157. elif filter_name is None:
  158. fourier_filter[:] = 1
  159. return fourier_filter[:, np.newaxis]
  160. def iradon(
  161. radon_image,
  162. theta=None,
  163. output_size=None,
  164. filter_name="ramp",
  165. interpolation="linear",
  166. circle=True,
  167. preserve_range=True,
  168. ):
  169. """Inverse radon transform.
  170. Reconstruct an image from the radon transform, using the filtered
  171. back projection algorithm.
  172. Parameters
  173. ----------
  174. radon_image : ndarray
  175. Image containing radon transform (sinogram). Each column of
  176. the image corresponds to a projection along a different
  177. angle. The tomography rotation axis should lie at the pixel
  178. index ``radon_image.shape[0] // 2`` along the 0th dimension of
  179. ``radon_image``.
  180. theta : array, optional
  181. Reconstruction angles (in degrees). Default: m angles evenly spaced
  182. between 0 and 180 (if the shape of `radon_image` is (N, M)).
  183. output_size : int, optional
  184. Number of rows and columns in the reconstruction.
  185. filter_name : str, optional
  186. Filter used in frequency domain filtering. Ramp filter used by default.
  187. Filters available: ramp, shepp-logan, cosine, hamming, hann.
  188. Assign None to use no filter.
  189. interpolation : str, optional
  190. Interpolation method used in reconstruction. Methods available:
  191. 'linear', 'nearest', and 'cubic' ('cubic' is slow).
  192. circle : bool, optional
  193. Assume the reconstructed image is zero outside the inscribed circle.
  194. Also changes the default output_size to match the behaviour of
  195. ``radon`` called with ``circle=True``.
  196. preserve_range : bool, optional
  197. Whether to keep the original range of values. Otherwise, the input
  198. image is converted according to the conventions of `img_as_float`.
  199. Also see https://scikit-image.org/docs/dev/user_guide/data_types.html
  200. Returns
  201. -------
  202. reconstructed : ndarray
  203. Reconstructed image. The rotation axis will be located in the pixel
  204. with indices
  205. ``(reconstructed.shape[0] // 2, reconstructed.shape[1] // 2)``.
  206. .. versionchanged:: 0.19
  207. In ``iradon``, ``filter`` argument is deprecated in favor of
  208. ``filter_name``.
  209. References
  210. ----------
  211. .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
  212. Imaging", IEEE Press 1988.
  213. .. [2] B.R. Ramesh, N. Srinivasa, K. Rajgopal, "An Algorithm for Computing
  214. the Discrete Radon Transform With Some Applications", Proceedings of
  215. the Fourth IEEE Region 10 International Conference, TENCON '89, 1989
  216. Notes
  217. -----
  218. It applies the Fourier slice theorem to reconstruct an image by
  219. multiplying the frequency domain of the filter with the FFT of the
  220. projection data. This algorithm is called filtered back projection.
  221. """
  222. if radon_image.ndim != 2:
  223. raise ValueError('The input image must be 2-D')
  224. if theta is None:
  225. theta = np.linspace(0, 180, radon_image.shape[1], endpoint=False)
  226. angles_count = len(theta)
  227. if angles_count != radon_image.shape[1]:
  228. raise ValueError(
  229. "The given ``theta`` does not match the number of "
  230. "projections in ``radon_image``."
  231. )
  232. interpolation_types = ('linear', 'nearest', 'cubic')
  233. if interpolation not in interpolation_types:
  234. raise ValueError(f"Unknown interpolation: {interpolation}")
  235. filter_types = ('ramp', 'shepp-logan', 'cosine', 'hamming', 'hann', None)
  236. if filter_name not in filter_types:
  237. raise ValueError(f"Unknown filter: {filter_name}")
  238. radon_image = convert_to_float(radon_image, preserve_range)
  239. dtype = radon_image.dtype
  240. img_shape = radon_image.shape[0]
  241. if output_size is None:
  242. # If output size not specified, estimate from input radon image
  243. if circle:
  244. output_size = img_shape
  245. else:
  246. output_size = int(np.floor(np.sqrt((img_shape) ** 2 / 2.0)))
  247. if circle:
  248. radon_image = _sinogram_circle_to_square(radon_image)
  249. img_shape = radon_image.shape[0]
  250. # Resize image to next power of two (but no less than 64) for
  251. # Fourier analysis; speeds up Fourier and lessens artifacts
  252. projection_size_padded = max(64, int(2 ** np.ceil(np.log2(2 * img_shape))))
  253. pad_width = ((0, projection_size_padded - img_shape), (0, 0))
  254. img = np.pad(radon_image, pad_width, mode='constant', constant_values=0)
  255. # Apply filter in Fourier domain
  256. fourier_filter = _get_fourier_filter(projection_size_padded, filter_name)
  257. projection = fft(img, axis=0) * fourier_filter
  258. radon_filtered = np.real(ifft(projection, axis=0)[:img_shape, :])
  259. # Reconstruct image by interpolation
  260. reconstructed = np.zeros((output_size, output_size), dtype=dtype)
  261. radius = output_size // 2
  262. xpr, ypr = np.mgrid[:output_size, :output_size] - radius
  263. x = np.arange(img_shape) - img_shape // 2
  264. for col, angle in zip(radon_filtered.T, np.deg2rad(theta)):
  265. t = ypr * np.cos(angle) - xpr * np.sin(angle)
  266. if interpolation == 'linear':
  267. interpolant = partial(np.interp, xp=x, fp=col, left=0, right=0)
  268. else:
  269. interpolant = interp1d(
  270. x, col, kind=interpolation, bounds_error=False, fill_value=0
  271. )
  272. reconstructed += interpolant(t)
  273. if circle:
  274. out_reconstruction_circle = (xpr**2 + ypr**2) > radius**2
  275. reconstructed[out_reconstruction_circle] = 0.0
  276. return reconstructed * np.pi / (2 * angles_count)
  277. def order_angles_golden_ratio(theta):
  278. """Order angles to reduce the amount of correlated information in
  279. subsequent projections.
  280. Parameters
  281. ----------
  282. theta : array of floats, shape (M,)
  283. Projection angles in degrees. Duplicate angles are not allowed.
  284. Returns
  285. -------
  286. indices_generator : generator yielding unsigned integers
  287. The returned generator yields indices into ``theta`` such that
  288. ``theta[indices]`` gives the approximate golden ratio ordering
  289. of the projections. In total, ``len(theta)`` indices are yielded.
  290. All non-negative integers < ``len(theta)`` are yielded exactly once.
  291. Notes
  292. -----
  293. The method used here is that of the golden ratio introduced
  294. by T. Kohler.
  295. References
  296. ----------
  297. .. [1] Kohler, T. "A projection access scheme for iterative
  298. reconstruction based on the golden section." Nuclear Science
  299. Symposium Conference Record, 2004 IEEE. Vol. 6. IEEE, 2004.
  300. .. [2] Winkelmann, Stefanie, et al. "An optimal radial profile order
  301. based on the Golden Ratio for time-resolved MRI."
  302. Medical Imaging, IEEE Transactions on 26.1 (2007): 68-76.
  303. """
  304. interval = 180
  305. remaining_indices = list(np.argsort(theta)) # indices into theta
  306. # yield an arbitrary angle to start things off
  307. angle = theta[remaining_indices[0]]
  308. yield remaining_indices.pop(0)
  309. # determine subsequent angles using the golden ratio method
  310. angle_increment = interval / golden_ratio**2
  311. while remaining_indices:
  312. remaining_angles = theta[remaining_indices]
  313. angle = (angle + angle_increment) % interval
  314. index_above = np.searchsorted(remaining_angles, angle)
  315. index_below = index_above - 1
  316. index_above %= len(remaining_indices)
  317. diff_below = abs(angle - remaining_angles[index_below])
  318. distance_below = min(diff_below % interval, diff_below % -interval)
  319. diff_above = abs(angle - remaining_angles[index_above])
  320. distance_above = min(diff_above % interval, diff_above % -interval)
  321. if distance_below < distance_above:
  322. yield remaining_indices.pop(index_below)
  323. else:
  324. yield remaining_indices.pop(index_above)
  325. def iradon_sart(
  326. radon_image,
  327. theta=None,
  328. image=None,
  329. projection_shifts=None,
  330. clip=None,
  331. relaxation=0.15,
  332. dtype=None,
  333. ):
  334. """Inverse radon transform.
  335. Reconstruct an image from the radon transform, using a single iteration of
  336. the Simultaneous Algebraic Reconstruction Technique (SART) algorithm.
  337. Parameters
  338. ----------
  339. radon_image : ndarray, shape (M, N)
  340. Image containing radon transform (sinogram). Each column of
  341. the image corresponds to a projection along a different angle. The
  342. tomography rotation axis should lie at the pixel index
  343. ``radon_image.shape[0] // 2`` along the 0th dimension of
  344. ``radon_image``.
  345. theta : array, shape (N,), optional
  346. Reconstruction angles (in degrees). Default: m angles evenly spaced
  347. between 0 and 180 (if the shape of `radon_image` is (N, M)).
  348. image : ndarray, shape (M, M), optional
  349. Image containing an initial reconstruction estimate. Default is an array of zeros.
  350. projection_shifts : array, shape (N,), optional
  351. Shift the projections contained in ``radon_image`` (the sinogram) by
  352. this many pixels before reconstructing the image. The i'th value
  353. defines the shift of the i'th column of ``radon_image``.
  354. clip : length-2 sequence of floats, optional
  355. Force all values in the reconstructed tomogram to lie in the range
  356. ``[clip[0], clip[1]]``
  357. relaxation : float, optional
  358. Relaxation parameter for the update step. A higher value can
  359. improve the convergence rate, but one runs the risk of instabilities.
  360. Values close to or higher than 1 are not recommended.
  361. dtype : dtype, optional
  362. Output data type, must be floating point. By default, if input
  363. data type is not float, input is cast to double, otherwise
  364. dtype is set to input data type.
  365. Returns
  366. -------
  367. reconstructed : ndarray
  368. Reconstructed image. The rotation axis will be located in the pixel
  369. with indices
  370. ``(reconstructed.shape[0] // 2, reconstructed.shape[1] // 2)``.
  371. Notes
  372. -----
  373. Algebraic Reconstruction Techniques are based on formulating the tomography
  374. reconstruction problem as a set of linear equations. Along each ray,
  375. the projected value is the sum of all the values of the cross section along
  376. the ray. A typical feature of SART (and a few other variants of algebraic
  377. techniques) is that it samples the cross section at equidistant points
  378. along the ray, using linear interpolation between the pixel values of the
  379. cross section. The resulting set of linear equations are then solved using
  380. a slightly modified Kaczmarz method.
  381. When using SART, a single iteration is usually sufficient to obtain a good
  382. reconstruction. Further iterations will tend to enhance high-frequency
  383. information, but will also often increase the noise.
  384. References
  385. ----------
  386. .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
  387. Imaging", IEEE Press 1988.
  388. .. [2] AH Andersen, AC Kak, "Simultaneous algebraic reconstruction
  389. technique (SART): a superior implementation of the ART algorithm",
  390. Ultrasonic Imaging 6 pp 81--94 (1984)
  391. .. [3] S Kaczmarz, "Angenäherte auflösung von systemen linearer
  392. gleichungen", Bulletin International de l’Academie Polonaise des
  393. Sciences et des Lettres 35 pp 355--357 (1937)
  394. .. [4] Kohler, T. "A projection access scheme for iterative
  395. reconstruction based on the golden section." Nuclear Science
  396. Symposium Conference Record, 2004 IEEE. Vol. 6. IEEE, 2004.
  397. .. [5] Kaczmarz' method, Wikipedia,
  398. https://en.wikipedia.org/wiki/Kaczmarz_method
  399. """
  400. if radon_image.ndim != 2:
  401. raise ValueError('radon_image must be two dimensional')
  402. if dtype is None:
  403. if radon_image.dtype.char in 'fd':
  404. dtype = radon_image.dtype
  405. else:
  406. warn(
  407. "Only floating point data type are valid for SART inverse "
  408. "radon transform. Input data is cast to float. To disable "
  409. "this warning, please cast image_radon to float."
  410. )
  411. dtype = np.dtype(float)
  412. elif np.dtype(dtype).char not in 'fd':
  413. raise ValueError(
  414. "Only floating point data type are valid for inverse " "radon transform."
  415. )
  416. dtype = np.dtype(dtype)
  417. radon_image = radon_image.astype(dtype, copy=False)
  418. reconstructed_shape = (radon_image.shape[0], radon_image.shape[0])
  419. if theta is None:
  420. theta = np.linspace(0, 180, radon_image.shape[1], endpoint=False, dtype=dtype)
  421. elif len(theta) != radon_image.shape[1]:
  422. raise ValueError(
  423. f'Shape of theta ({len(theta)}) does not match the '
  424. f'number of projections ({radon_image.shape[1]})'
  425. )
  426. else:
  427. theta = np.asarray(theta, dtype=dtype)
  428. if image is None:
  429. image = np.zeros(reconstructed_shape, dtype=dtype)
  430. elif image.shape != reconstructed_shape:
  431. raise ValueError(
  432. f'Shape of image ({image.shape}) does not match first dimension '
  433. f'of radon_image ({reconstructed_shape})'
  434. )
  435. elif image.dtype != dtype:
  436. warn(f'image dtype does not match output dtype: ' f'image is cast to {dtype}')
  437. image = np.asarray(image, dtype=dtype)
  438. if projection_shifts is None:
  439. projection_shifts = np.zeros((radon_image.shape[1],), dtype=dtype)
  440. elif len(projection_shifts) != radon_image.shape[1]:
  441. raise ValueError(
  442. f'Shape of projection_shifts ({len(projection_shifts)}) does not match the '
  443. f'number of projections ({radon_image.shape[1]})'
  444. )
  445. else:
  446. projection_shifts = np.asarray(projection_shifts, dtype=dtype)
  447. if clip is not None:
  448. if len(clip) != 2:
  449. raise ValueError('clip must be a length-2 sequence')
  450. clip = np.asarray(clip, dtype=dtype)
  451. for angle_index in order_angles_golden_ratio(theta):
  452. image_update = sart_projection_update(
  453. image,
  454. theta[angle_index],
  455. radon_image[:, angle_index],
  456. projection_shifts[angle_index],
  457. )
  458. image += relaxation * image_update
  459. if clip is not None:
  460. image = np.clip(image, clip[0], clip[1])
  461. return image