7.5.3. algotom.util.utility
¶
Module of utility methods:
Methods for parallel computing, geometric transformation, masking.
- Methods for customizing stripe/ring removal methods
sort_forward
sort_backward
separate_frequency_component
generate_fitted_image
detect_stripe
calculate_regularization_coefficient
make_2d_butterworth_window
make_2d_damping_window
apply_wavelet_decomposition
apply_wavelet_reconstruction
apply_filter_to_wavelet_component
interpolate_inside_stripe
transform_slice_forward
transform_slice_backward
- Customized smoothing filters:
apply_gaussian_filter (in the Fourier space)
apply_regularization_filter
- Methods for grid scans:
detect_sample
fix_non_sample_areas
locate_slice
locate_slice_chunk
- Methods for speckle-based tomography
generate_spiral_positions
- Method for finding the center of rotation by visual inspection.
find_center_visual_sinograms
Functions:
|
Apply a processing method to slices of a 3D array in parallel. |
|
Apply a geometric transformation to a 2D array |
|
Create a circle mask. |
|
Sort gray-scales of an image along an axis. |
|
Sort gray-scales of an image using an index array provided. |
|
Separate low and high frequency components of an image along an axis. |
|
Apply a polynomial fitting along an axis of an image. |
|
Locate stripe positions using Algorithm 4 in Ref. |
|
Calculate coefficients used for the regularization-based method. |
|
Create a 2d window from the 1D Butterworth window. |
|
Make 2D damping window from a list of 1D window for a Fourier-space filter, i.e. a high-pass filter. |
|
Apply 2D wavelet decomposition. |
|
Apply 2D wavelet reconstruction. |
|
Apply a filter to a component of the wavelet decomposition of an image. |
|
Interpolate gray-scales inside vertical stripes of an image. |
|
Generate coordinates of a rectangular grid from polar coordinates. |
|
Generate polar coordinates from grid coordinates. |
|
Transform a reconstructed image into polar coordinates. |
|
Transform a reconstructed image in polar coordinates back to rectangular coordinates. |
|
Create a 2D Gaussian window. |
|
Filtering an image in the Fourier domain using a 2D Gaussian window. |
|
Apply a regularization filter using the method in Ref. |
|
Transform a 1d-window to 2d-window. |
|
To check if there is a sample in a sinogram using the "double-wedge" property of the Fourier transform of the sinogram (Ref. |
|
Used to fix overlap values of grid-cells without sample by copying from its neighbours. |
|
Locate slice indices in grid-rows given a slice index of the reconstruction data as a whole. |
|
Locate slice indices in grid-rows given slice indices of the reconstruction data as a whole. |
|
Generate Fermat spiral positions. |
|
For visually finding the center-of-rotation (COR) using converted 360-degree sinograms from a 180-degree sinogram at different CORs (Ref. |
- algotom.util.utility.parallel_process_slices(data, method, parameters, axis=1, ncore=None, prefer='threads', **kwargs)[source]¶
Apply a processing method to slices of a 3D array in parallel.
- Parameters
data (array_like or hdf object) – 3D data array or HDF dataset.
method (callable) – Function to apply to each slice.
parameters (list) – List of positional parameters for the method.
axis (int) – Axis along which the method is applied.
ncore (int, optional) – Number of cpu-cores used for computing. Automatically selected if None.
prefer ({“threads”, “processes”}) – Preferred parallel backend (“threads” for I/O bound tasks or “processes” for CPU-bound tasks).
**kwargs (dict) – Additional keyword parameters to pass to the method.
- Returns
array_like – Same axis-definition as the input.
- algotom.util.utility.mapping(mat, x_mat, y_mat, order=1, mode='reflect')[source]¶
Apply a geometric transformation to a 2D array
- Parameters
mat (array_like) – 2D array.
x_mat (array_like) – 2D array of the x-coordinates.
y_mat (array_like) – 2D array of the y-coordinates.
order (int, optional) – The order of the spline interpolation, default is 1. The order has to be in the range 0-5.
mode ({‘reflect’, ‘constant’, ‘nearest’, ‘mirror’, ‘wrap’}, optional) – The mode parameter determines how the input array is extended beyond its boundaries. Default is ‘reflect’.
- Returns
array_like – 2D array.
- algotom.util.utility.make_circle_mask(width, ratio)[source]¶
Create a circle mask.
- Parameters
width (int) – Width of a square array.
ratio (float) – Ratio between the diameter of the mask and the width of the array.
- Returns
array_like – Square array.
- algotom.util.utility.sort_forward(mat, axis=0)[source]¶
Sort gray-scales of an image along an axis. e.g. axis=0 is to sort along each column.
- Parameters
mat (array_like) – 2D array.
axis (int) – Axis along which to sort.
- Returns
mat_sort (array_like) – 2D array. Sorted image.
mat_index (array_like) – 2D array. Index array used for sorting backward.
- algotom.util.utility.sort_backward(mat, mat_index, axis=0)[source]¶
Sort gray-scales of an image using an index array provided. e.g. axis=0 is to sort each column.
- Parameters
mat (array_like) – 2D array.
mat_index (array_like) – 2D array. Index array used for sorting.
axis (int) – Axis along which to sort.
- Returns
mat_sort (array_like) – 2D array. Sorted image.
- algotom.util.utility.separate_frequency_component(mat, axis=0, window=None)[source]¶
Separate low and high frequency components of an image along an axis. e.g. axis=0 is to apply the separation to each column.
- Parameters
mat (array_like) – 2D array.
axis (int) – Axis along which to apply the filter.
window (array_like or dict) – 1D array or a dictionary which given the name of a window in the scipy_window list and its parameters (without window-length). E.g window={“name”: “gaussian”, “sigma”: 5}
- Returns
mat_low (array_like) – 2D array. Low-frequency image.
mat_high (array_like) – 2D array. High-frequency image.
- algotom.util.utility.generate_fitted_image(mat, order, axis=0, num_chunk=1)[source]¶
Apply a polynomial fitting along an axis of an image. e.g. axis=0 is to apply the fitting to each column.
- Parameters
mat (array_like) – 2D array.
order (int) – Order of the polynomial used to fit.
axis (int) – Axis along which to apply the filter.
num_chunk (int) – Number of chunks of rows or columns to apply the fitting.
- Returns
mat_fit (array_like)
- algotom.util.utility.detect_stripe(list_data, snr)[source]¶
Locate stripe positions using Algorithm 4 in Ref. [1]
- Parameters
list_data (array_like) – 1D array. Normalized data.
snr (float) – Ratio (>1.0) for stripe detection. Greater is less sensitive.
- Returns
array_like – 1D binary mask.
References
- algotom.util.utility.calculate_regularization_coefficient(width, alpha)[source]¶
Calculate coefficients used for the regularization-based method. Eq. (7) in Ref. [1].
- Parameters
width (int) – Width of a square array.
alpha (float) – Regularization parameter.
- Returns
float – Square array.
References
- algotom.util.utility.make_2d_butterworth_window(width, height, u, v, n)[source]¶
Create a 2d window from the 1D Butterworth window.
- Parameters
height (int) – Height of the window.
width (int) – Width of the window.
u (int) – Cutoff frequency.
n (int) – Filter order.
v (int) – Number of rows (= 2*v) around the height middle are the 1D Butterworth windows.
- Returns
array_like – 2D array.
- algotom.util.utility.make_2d_damping_window(width, height, size, window_name='gaussian')[source]¶
Make 2D damping window from a list of 1D window for a Fourier-space filter, i.e. a high-pass filter.
- Parameters
height (int) – Height of the window.
width (int) – Width of the window.
size (int) – Sigma of a Gaussian window or cutoff frequency of a Butterworth window.
window_name (str, optional) – Two options: “gaussian” or “butter”.
- Returns
array_like – 2D array of the window.
- algotom.util.utility.apply_wavelet_decomposition(mat, wavelet_name, level=None)[source]¶
Apply 2D wavelet decomposition.
- Parameters
mat (array_like) – 2D array.
wavelet_name (str) – Name of a wavelet. E.g. “db5”
level (int, optional) – Decomposition level. It is constrained to return an array with a minimum size of larger than 16 pixels.
- Returns
list – The first element is an 2D-array, next elements are tuples of three 2D-arrays. i.e [mat_n, (cH_level_n, cV_level_n, cD_level_n), …, (cH_level_1, cV_level_1, cD_level_1)]
- algotom.util.utility.apply_wavelet_reconstruction(data, wavelet_name, ignore_level=None)[source]¶
Apply 2D wavelet reconstruction.
- Parameters
data (list or tuple) – The first element is an 2D-array, next elements are tuples of three 2D-arrays. i.e [mat_n, (cH_level_n, cV_level_n, cD_level_n), …, (cH_level_1, cV_level_1, cD_level_1)].
wavelet_name (str) – Name of a wavelet. E.g. “db5”
ignore_level (int, optional) – Decomposition level to be ignored for reconstruction.
- Returns
array_like – 2D array. Note that the sizes of the array are always even numbers.
- algotom.util.utility.check_level(level, n_level)[source]¶
Supplementary method for the method of “apply_filter_to_wavelet_component”. To check if the provided level is in the right format.
- algotom.util.utility.apply_filter_to_wavelet_component(data, level=None, order=1, method=None, para=None)[source]¶
Apply a filter to a component of the wavelet decomposition of an image.
- Parameters
data (list or tuple) – The first element is an 2D-array, next elements are tuples of three 2D-arrays. i.e [mat_n, (cH_level_n, cV_level_n, cD_level_n), …, (cH_level_1, cV_level_1, cD_level_1)].
level (int, list of int, or None) – Decomposition level to be applied the filter.
order ({0, 1, 2}) – Specify which component in a tuple, (cH_level_n, cV_level_n, cD_level_n), to be filtered.
method (str) – Name of the filter in the namespace. E.g. method=”gaussian_filter”
para (list or tuple) – Parameters of the filter. E.g para=[(1,11)]
- Returns
list or tuple – The first element is an 2D-array, next elements are tuples of three 2D-arrays. i.e [mat_n, (cH_level_n, cV_level_n, cD_level_n), …, (cH_level_1, cV_level_1, cD_level_1)].
- algotom.util.utility.interpolate_inside_stripe(mat, list_mask, kind='linear')[source]¶
Interpolate gray-scales inside vertical stripes of an image. Stripe locations given by a binary 1D-mask.
- Parameters
mat (array_like) – 2D array.
list_mask (array_like) – 1D array. Must equal the width of an image.
kind ({‘linear’, ‘cubic’, ‘quintic’}, optional) – The kind of spline interpolation to use. Default is ‘linear’.
- Returns
array_like
- algotom.util.utility.rectangular_from_polar(width_reg, height_reg, width_pol, height_pol)[source]¶
Generate coordinates of a rectangular grid from polar coordinates.
- Parameters
width_reg (int) – Width of an image in the Cartesian coordinate system.
height_reg (int) – Height of an image in the Cartesian coordinate system.
width_pol (int) – Width of an image in the polar coordinate system.
height_pol (int) – Height of an image in the polar coordinate system.
- Returns
x_mat (array_like) – 2D array. Broadcast of the x-coordinates.
y_mat (array_like) – 2D array. Broadcast of the y-coordinates.
- algotom.util.utility.polar_from_rectangular(width_pol, height_pol, width_reg, height_reg)[source]¶
Generate polar coordinates from grid coordinates.
- Parameters
width_pol (int) – Width of an image in the polar coordinate system.
height_pol (int) – Height of an image in the polar coordinate system.
width_reg (int) – Width of an image in the Cartesian coordinate system.
height_reg (int) – Height of an image in the Cartesian coordinate system.
- Returns
r_mat (array_like) – 2D array. Broadcast of the r-coordinates.
theta_mat (array_like) – 2D array. Broadcast of the theta-coordinates.
- algotom.util.utility.transform_slice_forward(mat, coord_mat=None)[source]¶
Transform a reconstructed image into polar coordinates.
- Parameters
mat (array_like) – Square array. Reconstructed image.
coord_mat (tuple of array_like, optional) – (Square array of x-coordinates , square array of y-coordinates) or generated if None.
- Returns
array_like – Transformed image.
- algotom.util.utility.transform_slice_backward(mat, coord_mat=None)[source]¶
Transform a reconstructed image in polar coordinates back to rectangular coordinates.
- Parameters
mat (array_like) – Square array. Reconstructed image in polar coordinates.
coord_mat (tuple of array_like, optional) – (Square array of r-coordinates , square array of theta-coordinates) or generated if None.
- Returns
array_like – Transformed image.
- algotom.util.utility.make_2d_gaussian_window(height, width, sigma_x, sigma_y)[source]¶
Create a 2D Gaussian window.
- Parameters
height (int) – Height of the image.
width (int) – Width of the image.
sigma_x (int) – Sigma in the x-direction.
sigma_y (int) – Sigma in the y-direction.
- Returns
array_like – 2D array.
- algotom.util.utility.apply_gaussian_filter(mat, sigma_x, sigma_y, pad=None, mode=None)[source]¶
Filtering an image in the Fourier domain using a 2D Gaussian window. Smaller is stronger.
- Parameters
mat (array_like) – 2D array.
sigma_x (int) – Sigma in the x-direction.
sigma_y (int) – Sigma in the y-direction.
pad (int or None) – Padding for the Fourier transform.
mode (str, list of str, or tuple of str) – Padding method. One of options : ‘reflect’, ‘edge’, ‘constant’. Full list is at: https://numpy.org/doc/stable/reference/generated/numpy.pad.html
- Returns
array_like – 2D array. Filtered image.
- algotom.util.utility.apply_1d_regularizer(list_data, sijmat)[source]¶
Supplementary method for the method of “apply_regularization_filter”. To apply a regularizer to an 1D-array.
- algotom.util.utility.apply_regularization_filter(mat, alpha, axis=1, ncore=None)[source]¶
Apply a regularization filter using the method in Ref. [1]. Note that it’s computationally costly.
- Parameters
mat (array_like) – 2D array
alpha (float) – Regularization parameter, e.g. 0.001. Smaller is stronger.
axis (int) – Axis along which to apply the filter.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
array_like – 2D array. Smoothed image.
References
- algotom.util.utility.transform_1d_window_to_2d(win_1d, order=1, mode='reflect')[source]¶
Transform a 1d-window to 2d-window. Useful for designing a Fourier filter.
- Parameters
win_1d (array_like) – 1D array.
order (int, optional) – The order of the spline interpolation, default is 1. The order has to be in the range 0-5.
mode ({‘reflect’, ‘constant’, ‘nearest’, ‘mirror’, ‘wrap’}, optional) – The mode parameter determines how the input array is extended beyond its boundaries. Default is ‘reflect’.
- Returns
win_2d (array_like) – Square array, a 2D version of the 1d-window.
- algotom.util.utility.detect_sample(sinogram, sino_type='180')[source]¶
To check if there is a sample in a sinogram using the “double-wedge” property of the Fourier transform of the sinogram (Ref. [1]).
- Parameters
sinogram (array_like) – 2D array. Sinogram image
sino_type ({“180”, “360”}) – Sinogram type : 180-degree or 360-degree.
- Returns
bool – True if there is a sample.
References
- algotom.util.utility.fix_non_sample_areas(overlap_metadata, direction='horizontal')[source]¶
Used to fix overlap values of grid-cells without sample by copying from its neighbours. Input is a 3d-array of overlapping values for each grid cell. For example, to a 5 x 3 (n_row x n_column) grid scans, the shape for overlapping values in the horizontal direction is: 5 x 2 (n_column - 1) x 2 (overlap, side). The shape for overlapping values in the vertical direction is: 4 (n_row - 1) x 3 x 2 (overlap, side). The order of calculating overlapping values in a grid is left-to-right, top-to-bottom.
- Parameters
overlap_metadata (array_like) – A matrix of overlap values of each grid-cell where each element is a list of [overlap, side].
direction ({“horizontal”, “vertical”}) – Direction of overlapping calculation.
- Returns
metadata (array_like)
- algotom.util.utility.locate_slice(slice_idx, height, overlap_metadata)[source]¶
Locate slice indices in grid-rows given a slice index of the reconstruction data as a whole.
- Parameters
slice_idx (int) – Slice index of full reconstruction data.
height (int) – Height of a projection image of each grid-cell.
overlap_metadata (array_like) – A matrix of overlap values of each grid-row where each element is a list of [overlap, side]. Used to stitch the grid-data along the row-direction.
- Returns
list of int and float – If the slice is not in the overlapping area between two grid-rows, the result is a list of [grid_row_index, slice_index, weight_factor]. If the slice is in the overlapping area between two grid-rows, the result is a list of [[grid_row_index_0, slice_index_0, weight_factor_0], [grid_row_index_1, slice_index_1, weight_factor_1]]
- algotom.util.utility.locate_slice_chunk(slice_start, slice_stop, height, overlap_metadata)[source]¶
Locate slice indices in grid-rows given slice indices of the reconstruction data as a whole.
- Parameters
slice_start (int) – Starting index of full reconstruction data.
slice_stop (int) – Stopping index of full reconstruction data.
height (int) – Height of a projection image of each grid-cell.
overlap_metadata (array_like) – A matrix of overlap values of each grid-row where each element is a list of [overlap, side]. Used to stitch the grid-data along the row-direction.
- Returns
list of list of int and float – List of results for each slice index. If a slice is not in the overlapping area between two grid-rows, the result is a list of [grid_row_index, slice_index, weight_factor]. If a slice is in the overlapping area between two grid-rows, the result is a list of [[grid_row_index_0, slice_index_0, weight_factor_0], [grid_row_index_1, slice_index_1, weight_factor_1]].
- algotom.util.utility.generate_spiral_positions(step, num_pos, height, width, spiral_shape=1.0)[source]¶
Generate Fermat spiral positions. Unit is pixel.
- Parameters
step (int) – Step size in pixel. ~ 20-> 40
num_pos (int) – Number of positions.
height (int) – Height of the field of view (in pixel).
width (int) – Width of the field of view (in pixel).
spiral_shape (float, optional) – To define the spiral shape.
- Returns
array_like – 2D array. List of (x,y) positions
- algotom.util.utility.find_center_visual_sinograms(sino_180, output, start, stop, step=1, zoom=1.0, display=False)[source]¶
For visually finding the center-of-rotation (COR) using converted 360-degree sinograms from a 180-degree sinogram at different CORs (Ref. [1]).
- Parameters
sino_180 (array_like) – 2D array. 180-degree sinogram.
output (str) – Base folder for saving converted 360-degree sinograms.
start (float) – Starting point for searching CoR.
stop (float) – Ending point for searching CoR.
step (float) – Searching step.
zoom (float) – To resize output images. For example, 0.5 <=> reduce the size of output images by half.
display (bool) – Print the output if True.
- Returns
str – Folder path to tif images.
References