nabu.reconstruction.rings module

class nabu.reconstruction.rings.MunchDeringer(sigma, sinos_shape, levels=None, wname='db15', padding=None, padding_mode='edge')[source]

Bases: object

Initialize a “Munch Et Al” sinogram deringer. See References for more information.

Parameters:
  • sigma (float) – Standard deviation of the damping parameter. The higher value of sigma, the more important the filtering effect on the rings.

  • levels (int, optional) – Number of wavelets decomposition levels. By default (None), the maximum number of decomposition levels is used.

  • wname (str, optional) – Default is “db15” (Daubechies, 15 vanishing moments)

  • sinos_shape (tuple, optional) – Shape of the sinogram (or sinograms stack).

  • padding (tuple of two int, optional) – Horizontal padding to use for reducing the aliasing artefacts

References

B. Munch, P. Trtik, F. Marone, M. Stampanoni, Stripe and ring artifact removal with combined wavelet-Fourier filtering, Optics Express 17(10):8567-8591, 2009.

remove_rings(sinos, output=None)[source]

Main function to performs rings artefacts removal on sinogram(s). CAUTION: this function defaults to in-place processing, meaning that the sinogram(s) you pass will be overwritten.

Parameters:
  • sinos (numpy.ndarray) – Sinogram or stack of sinograms.

  • output (numpy.ndarray, optional) – Output array. If set to None (default), the output overwrites the input.

class nabu.reconstruction.rings.VoDeringer(sinos_shape, **remove_all_stripe_options)[source]

Bases: object

An interface to Nghia Vo’s “remove_all_stripe”. Needs algotom to run.

remove_rings_sinogram(sino, output=None)[source]
remove_rings_sinograms(sinos, output=None)[source]
remove_rings_radios(radios)[source]
remove_rings(sinos, output=None)
class nabu.reconstruction.rings.SinoMeanDeringer(sinos_shape, mode='subtract', filter_cutoff=None, padding_mode='edge', fft_num_threads=None)[source]

Bases: object

Rings correction with mean subtraction/division. The principle of this method is to subtract (or divide) the sinogram by its mean along a certain axis. In short:

sinogram -= filt(sinogram.mean(axis=0))

where filt is some bandpass filter.

Parameters:
  • sinos_shape (tuple of int) – Sinograms shape, in the form (n_angles, n_x) or (n_sinos, n_angles, n_x)

  • mode (str, optional) – Operation to do on the sinogram, either “subtract” or “divide”

  • filter_cutoff (tuple, optional) –

    Cut-off of the bandpass filter applied on the sinogram profiles. Empty (default) means no filtering. Possible values forms are:

    • (sigma_low, sigma_high): two float values defining the standard deviation of gaussian(sigma_low) * (1 - gaussian(sigma_high)). High values of sigma mean stronger effect of associated filters.

    • ((cutoff_low, transition_low), (cutoff_high, transition_high)) where “cutoff” is in normalized Nyquist frequency (0.5 is the maximum frequency), and “transition” is the width of filter decay in fraction of the cutoff frequency

  • padding_mode (str, optional) – Padding mode when filtering the sinogram profile. Should be “constant” (i.e “zeros”) for mathematical correctness, but in practice this yields a Gibbs effect when replicating the sinogram, so “edges” is recommended.

  • fft_num_threads (int, optional) – How many threads to use for computing the fast Fourier transform when filtering the sinogram profile. Defaut is all the available threads.

supported_modes = ['subtract', 'divide']
remove_rings_sinogram(sino, output=None)[source]
remove_rings_sinograms(sinos, output=None)[source]
remove_rings(sinos, output=None)