tomopy.misc.corr
¶
Module for data correction and masking functions.
Functions:
|
Change dynamic range of values in an array. |
|
Apply circular mask to a 3D array. |
|
Apply Gaussian filter to 3D array along specified axis. |
|
Apply median filter to 3D array along specified axis. |
|
Apply median filter to 3D array along 0 axis with GPU support. |
|
Remove nonfinite values from a 3D array using an in-place 2D median filter. |
|
Apply Sobel filter to 3D array along specified axis. |
|
Replace NaN values in array with a given value. |
|
Replace negative values in array with a given value. |
|
Remove high intensity bright spots from a N-dimensional array by chunking along the specified dimension, and performing (N-1)-dimensional median filtering along the other dimensions. |
|
Remove high intensity bright spots from a 3D array along axis 0 dimension using GPU. |
|
Remove ring artifacts from images in the reconstructed domain. |
- tomopy.misc.corr.adjust_range(arr, dmin=None, dmax=None)[source]¶
Change dynamic range of values in an array.
- Parameters
arr (ndarray) – Input array.
dmin, dmax (float, optional) – Mininum and maximum values to rescale data.
- Returns
ndarray – Output array.
- tomopy.misc.corr.circ_mask(arr, axis, ratio=1, val=0.0, ncore=None)[source]¶
Apply circular mask to a 3D array.
- Parameters
arr (ndarray) – Arbitrary 3D array.
axis (int) – Axis along which mask will be performed.
ratio (int, optional) – Ratio of the mask’s diameter in pixels to the smallest edge size along given axis.
val (int, optional) – Value for the masked region.
- Returns
ndarray – Masked array.
- tomopy.misc.corr.enhance_projs_aps_1id(imgstack, median_ks=5, ncore=None)[source]¶
Enhance the projection images with weak contrast collected at APS 1ID
This filter uses a median fileter (will be switched to enhanced recursive median fileter, ERMF, in the future) for denoising, and a histogram equalization for dynamic range adjustment to bring out the details.
- Parameters
imgstack (np.ndarray) – tomopy images stacks (axis_0 is the oemga direction)
median_ks (int, optional) – 2D median filter kernel size for local noise suppresion
ncore (int, optional) – number of cores used for speed up
- Returns
ndarray – 3D enhanced image stacks.
- tomopy.misc.corr.gaussian_filter(arr, sigma=3, order=0, axis=0, ncore=None)[source]¶
Apply Gaussian filter to 3D array along specified axis.
- Parameters
arr (ndarray) – Input array.
sigma (scalar or sequence of scalars) – Standard deviation for Gaussian kernel. The standard deviations of the Gaussian filter are given for each axis as a sequence, or as a single number, in which case it is equal for all axes.
order ({0, 1, 2, 3} or sequence from same set, optional) – Order of the filter along each axis is given as a sequence of integers, or as a single number. An order of 0 corresponds to convolution with a Gaussian kernel. An order of 1, 2, or 3 corresponds to convolution with the first, second or third derivatives of a Gaussian. Higher order derivatives are not implemented
axis (int, optional) – Axis along which median filtering is performed.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
ndarray – 3D array of same shape as input.
- tomopy.misc.corr.median_filter(arr, size=3, axis=0, ncore=None)[source]¶
Apply median filter to 3D array along specified axis.
- Parameters
arr (ndarray) – Input array.
size (int, optional) – The size of the filter.
axis (int, optional) – Axis along which median filtering is performed.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
ndarray – Median filtered 3D array.
- tomopy.misc.corr.median_filter_cuda(arr, size=3, axis=0)[source]¶
Apply median filter to 3D array along 0 axis with GPU support. The winAllow is for A6000, Tian X support 3 to 8
- Parameters
arr (ndarray) – Input array.
size (int, optional) – The size of the filter.
axis (int, optional) – Axis along which median filtering is performed.
- Returns
ndarray – Median filtered 3D array.
Example
import tomocuda tomocuda.remove_outlier_cuda(arr, dif, 5)
For more information regarding install and using tomocuda, check https://github.com/kyuepublic/tomocuda for more information
- tomopy.misc.corr.median_filter_nonfinite(arr, size=3, callback=None)[source]¶
Remove nonfinite values from a 3D array using an in-place 2D median filter.
The 2D selective median filter is applied along the last two axes of the array.
New in version 1.11.
- Parameters
arr (ndarray) – The 3D array with nonfinite values in it.
size (int, optional) – The size of the filter.
callback (func(total, description, unit)) – A function called after every internal loop iteration. total is number of loop iterations. description is ‘Nonfinite median filter’. unit is ‘ prjs’.
- Returns
ndarray – The corrected 3D array with all nonfinite values removed based upon the local median value defined by the kernel size.
- Raises
ValueError – If the filter comes across a kernel only containing non-finite values a ValueError is raised for the user to increase their kernel size.
- tomopy.misc.corr.remove_nan(arr, val=0.0, ncore=None)[source]¶
Replace NaN values in array with a given value.
- Parameters
arr (ndarray) – Input array.
val (float, optional) – Values to be replaced with NaN values in array.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
ndarray – Corrected array.
- tomopy.misc.corr.remove_neg(arr, val=0.0, ncore=None)[source]¶
Replace negative values in array with a given value.
- Parameters
arr (ndarray) – Input array.
val (float, optional) – Values to be replaced with negative values in array.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
ndarray – Corrected array.
- tomopy.misc.corr.remove_outlier(arr, dif, size=3, axis=0, ncore=None, out=None)[source]¶
Remove high intensity bright spots from a N-dimensional array by chunking along the specified dimension, and performing (N-1)-dimensional median filtering along the other dimensions.
- Parameters
arr (ndarray) – Input array.
dif (float) – Expected difference value between outlier value and the median value of the array.
size (int) – Size of the median filter.
axis (int, optional) – Axis along which to chunk.
ncore (int, optional) – Number of cores that will be assigned to jobs.
out (ndarray, optional) – Output array for result. If same as arr, process will be done in-place.
- Returns
ndarray – Corrected array.
- tomopy.misc.corr.remove_outlier1d(arr, dif, size=3, axis=0, ncore=None, out=None)[source]¶
Remove high intensity bright spots from an array, using a one-dimensional median filter along the specified axis.
- Parameters
arr (ndarray) – Input array.
dif (float) – Expected difference value between outlier value and the median value of the array.
size (int) – Size of the median filter.
axis (int, optional) – Axis along which median filtering is performed.
ncore (int, optional) – Number of cores that will be assigned to jobs.
out (ndarray, optional) – Output array for result. If same as arr, process will be done in-place.
- Returns
ndarray – Corrected array.
- tomopy.misc.corr.remove_outlier_cuda(arr, dif, size=3, axis=0)[source]¶
Remove high intensity bright spots from a 3D array along axis 0 dimension using GPU.
- Parameters
arr (ndarray) – Input array.
dif (float) – Expected difference value between outlier value and the median value of the array.
size (int) – Size of the median filter.
axis (int, optional) – Axis along which outlier removal is performed.
- Returns
ndarray – Corrected array.
Example
>>> import tomocuda >>> tomocuda.remove_outlier_cuda(arr, dif, 5)
For more information regarding install and using tomocuda, check https://github.com/kyuepublic/tomocuda for more information
- tomopy.misc.corr.remove_ring(rec, center_x=None, center_y=None, thresh=300.0, thresh_max=300.0, thresh_min=- 100.0, theta_min=30, rwidth=30, int_mode='WRAP', ncore=None, nchunk=None, out=None)[source]¶
Remove ring artifacts from images in the reconstructed domain. Descriptions of parameters need to be more clear for sure.
- Parameters
arr (ndarray) – Array of reconstruction data
center_x (float, optional) – abscissa location of center of rotation
center_y (float, optional) – ordinate location of center of rotation
thresh (float, optional) – maximum value of an offset due to a ring artifact
thresh_max (float, optional) – max value for portion of image to filter
thresh_min (float, optional) – min value for portion of image to filer
theta_min (int, optional) – Features larger than twice this angle (degrees) will be considered a ring artifact. Must be less than 180 degrees.
rwidth (int, optional) – Maximum width of the rings to be filtered in pixels
int_mode (str, optional) – ‘WRAP’ for wrapping at 0 and 360 degrees, ‘REFLECT’ for reflective boundaries at 0 and 180 degrees.
ncore (int, optional) – Number of cores that will be assigned to jobs.
nchunk (int, optional) – Chunk size for each core.
out (ndarray, optional) – Output array for result. If same as arr, process will be done in-place.
- Returns
ndarray – Corrected reconstruction data
- tomopy.misc.corr.sobel_filter(arr, axis=0, ncore=None)[source]¶
Apply Sobel filter to 3D array along specified axis.
- Parameters
arr (ndarray) – Input array.
axis (int, optional) – Axis along which sobel filtering is performed.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
ndarray – 3D array of same shape as input.