tomopy.prep.stripe

Module for pre-processing tasks.

Functions:

remove_stripe_fw(tomo[, level, wname, …])

Remove horizontal stripes from sinogram using the Fourier-Wavelet (FW) based method [B4].

remove_stripe_ti(tomo[, nblock, alpha, …])

Remove horizontal stripes from sinogram using Titarenko’s approach [B12].

remove_stripe_sf(tomo[, size, ncore, nchunk])

Normalize raw projection data using a smoothing filter approach.

remove_stripe_based_sorting(tomo[, size, …])

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 3 in the paper.

remove_stripe_based_filtering(tomo[, sigma, …])

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 2 in the paper.

remove_stripe_based_fitting(tomo[, order, …])

Remove horizontal stripes from sinogram using Nghia Vo’s approach [B23] Algorithm 1 in the paper.

remove_large_stripe(tomo[, snr, size, …])

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 5 in the paper.

remove_dead_stripe(tomo[, snr, size, ncore, …])

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 6 in the paper.

remove_all_stripe(tomo[, snr, la_size, …])

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Combine algorithms 6,5,4,3 to remove all types of stripes.

tomopy.prep.stripe.remove_all_stripe(tomo, snr=3, la_size=61, sm_size=21, ncore=None, nchunk=None)[source]

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Combine algorithms 6,5,4,3 to remove all types of stripes.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • snr (float) – Ratio used to locate large stripes. Greater is less sensitive.

  • la_size (int) – Window size of the median filter to remove large stripes.

  • sm_size (int) – Window size of the median filter to remove small-to-medium stripes.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_dead_stripe(tomo, snr=3, size=51, ncore=None, nchunk=None)[source]

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 6 in the paper. Remove unresponsive and fluctuating stripes.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • snr (float) – Ratio used to detect locations of large stripes. Greater is less sensitive.

  • size (int) – Window size of the median filter.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_large_stripe(tomo, snr=3, size=51, ncore=None, nchunk=None)[source]

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 5 in the paper. Remove large stripes.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • snr (float) – Ratio used to locate of large stripes. Greater is less sensitive.

  • size (int) – Window size of the median filter.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_based_filtering(tomo, sigma=3, size=None, ncore=None, nchunk=None)[source]

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 2 in the paper. Remove stripes using the filtering technique.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • sigma (float) – Sigma of the Gaussian window which is used to separate the low-pass and high-pass components of the intensity profiles of each column.

  • size (int) – Window size of the median filter.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_based_fitting(tomo, order=3, sigma=5, 20, ncore=None, nchunk=None)[source]

Remove horizontal stripes from sinogram using Nghia Vo’s approach [B23] Algorithm 1 in the paper. Remove stripes using the fitting technique.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • order (int) – Polynomial fit order.

  • sigma (tuple of 2 floats) – Sigmas of a 2D Gaussian window in x and y direction.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_based_sorting(tomo, size=None, ncore=None, nchunk=None)[source]

Remove stripe artifacts from sinogram using Nghia Vo’s approach [B23] Algorithm 3 in the paper. Remove stripes using the sorting technique. Work particularly well for removing partial stripes.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • size (int) – Window size of the median filter.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_fw(tomo, level=None, wname='db5', sigma=2, pad=True, ncore=None, nchunk=None)[source]

Remove horizontal stripes from sinogram using the Fourier-Wavelet (FW) based method [B4].

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • level (int, optional) – Number of discrete wavelet transform levels.

  • wname (str, optional) – Type of the wavelet filter. ‘haar’, ‘db5’, sym5’, etc.

  • sigma (float, optional) – Damping parameter in Fourier space.

  • pad (bool, optional) – If True, extend the size of the sinogram by padding with zeros.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_sf(tomo, size=5, ncore=None, nchunk=None)[source]

Normalize raw projection data using a smoothing filter approach.

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • size (int, optional) – Size of the smoothing filter.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.

tomopy.prep.stripe.remove_stripe_ti(tomo, nblock=0, alpha=1.5, ncore=None, nchunk=None)[source]

Remove horizontal stripes from sinogram using Titarenko’s approach [B12].

Parameters
  • tomo (ndarray) – 3D tomographic data.

  • nblock (int, optional) – Number of blocks.

  • alpha (int, optional) – Damping factor.

  • ncore (int, optional) – Number of cores that will be assigned to jobs.

  • nchunk (int, optional) – Chunk size for each core.

Returns

ndarray – Corrected 3D tomographic data.