| 123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536 |
- import numpy as np
- from scipy.interpolate import interp1d
- from scipy.constants import golden_ratio
- from scipy.fft import fft, ifft, fftfreq, fftshift
- from ._warps import warp
- from ._radon_transform import sart_projection_update
- from .._shared.utils import convert_to_float
- from warnings import warn
- from functools import partial
- __all__ = ['radon', 'order_angles_golden_ratio', 'iradon', 'iradon_sart']
- def radon(image, theta=None, circle=True, *, preserve_range=False):
- """
- Calculates the radon transform of an image given specified
- projection angles.
- Parameters
- ----------
- image : ndarray
- Input image. The rotation axis will be located in the pixel with
- indices ``(image.shape[0] // 2, image.shape[1] // 2)``.
- theta : array, optional
- Projection angles (in degrees). If `None`, the value is set to
- np.arange(180).
- circle : bool, optional
- Assume image is zero outside the inscribed circle, making the
- width of each projection (the first dimension of the sinogram)
- equal to ``min(image.shape)``.
- preserve_range : bool, optional
- Whether to keep the original range of values. Otherwise, the input
- image is converted according to the conventions of `img_as_float`.
- Also see https://scikit-image.org/docs/dev/user_guide/data_types.html
- Returns
- -------
- radon_image : ndarray
- Radon transform (sinogram). The tomography rotation axis will lie
- at the pixel index ``radon_image.shape[0] // 2`` along the 0th
- dimension of ``radon_image``.
- References
- ----------
- .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
- Imaging", IEEE Press 1988.
- .. [2] B.R. Ramesh, N. Srinivasa, K. Rajgopal, "An Algorithm for Computing
- the Discrete Radon Transform With Some Applications", Proceedings of
- the Fourth IEEE Region 10 International Conference, TENCON '89, 1989
- Notes
- -----
- Based on code of Justin K. Romberg
- (https://www.clear.rice.edu/elec431/projects96/DSP/bpanalysis.html)
- """
- if image.ndim != 2:
- raise ValueError('The input image must be 2-D')
- if theta is None:
- theta = np.arange(180)
- image = convert_to_float(image, preserve_range)
- if circle:
- shape_min = min(image.shape)
- radius = shape_min // 2
- img_shape = np.array(image.shape)
- coords = np.array(np.ogrid[: image.shape[0], : image.shape[1]], dtype=object)
- dist = ((coords - img_shape // 2) ** 2).sum(0)
- outside_reconstruction_circle = dist > radius**2
- if np.any(image[outside_reconstruction_circle]):
- warn(
- 'Radon transform: image must be zero outside the '
- 'reconstruction circle'
- )
- # Crop image to make it square
- slices = tuple(
- (
- slice(int(np.ceil(excess / 2)), int(np.ceil(excess / 2) + shape_min))
- if excess > 0
- else slice(None)
- )
- for excess in (img_shape - shape_min)
- )
- padded_image = image[slices]
- else:
- diagonal = np.sqrt(2) * max(image.shape)
- pad = [int(np.ceil(diagonal - s)) for s in image.shape]
- new_center = [(s + p) // 2 for s, p in zip(image.shape, pad)]
- old_center = [s // 2 for s in image.shape]
- pad_before = [nc - oc for oc, nc in zip(old_center, new_center)]
- pad_width = [(pb, p - pb) for pb, p in zip(pad_before, pad)]
- padded_image = np.pad(image, pad_width, mode='constant', constant_values=0)
- # padded_image is always square
- if padded_image.shape[0] != padded_image.shape[1]:
- raise ValueError('padded_image must be a square')
- center = padded_image.shape[0] // 2
- radon_image = np.zeros((padded_image.shape[0], len(theta)), dtype=image.dtype)
- for i, angle in enumerate(np.deg2rad(theta)):
- cos_a, sin_a = np.cos(angle), np.sin(angle)
- R = np.array(
- [
- [cos_a, sin_a, -center * (cos_a + sin_a - 1)],
- [-sin_a, cos_a, -center * (cos_a - sin_a - 1)],
- [0, 0, 1],
- ]
- )
- rotated = warp(padded_image, R, clip=False)
- radon_image[:, i] = rotated.sum(0)
- return radon_image
- def _sinogram_circle_to_square(sinogram):
- diagonal = int(np.ceil(np.sqrt(2) * sinogram.shape[0]))
- pad = diagonal - sinogram.shape[0]
- old_center = sinogram.shape[0] // 2
- new_center = diagonal // 2
- pad_before = new_center - old_center
- pad_width = ((pad_before, pad - pad_before), (0, 0))
- return np.pad(sinogram, pad_width, mode='constant', constant_values=0)
- def _get_fourier_filter(size, filter_name):
- """Construct the Fourier filter.
- This computation lessens artifacts and removes a small bias as
- explained in [1], Chap 3. Equation 61.
- Parameters
- ----------
- size : int
- filter size. Must be even.
- filter_name : str
- Filter used in frequency domain filtering. Filters available:
- ramp, shepp-logan, cosine, hamming, hann. Assign None to use
- no filter.
- Returns
- -------
- fourier_filter: ndarray
- The computed Fourier filter.
- References
- ----------
- .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
- Imaging", IEEE Press 1988.
- """
- n = np.concatenate(
- (
- np.arange(1, size / 2 + 1, 2, dtype=int),
- np.arange(size / 2 - 1, 0, -2, dtype=int),
- )
- )
- f = np.zeros(size)
- f[0] = 0.25
- f[1::2] = -1 / (np.pi * n) ** 2
- # Computing the ramp filter from the fourier transform of its
- # frequency domain representation lessens artifacts and removes a
- # small bias as explained in [1], Chap 3. Equation 61
- fourier_filter = 2 * np.real(fft(f)) # ramp filter
- if filter_name == "ramp":
- pass
- elif filter_name == "shepp-logan":
- # Start from first element to avoid divide by zero
- omega = np.pi * fftfreq(size)[1:]
- fourier_filter[1:] *= np.sin(omega) / omega
- elif filter_name == "cosine":
- freq = np.linspace(0, np.pi, size, endpoint=False)
- cosine_filter = fftshift(np.sin(freq))
- fourier_filter *= cosine_filter
- elif filter_name == "hamming":
- fourier_filter *= fftshift(np.hamming(size))
- elif filter_name == "hann":
- fourier_filter *= fftshift(np.hanning(size))
- elif filter_name is None:
- fourier_filter[:] = 1
- return fourier_filter[:, np.newaxis]
- def iradon(
- radon_image,
- theta=None,
- output_size=None,
- filter_name="ramp",
- interpolation="linear",
- circle=True,
- preserve_range=True,
- ):
- """Inverse radon transform.
- Reconstruct an image from the radon transform, using the filtered
- back projection algorithm.
- Parameters
- ----------
- radon_image : ndarray
- Image containing radon transform (sinogram). Each column of
- the image corresponds to a projection along a different
- angle. The tomography rotation axis should lie at the pixel
- index ``radon_image.shape[0] // 2`` along the 0th dimension of
- ``radon_image``.
- theta : array, optional
- Reconstruction angles (in degrees). Default: m angles evenly spaced
- between 0 and 180 (if the shape of `radon_image` is (N, M)).
- output_size : int, optional
- Number of rows and columns in the reconstruction.
- filter_name : str, optional
- Filter used in frequency domain filtering. Ramp filter used by default.
- Filters available: ramp, shepp-logan, cosine, hamming, hann.
- Assign None to use no filter.
- interpolation : str, optional
- Interpolation method used in reconstruction. Methods available:
- 'linear', 'nearest', and 'cubic' ('cubic' is slow).
- circle : bool, optional
- Assume the reconstructed image is zero outside the inscribed circle.
- Also changes the default output_size to match the behaviour of
- ``radon`` called with ``circle=True``.
- preserve_range : bool, optional
- Whether to keep the original range of values. Otherwise, the input
- image is converted according to the conventions of `img_as_float`.
- Also see https://scikit-image.org/docs/dev/user_guide/data_types.html
- Returns
- -------
- reconstructed : ndarray
- Reconstructed image. The rotation axis will be located in the pixel
- with indices
- ``(reconstructed.shape[0] // 2, reconstructed.shape[1] // 2)``.
- .. versionchanged:: 0.19
- In ``iradon``, ``filter`` argument is deprecated in favor of
- ``filter_name``.
- References
- ----------
- .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
- Imaging", IEEE Press 1988.
- .. [2] B.R. Ramesh, N. Srinivasa, K. Rajgopal, "An Algorithm for Computing
- the Discrete Radon Transform With Some Applications", Proceedings of
- the Fourth IEEE Region 10 International Conference, TENCON '89, 1989
- Notes
- -----
- It applies the Fourier slice theorem to reconstruct an image by
- multiplying the frequency domain of the filter with the FFT of the
- projection data. This algorithm is called filtered back projection.
- """
- if radon_image.ndim != 2:
- raise ValueError('The input image must be 2-D')
- if theta is None:
- theta = np.linspace(0, 180, radon_image.shape[1], endpoint=False)
- angles_count = len(theta)
- if angles_count != radon_image.shape[1]:
- raise ValueError(
- "The given ``theta`` does not match the number of "
- "projections in ``radon_image``."
- )
- interpolation_types = ('linear', 'nearest', 'cubic')
- if interpolation not in interpolation_types:
- raise ValueError(f"Unknown interpolation: {interpolation}")
- filter_types = ('ramp', 'shepp-logan', 'cosine', 'hamming', 'hann', None)
- if filter_name not in filter_types:
- raise ValueError(f"Unknown filter: {filter_name}")
- radon_image = convert_to_float(radon_image, preserve_range)
- dtype = radon_image.dtype
- img_shape = radon_image.shape[0]
- if output_size is None:
- # If output size not specified, estimate from input radon image
- if circle:
- output_size = img_shape
- else:
- output_size = int(np.floor(np.sqrt((img_shape) ** 2 / 2.0)))
- if circle:
- radon_image = _sinogram_circle_to_square(radon_image)
- img_shape = radon_image.shape[0]
- # Resize image to next power of two (but no less than 64) for
- # Fourier analysis; speeds up Fourier and lessens artifacts
- projection_size_padded = max(64, int(2 ** np.ceil(np.log2(2 * img_shape))))
- pad_width = ((0, projection_size_padded - img_shape), (0, 0))
- img = np.pad(radon_image, pad_width, mode='constant', constant_values=0)
- # Apply filter in Fourier domain
- fourier_filter = _get_fourier_filter(projection_size_padded, filter_name)
- projection = fft(img, axis=0) * fourier_filter
- radon_filtered = np.real(ifft(projection, axis=0)[:img_shape, :])
- # Reconstruct image by interpolation
- reconstructed = np.zeros((output_size, output_size), dtype=dtype)
- radius = output_size // 2
- xpr, ypr = np.mgrid[:output_size, :output_size] - radius
- x = np.arange(img_shape) - img_shape // 2
- for col, angle in zip(radon_filtered.T, np.deg2rad(theta)):
- t = ypr * np.cos(angle) - xpr * np.sin(angle)
- if interpolation == 'linear':
- interpolant = partial(np.interp, xp=x, fp=col, left=0, right=0)
- else:
- interpolant = interp1d(
- x, col, kind=interpolation, bounds_error=False, fill_value=0
- )
- reconstructed += interpolant(t)
- if circle:
- out_reconstruction_circle = (xpr**2 + ypr**2) > radius**2
- reconstructed[out_reconstruction_circle] = 0.0
- return reconstructed * np.pi / (2 * angles_count)
- def order_angles_golden_ratio(theta):
- """Order angles to reduce the amount of correlated information in
- subsequent projections.
- Parameters
- ----------
- theta : array of floats, shape (M,)
- Projection angles in degrees. Duplicate angles are not allowed.
- Returns
- -------
- indices_generator : generator yielding unsigned integers
- The returned generator yields indices into ``theta`` such that
- ``theta[indices]`` gives the approximate golden ratio ordering
- of the projections. In total, ``len(theta)`` indices are yielded.
- All non-negative integers < ``len(theta)`` are yielded exactly once.
- Notes
- -----
- The method used here is that of the golden ratio introduced
- by T. Kohler.
- References
- ----------
- .. [1] Kohler, T. "A projection access scheme for iterative
- reconstruction based on the golden section." Nuclear Science
- Symposium Conference Record, 2004 IEEE. Vol. 6. IEEE, 2004.
- .. [2] Winkelmann, Stefanie, et al. "An optimal radial profile order
- based on the Golden Ratio for time-resolved MRI."
- Medical Imaging, IEEE Transactions on 26.1 (2007): 68-76.
- """
- interval = 180
- remaining_indices = list(np.argsort(theta)) # indices into theta
- # yield an arbitrary angle to start things off
- angle = theta[remaining_indices[0]]
- yield remaining_indices.pop(0)
- # determine subsequent angles using the golden ratio method
- angle_increment = interval / golden_ratio**2
- while remaining_indices:
- remaining_angles = theta[remaining_indices]
- angle = (angle + angle_increment) % interval
- index_above = np.searchsorted(remaining_angles, angle)
- index_below = index_above - 1
- index_above %= len(remaining_indices)
- diff_below = abs(angle - remaining_angles[index_below])
- distance_below = min(diff_below % interval, diff_below % -interval)
- diff_above = abs(angle - remaining_angles[index_above])
- distance_above = min(diff_above % interval, diff_above % -interval)
- if distance_below < distance_above:
- yield remaining_indices.pop(index_below)
- else:
- yield remaining_indices.pop(index_above)
- def iradon_sart(
- radon_image,
- theta=None,
- image=None,
- projection_shifts=None,
- clip=None,
- relaxation=0.15,
- dtype=None,
- ):
- """Inverse radon transform.
- Reconstruct an image from the radon transform, using a single iteration of
- the Simultaneous Algebraic Reconstruction Technique (SART) algorithm.
- Parameters
- ----------
- radon_image : ndarray, shape (M, N)
- Image containing radon transform (sinogram). Each column of
- the image corresponds to a projection along a different angle. The
- tomography rotation axis should lie at the pixel index
- ``radon_image.shape[0] // 2`` along the 0th dimension of
- ``radon_image``.
- theta : array, shape (N,), optional
- Reconstruction angles (in degrees). Default: m angles evenly spaced
- between 0 and 180 (if the shape of `radon_image` is (N, M)).
- image : ndarray, shape (M, M), optional
- Image containing an initial reconstruction estimate. Default is an array of zeros.
- projection_shifts : array, shape (N,), optional
- Shift the projections contained in ``radon_image`` (the sinogram) by
- this many pixels before reconstructing the image. The i'th value
- defines the shift of the i'th column of ``radon_image``.
- clip : length-2 sequence of floats, optional
- Force all values in the reconstructed tomogram to lie in the range
- ``[clip[0], clip[1]]``
- relaxation : float, optional
- Relaxation parameter for the update step. A higher value can
- improve the convergence rate, but one runs the risk of instabilities.
- Values close to or higher than 1 are not recommended.
- dtype : dtype, optional
- Output data type, must be floating point. By default, if input
- data type is not float, input is cast to double, otherwise
- dtype is set to input data type.
- Returns
- -------
- reconstructed : ndarray
- Reconstructed image. The rotation axis will be located in the pixel
- with indices
- ``(reconstructed.shape[0] // 2, reconstructed.shape[1] // 2)``.
- Notes
- -----
- Algebraic Reconstruction Techniques are based on formulating the tomography
- reconstruction problem as a set of linear equations. Along each ray,
- the projected value is the sum of all the values of the cross section along
- the ray. A typical feature of SART (and a few other variants of algebraic
- techniques) is that it samples the cross section at equidistant points
- along the ray, using linear interpolation between the pixel values of the
- cross section. The resulting set of linear equations are then solved using
- a slightly modified Kaczmarz method.
- When using SART, a single iteration is usually sufficient to obtain a good
- reconstruction. Further iterations will tend to enhance high-frequency
- information, but will also often increase the noise.
- References
- ----------
- .. [1] AC Kak, M Slaney, "Principles of Computerized Tomographic
- Imaging", IEEE Press 1988.
- .. [2] AH Andersen, AC Kak, "Simultaneous algebraic reconstruction
- technique (SART): a superior implementation of the ART algorithm",
- Ultrasonic Imaging 6 pp 81--94 (1984)
- .. [3] S Kaczmarz, "Angenäherte auflösung von systemen linearer
- gleichungen", Bulletin International de l’Academie Polonaise des
- Sciences et des Lettres 35 pp 355--357 (1937)
- .. [4] Kohler, T. "A projection access scheme for iterative
- reconstruction based on the golden section." Nuclear Science
- Symposium Conference Record, 2004 IEEE. Vol. 6. IEEE, 2004.
- .. [5] Kaczmarz' method, Wikipedia,
- https://en.wikipedia.org/wiki/Kaczmarz_method
- """
- if radon_image.ndim != 2:
- raise ValueError('radon_image must be two dimensional')
- if dtype is None:
- if radon_image.dtype.char in 'fd':
- dtype = radon_image.dtype
- else:
- warn(
- "Only floating point data type are valid for SART inverse "
- "radon transform. Input data is cast to float. To disable "
- "this warning, please cast image_radon to float."
- )
- dtype = np.dtype(float)
- elif np.dtype(dtype).char not in 'fd':
- raise ValueError(
- "Only floating point data type are valid for inverse " "radon transform."
- )
- dtype = np.dtype(dtype)
- radon_image = radon_image.astype(dtype, copy=False)
- reconstructed_shape = (radon_image.shape[0], radon_image.shape[0])
- if theta is None:
- theta = np.linspace(0, 180, radon_image.shape[1], endpoint=False, dtype=dtype)
- elif len(theta) != radon_image.shape[1]:
- raise ValueError(
- f'Shape of theta ({len(theta)}) does not match the '
- f'number of projections ({radon_image.shape[1]})'
- )
- else:
- theta = np.asarray(theta, dtype=dtype)
- if image is None:
- image = np.zeros(reconstructed_shape, dtype=dtype)
- elif image.shape != reconstructed_shape:
- raise ValueError(
- f'Shape of image ({image.shape}) does not match first dimension '
- f'of radon_image ({reconstructed_shape})'
- )
- elif image.dtype != dtype:
- warn(f'image dtype does not match output dtype: ' f'image is cast to {dtype}')
- image = np.asarray(image, dtype=dtype)
- if projection_shifts is None:
- projection_shifts = np.zeros((radon_image.shape[1],), dtype=dtype)
- elif len(projection_shifts) != radon_image.shape[1]:
- raise ValueError(
- f'Shape of projection_shifts ({len(projection_shifts)}) does not match the '
- f'number of projections ({radon_image.shape[1]})'
- )
- else:
- projection_shifts = np.asarray(projection_shifts, dtype=dtype)
- if clip is not None:
- if len(clip) != 2:
- raise ValueError('clip must be a length-2 sequence')
- clip = np.asarray(clip, dtype=dtype)
- for angle_index in order_angles_golden_ratio(theta):
- image_update = sart_projection_update(
- image,
- theta[angle_index],
- radon_image[:, angle_index],
- projection_shifts[angle_index],
- )
- image += relaxation * image_update
- if clip is not None:
- image = np.clip(image, clip[0], clip[1])
- return image
|