7.3.1. algotom.rec.reconstruction

Module of reconstruction methods:

  • Filtered back-projection (FBP) method for GPU and CPU.

  • Back-projection filtering (BPF) method for GPU and CPU.

  • Direct Fourier inversion (DFI) method.

  • Wrapper for Astra-Toolbox reconstruction methods (optional)

  • Wrapper for Tomopy-gridrec reconstruction method (optional)

  • Automatic determination of the center of rotation.

  • Tool to assist in manual determination of the center of rotation.

Functions:

fbp_reconstruction(sinogram, center[, ...])

Apply the FBP (filtered back-projection) reconstruction method to a sinogram-image or a chunk of sinogram-images.

bpf_reconstruction(sinogram, center[, ...])

Apply the BPF (back-projection filtering) reconstruction method to a sinogram-image or a chunk of sinogram-images.

dfi_reconstruction(sinogram, center[, ...])

Apply the DFI (direct Fourier inversion) reconstruction method (Ref.

gridrec_reconstruction(sinogram, center[, ...])

Apply the gridrec method to a sinogram-image or a chunk of sinogram-images.

astra_reconstruction(sinogram, center[, ...])

Wrapper of reconstruction methods implemented in the astra toolbox package.

find_center_based_slice_metric(sinogram, ...)

Find the center-of-rotation (COR) using metrics of reconstructed slices at different CORs.

find_center_visual_slices(sinogram, output, ...)

For visually finding the center-of-rotation (COR) using reconstructed slices at different CORs.

algotom.rec.reconstruction.make_smoothing_window(filter_name, width)[source]

Make a 1d smoothing window.

Parameters
  • filter_name ({“hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Window function used for filtering.

  • width (int) – Width of the window.

Returns

array_like – 1D array.

algotom.rec.reconstruction.make_2d_ramp_window(height, width, filter_name=None)[source]

Make the 2d ramp window (in the Fourier space) by repeating the 1d ramp window with the option of adding a smoothing window.

Parameters
  • height (int) – Height of the window.

  • width (int) – Width of the window.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Name of a smoothing window used.

Returns

complex ndarray – 2D array.

algotom.rec.reconstruction.apply_ramp_filter(sinogram, ramp_win=None, filter_name=None, pad=None, pad_mode='edge')[source]

Apply the ramp filter to a sinogram with the option of adding a smoothing filter.

Parameters
  • sinogram (array_like) – 2D array. Sinogram image.

  • ramp_win (complex ndarray or None) – Ramp window in the Fourier space.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Name of a smoothing window used.

  • pad (int or None) – To apply padding before the FFT. The value is set to 10% of the image width if None is given.

  • pad_mode (str) – Padding method. Full list can be found at numpy_pad documentation.

Returns

array_like – Filtered sinogram.

algotom.rec.reconstruction.back_projection_gpu(sinogram, angles, center, block=(16, 16), edge_pad=False)[source]

Implement the back-projection algorithm using GPU.

Parameters
  • sinogram (array_like) – 2D array. Sinogram image.

  • angles (array_like) – 1D array. Angles (radian) corresponding to the sinogram.

  • center (float) – Center of rotation.

  • edge_pad (bool) – Enable/disable edge padding.

  • block (tuple of int, optional) – Size of a GPU block. E.g. (8, 8), (16, 16), (32, 32), …

Returns

recon (array_like) – Back-projected image.

algotom.rec.reconstruction.back_projection_gpu_chunk(sinograms, angles, center, block=(16, 16), edge_pad=False)[source]

Implement the back-projection algorithm for a chunk of sinograms using GPU. Axis of a sinogram/slice in the 3D array is 1.

Parameters
  • sinograms (array_like) – 3D array. Sinogram images.

  • angles (array_like) – 1D array. Angles (radian) corresponding to a sinogram.

  • center (float) – Center of rotation.

  • edge_pad (bool) – Enable/disable edge padding.

  • block (tuple of int, optional) – Size of a GPU block. E.g. (8, 8), (16, 16), (32, 32), …

Returns

recons (array_like) – Back-projected images.

algotom.rec.reconstruction.back_projection_cpu(sinogram, angles, center, edge_pad=False)

Implement the back-projection algorithm using CPU.

Parameters
  • sinogram (array_like) – 2D array. (Filtered) sinogram image.

  • angles (array_like) – 1D array. Angles (radian) corresponding to the sinogram.

  • center (float) – Center of rotation.

  • edge_pad (bool) – Enable/disable edge padding.

Returns

recon (array_like) – Square array, back-projected image.

algotom.rec.reconstruction.fbp_reconstruction(sinogram, center, angles=None, ratio=1.0, ramp_win=None, filter_name='hann', pad=None, pad_mode='edge', apply_log=True, gpu=True, block=(16, 16), ncore=None)[source]

Apply the FBP (filtered back-projection) reconstruction method to a sinogram-image or a chunk of sinogram-images. Angular axis is 0. If input is 3D array, the slicing axis of sinograms must be 1, e.g. data[:, index, :].

Parameters
  • sinogram (array_like) – 2D/3D array. Sinogram image.

  • center (float) – Center of rotation.

  • angles (array_like, optional) – 1D array. List of angles (in radian) corresponding to the sinogram.

  • ratio (float, optional) – To apply a circle mask to the reconstructed image.

  • ramp_win (complex ndarray, optional) – Ramp window in the Fourier space. Generated if None.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Apply a smoothing filter.

  • pad (int, optional) – To apply padding before the FFT. The value is set to 10% of the image width if None is given.

  • pad_mode (str, optional) – Padding method. Full list can be found at numpy_pad documentation.

  • apply_log (bool, optional) – Apply the logarithm function to the sinogram before reconstruction.

  • gpu (bool, optional) – Use GPU for computing if True.

  • block (tuple of two integer-values, optional) – Size of a GPU block. E.g. (8, 8), (16, 16), (32, 32), …

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

Returns

array_like – Square array. Reconstructed image.

algotom.rec.reconstruction.make_circular_ramp_window(width, filter_name=None)[source]

Make a circular ramp window (2d) with the option of adding a smoothing window.

Parameters
  • width (int) – Width of the window.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Name of a smoothing window used.

Returns

array_like – Square array, size of (width, width)

algotom.rec.reconstruction.apply_circular_ramp_filter(rec_img, ramp_win=None, filter_name=None, pad=None, pad_mode='edge')[source]

Apply the circular ramp filter to a back-projected image.

Parameters
  • rec_img (array_like) – Square array. back-projected image.

  • ramp_win (array_like) – 2d circular ramp window, generated if None given.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Name of a smoothing window used.

  • pad (int or None) – To apply padding before the FFT. The value is set to 10% of the image width if None is given.

  • pad_mode (str) – Padding method. Full list can be found at numpy_pad documentation.

Returns

array_like – Square array.

algotom.rec.reconstruction.bpf_reconstruction(sinogram, center, angles=None, ratio=1.0, ramp_win=None, filter_name='hann', pad=None, pad_mode='edge', apply_log=True, gpu=True, block=(16, 16), ncore=None)[source]

Apply the BPF (back-projection filtering) reconstruction method to a sinogram-image or a chunk of sinogram-images. Angular axis is 0. If input is 3D array, the slicing axis of sinograms must be 1, e.g. data[:, index, :].

Parameters
  • sinogram (array_like) – 2D/3D array. Sinogram image.

  • center (float) – Center of rotation.

  • angles (array_like, optional) – 1D array. List of angles (in radian) corresponding to the sinogram.

  • ratio (float, optional) – Apply a circle mask to the reconstructed image.

  • ramp_win (complex ndarray, optional) – Circular ramp window, generated if None.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Apply a smoothing filter.

  • pad (int, optional) – Apply padding before the FFT. The value is set to 10% of the image width if None is given.

  • pad_mode (str, optional) – Padding method. Full list can be found at numpy_pad documentation.

  • apply_log (bool, optional) – Apply logarithm to sinogram before reconstruction.

  • gpu (bool, optional) – Use GPU for computing if True.

  • block (tuple of two integer-values, optional) – Size of a GPU block. E.g. (8, 8), (16, 16), (32, 32), …

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

Returns

array_like – Square array. Reconstructed image.

algotom.rec.reconstruction.generate_mapping_coordinate(width_sino, height_sino, width_rec, height_rec)[source]

Calculate coordinates in the sinogram space from coordinates in the reconstruction space (in the Fourier domain). They are used for the DFI (direct Fourier inversion) reconstruction method.

Parameters
  • width_sino (int) – Width of a sinogram image.

  • height_sino (int) – Height of a sinogram image.

  • width_rec (int) – Width of a reconstruction image.

  • height_rec (int) – Height of a reconstruction image.

Returns

  • r_mat (array_like) – 2D array. Broadcast of the r-coordinates.

  • theta_mat (array_like) – 2D array. Broadcast of the theta-coordinates.

algotom.rec.reconstruction.dfi_reconstruction(sinogram, center, angles=None, ratio=1.0, filter_name='hann', pad_rate=0.25, pad_mode='edge', apply_log=True, ncore=None)[source]

Apply the DFI (direct Fourier inversion) reconstruction method (Ref. [1]) to a sinogram-image or a chunk of sinogram-images. Angular axis is 0. If input is 3D array, the slicing axis of sinograms must be 1, e.g. data[:, index, :].

Parameters
  • sinogram (array_like) – 2D/3D array. Sinogram image.

  • center (float) – Center of rotation.

  • angles (array_like) – 1D array. List of angles (in radian) corresponding to the sinogram.

  • ratio (float) – To apply a circle mask to the reconstructed image.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Apply a smoothing filter.

  • pad_rate (float) – To apply padding before the FFT. The padding width equals to (pad_rate * image_width).

  • pad_mode (str) – Padding method. Full list can be found at numpy_pad documentation.

  • apply_log (bool) – Apply the logarithm function to the sinogram before reconstruction.

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

Returns

array_like – Square array. Reconstructed image.

References

[1] : https://doi.org/10.1071/PH560198

algotom.rec.reconstruction.gridrec_reconstruction(sinogram, center, angles=None, ratio=1.0, filter_name='shepp', apply_log=True, pad=100, filter_par=0.9, ncore=1)[source]

Apply the gridrec method to a sinogram-image or a chunk of sinogram-images. Angular axis is 0. If input is 3D array, the slicing axis of sinograms must be 1, e.g. data[:, index, :]. This is the wrapper of the gridrec method implemented in the Tomopy package: https://tomopy.readthedocs.io/en/latest/api/tomopy.recon.algorithm.html. Users must install Tomopy before using this function.

Parameters
  • sinogram (array_like) – 2D/3D array. Sinogram image.

  • center (float) – Center of rotation.

  • angles (array_like) – 1D array. List of angles (radian) corresponding to the sinogram.

  • ratio (float) – To apply a circle mask to the reconstructed image.

  • filter_name (str or None) – Apply a smoothing filter. Full list is at: https://github.com/tomopy/tomopy/blob/master/source/tomopy/recon/algorithm.py

  • filter_par (float) – Adjust the strength of the filter. Smaller is stronger.

  • apply_log (bool) – Apply the logarithm function to the sinogram before reconstruction.

  • pad (bool or int) – Apply edge padding to the nearest power of 2.

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

Returns

array_like – Square array.

algotom.rec.reconstruction.astra_reconstruction(sinogram, center, angles=None, ratio=1.0, method='FBP_CUDA', num_iter=1, filter_name='hann', pad=None, apply_log=True, ncore=1)[source]

Wrapper of reconstruction methods implemented in the astra toolbox package. https://www.astra-toolbox.com/docs/algs/index.html Users must install Astra Toolbox before using this function. Apply the method to a sinogram-image or a chunk of sinogram-images. Angular axis is 0. If input is 3D array, the slicing axis of sinograms must be 1, e.g. data[:, index, :]

Parameters
  • sinogram (array_like) – 2D/3D array. Sinogram image.

  • center (float) – Center of rotation.

  • angles (array_like) – 1D array. List of angles (radian) corresponding to the sinogram.

  • ratio (float) – To apply a circle mask to the reconstructed image.

  • method (str) – Reconstruction algorithms. For CPU: ‘FBP’, ‘SIRT’, ‘SART’, ‘ART’, and ‘CGLS’. For GPU: ‘FBP_CUDA’, ‘SIRT_CUDA’, ‘SART_CUDA’, and ‘CGLS_CUDA’.

  • num_iter (int) – Number of iterations if using iteration methods.

  • filter_name (str or None) – Apply filter if using FBP method. Options: ‘ram-lak’, ‘hamming’, ‘hann’, ‘lanczos’, ‘kaiser’, ‘parzen’,…

  • pad (int) – Padding to reduce the side effect of FFT.

  • apply_log (bool) – Apply the logarithm function to the sinogram before reconstruction.

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

Returns

array_like – Square array.

algotom.rec.reconstruction.find_center_based_slice_metric(sinogram, start, stop, step=0.5, metric='entropy', radius=2, zoom=1.0, method='fbp', gpu=True, angles=None, ratio=1.0, filter_name='hann', apply_log=True, ncore=None, sigma=2, invert_metric=False, metric_function=None, **kwargs)[source]

Find the center-of-rotation (COR) using metrics of reconstructed slices at different CORs. The entropy of histogram (Ref. [1]) is used by default. If customized metrics are used, the minimum value must be corresponding to the optimal center.

Parameters
  • sinogram (array_like) – 2D array. Sinogram image.

  • start (float) – Starting point for searching CoR.

  • stop (float) – Ending point for searching CoR.

  • step (float) – Sub-pixel searching step.

  • metric ({“entropy”, “sharpness”}) – Which metric to use.

  • radius (float) – Searching range with the sub-pixel step.

  • zoom (float) – To resize the sinogram for fast coarse-searching. For example, 0.5 <=> reduce the size of the image by half.

  • method ({“bpf”, “dfi”, “gridrec”, “fbp”, “astra”}) – To select a backend method for reconstruction.

  • gpu (bool, optional) – Use GPU for computing if True.

  • angles (array_like, optional) – 1D array. List of angles (in radian) corresponding to the sinogram.

  • ratio (float, optional) – To apply a circle mask to the reconstructed image.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Apply a smoothing filter before reconstruction.

  • apply_log (bool, optional) – Apply the logarithm function to the sinogram before reconstruction.

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

  • sigma (int) – Denoising the sinogram before reconstruction.

  • invert_metric (bool) – Invert the metric scale, used with a custom metric-function.

  • metric_function (obj) – Custom function to calculate metric, accepts keyword arguments (** kwargs).

Returns

float – Center-of-rotation.

References

[1] : https://doi.org/10.1364/JOSAA.23.001048

algotom.rec.reconstruction.find_center_visual_slices(sinogram, output, start, stop, step=1, zoom=0.5, method='fbp', gpu=False, angles=None, ratio=1.0, filter_name='hann', apply_log=True, ncore=None, display=False)[source]

For visually finding the center-of-rotation (COR) using reconstructed slices at different CORs.

Parameters
  • sinogram (array_like) – 2D array. Sinogram image.

  • output (str) – Base folder for saving reconstructed slices.

  • start (float) – Starting point for searching CoR.

  • stop (float) – Ending point for searching CoR.

  • step (float) – Searching step.

  • zoom (float) – To resize input and output images. For example, 0.5 <=> reduce the size of images by half.

  • method ({“bpf”, “dfi”, “gridrec”, “fbp”, “astra”}) – To select a backend method for reconstruction.

  • gpu (bool, optional) – Use GPU for computing if True.

  • angles (array_like, optional) – 1D array. List of angles (in radian) corresponding to the sinogram.

  • ratio (float, optional) – To apply a circle mask to the reconstructed image.

  • filter_name ({None, “hann”, “bartlett”, “blackman”, “hamming”, “nuttall”, “parzen”, “triang”}) – Apply a smoothing filter.

  • apply_log (bool, optional) – Apply the logarithm function to the sinogram before reconstruction.

  • ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.

  • display (bool) – Print the output if True.

Returns

str – Folder path to tif images.