# -*- coding: utf-8 -*- # deringers.py # Copyright (c) 2013, Elettra - Sincrotrone Trieste S.C.p.A. # All rights reserved. # # Redistribution and use in source and binary forms, with or without # modification, are permitted provided that the following conditions are met: # # * Redistributions of source code must retain the above copyright # notice, this list of conditions and the following disclaimer. # * Redistributions in binary form must reproduce the above copyright # notice, this list of conditions and the following disclaimer in the # documentation and/or other materials provided with the distribution. # * Neither the name of the copyright holders nor the names of any # contributors may be used to endorse or promote products derived # from this software without specific prior written permission. # # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" # AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE # IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE # ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE # LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR # CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF # SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS # INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN # CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) # ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE # POSSIBILITY OF SUCH DAMAGE. from numpy import real, float32, uint16, ComplexWarning, iinfo, finfo from numpy import exp, arange, floor, kron, ones from numpy import zeros, mean, median, var, vstack, copy from numpy.lib import pad from numpy.fft import fft, fftshift, ifft, ifftshift, fft2, ifft2 from scipy.signal import filtfilt, butter, medfilt from pywt import wavedec2, waverec2 from warnings import simplefilter """De-striping filters better suited for CT ring artifacts compensation. :Organization: Elettra - Sincrotrone Trieste S.C.p.A. :Version: 2013.05.01 Requirements ------------ * `Python 2.7 ` * `Numpy 1.7.1 ` * `Scipy 1.7 ` Notes ----- The API is not stable yet and might change in the future. References ---------- [1] F. Brun, A. Accardo, G. Kourousias, D. Dreossi, R. Pugliese. Effective implementation of ring artifacts removal filters for synchrotron radiation microtomographic images. Proc. of the 8th International Symposium on Image and Signal Processing (ISPA), pp. 672-676, Sept. 4-6, Trieste (Italy), 2013. """ def rivers_filter(im, n=11): """Process a sinogram image with the Rivers de-striping algorithm. Parameters ---------- im : array_like Image data as numpy array. n : int Size of the median low-pass filtering. Example (using tiffile.py) -------------------------- >>> im = imread('original.tif') >>> im = rivers_filter(im, 13) >>> imsave('filtered.tif', im) References ---------- M. Rivers, http://cars9.uchicago.edu/software/epics/tomoRecon.html. """ # Copy image: im_f = copy(im.astype(float32)) # Compute mean of each column: col = im_f.mean(axis=0) # Perform low pass filtering: flt_col = medfilt(col, n) # Apply compensation on each row: for i in range(0, im_f.shape[0]): im_f[i,] = im_f[i,] - (col - flt_col) # Return image according to input type: if (im.dtype == 'uint16'): # Check extrema for uint16 images: im_f [ im_f < iinfo(uint16).min ] = iinfo(uint16).min im_f [ im_f > iinfo(uint16).max ] = iinfo(uint16).max # Return image: return im_f.astype(uint16) else: return im_f def boinhaibel_filter(im, n=11): """Process a sinogram image with the Boin and Haibel de-striping algorithm. Parameters ---------- im : array_like Image data as numpy array. n : int Size of the median low-pass filtering. Example (using tiffile.py) -------------------------- >>> im = imread('original.tif') >>> im = boinhaibel_filter(im, 11) >>> imsave('filtered.tif', im) References ---------- M. Boin and A. Haibel, Compensation of ring artefacts in synchrotron tomographic images, Optics Express 14(25):12071-12075, 2006. """ # Copy image: im_f = copy(im.astype(float32)) # Compute sum of each column (avoid further division by zero): col = im_f.sum(axis=0) + finfo(float32).eps # Perform low pass filtering: flt_col = medfilt(col, n) # Apply compensation on each row: for i in range(0, im_f.shape[0]): im_f[i,] = im_f[i,] * (flt_col/col) # Return image according to input type: if (im.dtype == 'uint16'): # Check extrema for uint16 images: im_f [ im_f < iinfo(uint16).min ] = iinfo(uint16).min im_f [ im_f > iinfo(uint16).max ] = iinfo(uint16).max # Return image: return im_f.astype(uint16) else: return im_f def raven_filter(im, n, u0): """Process a sinogram image with the Raven de-striping algorithm. Parameters ---------- im : array_like Image data as numpy array. n : int Order of the Butterworth low-pass filtering. u0 : float Cutoff frequency of the Butterworth low-pass filtering. Example (using tiffile.py) -------------------------- >>> im = imread('original.tif') >>> im = raven_filter(im, 0.01) >>> imsave('filtered.tif', im) References ---------- C. Raven, Numerical removal of ring artifacts in microtomography, Review of Scientific Instruments 69(8):2978-2980, 1998. """ # Replicate padding of input image: padsize = int(im.shape[1]/2) im_f = pad(im.astype(float32), ( (0,0), (padsize,padsize) ), 'edge') # Compute FT: im_f = fft2(im_f) # Apply 1-D low-pass Butterworth filtering at zero frequency: b, a = butter(n, u0) im_f[0,:] = filtfilt(b, a, im_f[0,:]) # Compute inverse FFT of the filtered data: im_f = real(ifft2( im_f )) # Return image according to input type: if (im.dtype == 'uint16'): # Check extrema for uint16 images: im_f [ im_f < iinfo(uint16).min ] = iinfo(uint16).min im_f [ im_f > iinfo(uint16).max ] = iinfo(uint16).max # Return image: return im_f.astype(uint16) else: return im_f def munchetal_filter(im, wlevel, sigma, wname='db15'): """Process a sinogram image with the Munch et al. de-striping algorithm. Parameters ---------- im : array_like Image data as numpy array. wname : {'haar', 'db1'-'db20', 'sym2'-'sym20', 'coif1'-'coif5', 'dmey'} The wavelet transform to use. wlevel : int Levels of the wavelet decomposition. sigma : float Cutoff frequency of the Butterworth low-pass filtering. Example (using tiffile.py) -------------------------- >>> im = imread('original.tif') >>> im = munchetal_filter(im, 'db15', 4, 1.0) >>> imsave('filtered.tif', im) 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. """ # Disable a warning: simplefilter("ignore", ComplexWarning) # Wavelet decomposition: coeffs = wavedec2(im.astype(float32), wname, level=wlevel) coeffsFlt = [coeffs[0]] # FFT transform of horizontal frequency bands: for i in range(1, wlevel + 1): # FFT: fcV = fftshift(fft(coeffs[i][1], axis=0)) my, mx = fcV.shape # Damping of vertical stripes: damp = 1-exp(-(arange(-floor(my/2.),-floor(my/2.)+my)**2)/(2*(sigma**2))) dampprime = kron(ones((1,mx)), damp.reshape((damp.shape[0],1))) fcV = fcV * dampprime # Inverse FFT: fcVflt = ifft(ifftshift(fcV), axis=0) cVHDtup = (coeffs[i][0], fcVflt, coeffs[i][2]) coeffsFlt.append(cVHDtup) # Get wavelet reconstruction: im_f = real(waverec2(coeffsFlt, wname)) # Return image according to input type: if (im.dtype == 'uint16'): # Check extrema for uint16 images: im_f [ im_f < iinfo(uint16).min ] = iinfo(uint16).min im_f [ im_f > iinfo(uint16).max ] = iinfo(uint16).max # Return filtered image (an additional row and/or column might be present): return im_f[0:im.shape[0],0:im.shape[1]].astype(uint16) else: return im_f[0:im.shape[0],0:im.shape[1]] def sijberspostnov_filter(im, winsize, thresh): """Process a sinogram image with the Sijbers and Postnov de-striping algorithm. Parameters ---------- im : array_like Image data as numpy array. winsize : int Size of the local floating window used to look for homogeneity. thresh : float Image rows (within the floating window) having variance below this tresh will be corrected. Example (using tiffile.py) -------------------------- >>> im = imread('original.tif') >>> im = sijberspostnov_filter(im, 51, 0.001) >>> imsave('filtered.tif', im) References ---------- J. Sijbers and A. Postnov, Reduction of ring artifacts in high resolution micro-CT reconstructions, Physics in Medicine and Biology 49(14):247-253, 2004. """ # Initializations: im_f = copy(im).astype(float32) dimx = im.shape[1] dimy = im.shape[0] # Normalize thresh parameter: #thresh = thresh*65536.0 glob_art = zeros(dimx) prevsize = 0 # Within a sliding window: for i in range(0, dimx-winsize): ct = 0 matrix = zeros(winsize) # For each line of the current window: for j in range(0, dimy): # Compute the variance within current sliding window: v = im_f[j, i:(i+winsize) ] curr_var = var(v) # If variance is below threshold: if (curr_var < thresh): # Add current line with mean subtracted to a temporary matrix: v = v - mean(v) matrix = vstack([matrix,v]) ct = ct + 1 # Determine local artifact correction vector: if (ct > 1): ct = ct - 1 matrix = matrix[1:ct,:] loc_art = median(matrix, axis=0) else: if (ct == 1): loc_art = matrix[1,:] else: loc_art = zeros(winsize) # Determine global artifact correction vector: for k in range(0, winsize): if ( matrix.shape[0] > prevsize ): glob_art[k+i] = loc_art[k] prevsize = matrix.shape[0] # Correct each line of the input image: for i in range(0, im.shape[0]): im_f[i,:] = im_f[i,:] - glob_art # Return image according to input type: if (im.dtype == 'uint16'): # Check extrema for uint16 images: im_f [ im_f < iinfo(uint16).min ] = iinfo(uint16).min im_f [ im_f > iinfo(uint16).max ] = iinfo(uint16).max # Return image: return im_f.astype(uint16) else: return im_f