def test_sort_backward(self): mat = np.transpose( np.tile(np.arange(self.size - 1, -1, -1), (self.size, 1))) mat1, mat_idx = util.sort_forward(mat, axis=0) mat2 = util.sort_backward(mat1, mat_idx, axis=0) num = np.mean(np.abs(mat2 - mat)) self.assertTrue(num < self.eps)
def add_stripe_artifact(sinogram, size, position, strength_ratio=0.2, stripe_type="partial"): """ Add stripe artifacts to a sinogram. Parameters ---------- sinogram: array_like 2D array. Sinogram image. size : int Size of stripe artifact. position : int Position of the stripe. strength_ratio : float To define the strength of the artifact. The value is in the range of [0.0, 1.0]. stripe_type : {"partial", "full", "dead", "fluctuating"} Type of stripe as classified in Ref. [1]. Returns ------- array_like References ---------- .. [1] https://doi.org/10.1364/OE.26.028396 """ sinogram = np.copy(sinogram) (nrow, ncol) = sinogram.shape position = np.clip(position, 0, ncol - size - 1) strength_ratio = np.clip(strength_ratio, 0.0, 1.0) stripe = sinogram[:, position:position + size] if stripe_type == "partial": stripe_sort, mat_idx = util.sort_forward(stripe, axis=0) pos = int((1.0 - strength_ratio) * nrow) list_ratio = np.ones(nrow, dtype=np.float32) list_ratio[pos:nrow] = np.linspace(1.0, 1.0 - strength_ratio, nrow - pos) mat_ratio = np.tile(list_ratio, (size, 1)) stripe_sort = stripe_sort * np.transpose(mat_ratio) stripe_mod = util.sort_backward(stripe_sort, mat_idx, axis=0) elif stripe_type == "dead": stripe_mod = np.ones_like(stripe) * strength_ratio * np.max(sinogram) elif stripe_type == "fluctuating": std_dev = np.mean(sinogram[sinogram != 0.0]) * strength_ratio noise = np.random.normal(0.0, std_dev, size=stripe.shape) stripe_mod = stripe + noise else: list_ratio = (1 - strength_ratio) * np.ones(nrow, dtype=np.float32) mat_ratio = np.tile(list_ratio, (size, 1)) stripe_mod = stripe * np.transpose(mat_ratio) sinogram[:, position:position + size] = stripe_mod return sinogram
def remove_stripe_based_sorting(sinogram, size=21, dim=1, **options): """ Remove stripe artifacts in a sinogram using the sorting technique, algorithm 3 in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image. size : int Window size of the median filter. dim : {1, 2}, optional Dimension of the window. options : dict, optional Use another smoothing filter rather than the median filter. E.g. options={"method": "gaussian_filter", "para1": (1,21))} Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" sino_sort, sino_index = util.sort_forward(np.float32(sinogram), axis=0) if len(options) == 0: if dim == 2: sino_sort = ndi.median_filter(sino_sort, (size, size)) else: sino_sort = ndi.median_filter(sino_sort, (1, size)) else: if not isinstance(options, dict): raise ValueError(msg) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] para = tuple(opt.values())[1:] if method in dir(ndi): try: sino_sort = getattr(ndi, method)(sino_sort, *para) except: raise ValueError(msg) else: if method in dir(util): try: sino_sort = getattr(util, method)(sino_sort, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) return util.sort_backward(sino_sort, sino_index, axis=0)
def remove_stripe_based_fft(sinogram, u=20, n=8, v=1, sort=False): """ Remove stripes using the method in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image. u : int Cutoff frequency. n : int Filter order. v : int Number of rows (* 2) to be applied the filter. sort : bool, optional Apply sorting (Ref. [2]) if True. Returns ------- ndarray 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1063/1.1149043 .. [2] https://doi.org/10.1364/OE.26.028396 """ if sort is True: sinogram, sino_index = util.sort_forward(sinogram, axis=0) pad = min(150, int(0.1 * np.min(sinogram.shape))) sinogram = np.pad(sinogram, ((pad, pad), (0, 0)), mode='mean') sinogram = np.pad(sinogram, ((0, 0), (pad, pad)), mode='edge') (nrow, ncol) = sinogram.shape window_2d = util.make_2d_butterworth_window(ncol, nrow, u, v, n) sinogram = fft.ifft2( np.fft.ifftshift(np.fft.fftshift(fft.fft2(sinogram)) * window_2d)) sinogram = np.abs(sinogram[pad:nrow - pad, pad:ncol - pad]) if sort is True: sinogram = util.sort_backward(sinogram, sino_index, axis=0) return sinogram
def remove_stripe_based_filtering(sinogram, sigma=3, size=21, dim=1, sort=True, **options): """ Remove stripe artifacts in a sinogram using the filtering technique, algorithm 2 in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image sigma : int Sigma of the Gaussian window used to separate the low-pass and high-pass components of the intensity profile of each column. size : int Window size of the median filter. dim : {1, 2}, optional Dimension of the window. sort : bool, optional Apply sorting if True. options : dict, optional Use another smoothing filter rather than the median filter. E.g. options={"method": "gaussian_filter", "para1": (1,21))}. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" window = {"name": "gaussian", "sigma": sigma} sino_smooth, sino_sharp = util.separate_frequency_component( np.float32(sinogram), axis=0, window=window) if sort is True: sino_smooth, sino_index = util.sort_forward(sino_smooth, axis=0) if len(options) == 0: if dim == 2: sino_smooth = ndi.median_filter(sino_smooth, (size, size)) else: sino_smooth = ndi.median_filter(sino_smooth, (1, size)) else: if not isinstance(options, dict): raise ValueError(msg) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] if method in dir(ndi): para = tuple(opt.values())[1:] try: sino_smooth = getattr(ndi, method)(sino_smooth, *para) except: raise ValueError(msg) else: if method in dir(util): try: sino_smooth = getattr(util, method)(sino_smooth, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) if sort is True: sino_smooth = util.sort_backward(sino_smooth, sino_index, axis=0) return sino_smooth + sino_sharp
def remove_stripe_based_wavelet_fft(sinogram, level=5, size=1, wavelet_name="db9", window_name="gaussian", sort=False, **options): """ Remove stripes using the method in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image. level : int Wavelet decomposition level. size : int Damping parameter. Larger is stronger. wavelet_name : str Name of a wavelet. Search pywavelets API for a full list. window_name : str High-pass window. Two options: "gaussian" or "butter". sort : bool, optional Apply sorting (Ref. [2]) if True. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1364/OE.17.008567 .. [2] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" if sort is True: sinogram, sino_index = util.sort_forward(sinogram, axis=0) (nrow, ncol) = sinogram.shape pad = min(150, int(0.1 * min(nrow, ncol))) sinogram = np.pad(sinogram, ((pad, pad), (0, 0)), mode='mean') sinogram = np.pad(sinogram, ((0, 0), (pad, pad)), mode='edge') output_data = util.apply_wavelet_decomposition(sinogram, wavelet_name, level=level) output_data = [list(data) for data in output_data] n_level = len(output_data[1:]) for i in range(1, n_level + 1): if len(options) == 0: (height, width) = output_data[i][1].shape window = np.transpose( util.make_2d_damping_window(height, width, size, window_name)) output_data[i][1] = np.real( np.fft.ifft2( np.fft.ifftshift( np.fft.fftshift(np.fft.fft2(output_data[i][1])) * window))) else: if not isinstance(options, dict): raise ValueError(msg) mat_smooth = np.copy(output_data[i][1]) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] para = tuple(opt.values())[1:] if method in dir(ndi): try: mat_smooth = getattr(ndi, method)(mat_smooth, *para) except: raise ValueError(msg) elif method in dir(util): try: mat_smooth = getattr(util, method)(mat_smooth, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) output_data[i][1] = mat_smooth sinogram = util.apply_wavelet_reconstruction(output_data, wavelet_name) sinogram = sinogram[pad:nrow + pad, pad:ncol + pad] if sort is True: sinogram = util.sort_backward(sinogram, sino_index, axis=0) return sinogram
def remove_stripe_based_regularization(sinogram, alpha=0.0005, num_chunk=1, apply_log=True, sort=True): """ Remove stripes using the method in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image. alpha : float Regularization parameter, e.g. 0.0005. Smaller is stronger. num_chunk : int Number of chunks of rows. apply_log : bool Apply the logarithm function to the sinogram if True. sort : bool, optional Apply sorting (Ref. [2]) if True. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1016/j.aml.2010.08.022 .. [2] https://doi.org/10.1364/OE.26.028396 """ (nrow, ncol) = sinogram.shape sinogram = np.copy(sinogram) if sort is True: sinogram, sino_index = util.sort_forward(sinogram, axis=0) if apply_log is True: if np.any(sinogram <= 0.0): nmean = np.mean(np.abs(sinogram)) sinogram[sinogram <= 0.0] = nmean sinogram = -np.log(sinogram) else: sinogram = -np.log(sinogram) sijmat = util.calculate_regularization_coefficient(ncol, alpha) list_index = np.array_split(np.arange(nrow), num_chunk) list_grad = np.zeros(ncol, dtype=np.float32) mat_grad = np.zeros((ncol, ncol), dtype=np.float32) for pos in list_index: bindex = pos[0] eindex = pos[-1] + 1 list_mean = np.mean(sinogram[bindex:eindex], axis=0) list_grad[1:-1] = (-1) * np.diff(list_mean, 2) list_grad[0] = list_mean[0] - list_mean[1] list_grad[-1] = list_mean[-1] - list_mean[-2] mat_grad[:] = list_grad list_coe = np.sum(mat_grad * sijmat, axis=1) mat_coe = np.tile(list_coe, (eindex - bindex, 1)) sinogram[bindex:eindex, :] = sinogram[bindex:eindex, :] + mat_coe if sort is True: sinogram = util.sort_backward(sinogram, sino_index, axis=0) if apply_log is True: sinogram = np.exp(-sinogram) return sinogram
def remove_stripe_based_normalization(sinogram, sigma=15, num_chunk=1, sort=True, **options): """ Remove stripes using the method in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image. sigma : int Sigma of the Gaussian window. num_chunk : int Number of chunks of rows. sort : bool, optional Apply sorting (Ref. [2]) if True. options : dict, optional Use another smoothing 1D-filter rather than the Gaussian filter. E.g. options={"method": "median_filter", "para1": 21)}. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://www.mcs.anl.gov/research/projects/X-ray-cmt/rivers/ tutorial.html .. [2] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" \ "\n Note that the filter must be a 1D-filter." (nrow, _) = sinogram.shape sinogram = np.copy(sinogram) if sort is True: sinogram, sino_index = util.sort_forward(sinogram, axis=0) list_index = np.array_split(np.arange(nrow), num_chunk) for pos in list_index: bindex = pos[0] eindex = pos[-1] + 1 list_mean = np.mean(sinogram[bindex:eindex], axis=0) if len(options) == 0: list_filt = ndi.gaussian_filter(list_mean, sigma) else: if not isinstance(options, dict): raise ValueError(msg) list_filt = np.copy(list_mean) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] if method in dir(ndi): para = tuple(opt.values())[1:] try: list_filt = getattr(ndi, method)(list_filt, *para) except: raise ValueError(msg) else: if method in dir(util): try: list_filt = getattr(util, method)(list_filt, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) list_coe = list_filt - list_mean matcoe = np.tile(list_coe, (eindex - bindex, 1)) sinogram[bindex:eindex, :] = sinogram[bindex:eindex, :] + matcoe if sort is True: sinogram = util.sort_backward(sinogram, sino_index, axis=0) return sinogram
def remove_large_stripe(sinogram, snr=3.0, size=51, drop_ratio=0.1, norm=True, **options): """ Remove large stripe artifacts in a sinogram, algorithm 5 in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image snr : float Ratio (>1.0) used to detect stripe locations. Greater is less sensitive. size : int Window size of the median filter. drop_ratio : float, optional Ratio of pixels to be dropped, which is used to to reduce the possibility of the false detection of stripes. norm : bool, optional Apply normalization if True. options : dict, optional Use another smoothing filter rather than the median filter. E.g. options={"method": "gaussian_filter", "para1": (1,21))}. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" sinogram = np.copy(np.float32(sinogram)) drop_ratio = np.clip(drop_ratio, 0.0, 0.8) (nrow, ncol) = sinogram.shape ndrop = int(0.5 * drop_ratio * nrow) sino_sort, sino_index = util.sort_forward(sinogram, axis=0) if len(options) == 0: sino_smooth = ndi.median_filter(sino_sort, (1, size)) else: if not isinstance(options, dict): raise ValueError(msg) sino_smooth = np.copy(sino_sort) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] if method in dir(ndi): para = tuple(opt.values())[1:] try: sino_smooth = getattr(ndi, method)(sino_smooth, *para) except: raise ValueError(msg) else: if method in dir(util): try: sino_smooth = getattr(util, method)(sino_smooth, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) list1 = np.mean(sino_sort[ndrop:nrow - ndrop], axis=0) list2 = np.mean(sino_smooth[ndrop:nrow - ndrop], axis=0) list_fact = np.divide(list1, list2, out=np.ones_like(list1), where=list2 != 0) list_mask = util.detect_stripe(list_fact, snr) list_mask = np.float32(ndi.binary_dilation(list_mask, iterations=1)) if norm is True: sinogram = sinogram / np.tile(list_fact, (nrow, 1)) sino_corr = util.sort_backward(sino_smooth, sino_index, axis=0) xlist_miss = np.where(list_mask > 0.0)[0] sinogram[:, xlist_miss] = sino_corr[:, xlist_miss] return sinogram
def remove_stripe_based_fitting(sinogram, order=2, sigma=10, sort=False, num_chunk=1, **options): """ Remove stripe artifacts in a sinogram using the fitting technique, algorithm 1 in Ref. [1]. Angular direction is along the axis 0. Parameters ---------- sinogram : array_like 2D array. Sinogram image order : int Polynomial fit order. sigma : int Sigma of the Gaussian window in the x-direction. Smaller is stronger. sort : bool, optional Apply sorting if True. num_chunk : int Number of chunks of rows to apply the fitting. options : dict, optional Use another smoothing filter rather than the Fourier gaussian filter. E.g. options={"method": "gaussian_filter", "para1": (1,21))}. Returns ------- array_like 2D array. Stripe-removed sinogram. References ---------- .. [1] https://doi.org/10.1364/OE.26.028396 """ msg = "\n Please use the dictionary format: options={'method':" \ " 'filter_name', 'para1': parameter_1, 'para2': parameter_2}" (nrow, ncol) = sinogram.shape pad = min(150, int(0.1 * nrow)) sigmay = np.clip(min(60, int(0.1 * ncol)), 10, None) if sort is True: sinogram, sino_index = util.sort_forward(sinogram, axis=0) sino_fit = util.generate_fitted_image(sinogram, order, axis=0, num_chunk=num_chunk) if len(options) == 0: sino_filt = util.apply_gaussian_filter(sino_fit, sigma, sigmay, pad) else: if not isinstance(options, dict): raise ValueError(msg) sino_filt = np.copy(sino_fit) for opt_name in options: opt = options[opt_name] method = tuple(opt.values())[0] if method in dir(ndi): para = tuple(opt.values())[1:] try: sino_filt = getattr(ndi, method)(sino_filt, *para) except: raise ValueError(msg) else: if method in dir(util): try: sino_filt = getattr(util, method)(sino_filt, *para) except: raise ValueError(msg) else: raise ValueError("Can't find the method: '{}' in the" " namespace".format(method)) sino_filt = np.mean(np.abs(sino_fit)) * sino_filt / np.mean( np.abs(sino_filt)) sino_corr = ((sinogram / sino_fit) * sino_filt) if sort is True: sino_corr = util.sort_backward(sino_corr, sino_index, axis=0) return sino_corr