7.2.1. algotom.prep.calculation
¶
Module of calculation methods in the preprocessing stage:
Calculating the center-of-rotation (COR) using a 180-degree sinogram.
Determining the overlap-side and overlap-area between images.
Calculating the COR in a half-acquisition scan (360-degree scan with offset COR).
Using the similar technique as above to calculate the COR in a 180-degree scan from two projections.
Determining the relative translations between images using phase-correlation technique.
Calculating the COR using phase-correlation technique.
Functions:
|
Generate a double-wedge binary mask using Eq. |
|
Calculate a metric of an estimated center-of-rotation. |
|
Find the center-of-rotation (COR) using integer shifting. |
|
Find the center-of-rotation (COR) using sub-pixel shifting. |
|
Downsample an image by averaging. |
|
Find the center-of-rotation using the method described in Ref. |
|
Calculate the curvature of a fitted curve going through the minimum value of a metric list. |
|
Calculate the correlation metric. |
|
Calculate the correlation metrics between a rectangular region, defined by the window width, on the utmost left/right side of image 2 and the same size region in image 1 where the region is slided across image 1. |
|
Find the overlap area and overlap side between two images (Ref. |
|
Find the overlap-areas and overlap-sides of a list of images where the overlap side referring to the previous image. |
|
Find the center-of-rotation (COR) in a 360-degree scan with offset COR use the method presented in Ref. |
|
Return complex gradient of a 2D array. |
|
Find relative translation in x and y direction between images with haft-pixel accuracy (Ref. |
|
Find the center-of-rotation (COR) using projection images at 0-degree and 180-degree. |
|
Find the center-of-rotation (COR) using projection images at 0-degree and 180-degree based on a method in Ref. |
|
Calculate reconstructable height in a helical scan. |
|
Calculate the maximum index of a reconstructable slice in a helical scan. |
- algotom.prep.calculation.make_inverse_double_wedge_mask(height, width, radius, hor_drop=None, ver_drop=None)[source]¶
Generate a double-wedge binary mask using Eq. (3) in Ref. [1]. Values outside the double-wedge region correspond to 1.0.
- Parameters
height (int) – Image height.
width (int) – Image width.
radius (int) – Radius of an object, in pixel unit.
hor_drop (int or None, optional) – Number of rows (2 * hor_drop) around the middle of the mask with values set to zeros.
ver_drop (int or None, optional) – Number of columns (2 * ver_drop) around the middle of the mask with values set to zeros.
- Returns
array_like – 2D binary mask.
References
- algotom.prep.calculation.calculate_center_metric(center, sino_180, sino_flip, sino_comp, mask)[source]¶
Calculate a metric of an estimated center-of-rotation.
- Parameters
center (float) – Estimated center.
sino_180 (array_like) – 2D array. 180-degree sinogram.
sino_flip (array_like) – 2D array. Flip the 180-degree sinogram in the left/right direction.
sino_comp (array_like) – 2D array. Used to fill the gap left by image shifting.
mask (array_like) – 2D array. Used to select coefficients in the double-wedge region.
- Returns
float – Metric.
- algotom.prep.calculation.coarse_search_cor(sino_180, start, stop, ratio=0.5, denoise=True, ncore=None, hor_drop=None, ver_drop=None)[source]¶
Find the center-of-rotation (COR) using integer shifting.
- Parameters
sino_180 (array_like) – 2D array. 180-degree sinogram.
start (int) – Starting point for searching COR.
stop (int) – Ending point for searching COR.
ratio (float) – Ratio between a sample and the width of the sinogram.
denoise (bool, optional) – Apply a smoothing filter.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
hor_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
ver_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
- Returns
float – Center of rotation.
- algotom.prep.calculation.fine_search_cor(sino_180, start, radius, step, ratio=0.5, denoise=True, ncore=None, hor_drop=None, ver_drop=None)[source]¶
Find the center-of-rotation (COR) using sub-pixel shifting.
- Parameters
sino_180 (array_like) – 2D array. 180-degree sinogram.
start (float) – Starting point for searching COR.
radius (float) – Searching range: [start - radius; start + radius].
step (float) – Searching step.
ratio (float) – Ratio between a sample and the width of the sinogram.
denoise (bool, optional) – Apply a smoothing filter.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
hor_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
ver_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
- Returns
float – Center of rotation.
- algotom.prep.calculation.downsample_cor(image, dsp_fact0, dsp_fact1)[source]¶
Downsample an image by averaging.
- Parameters
image (array_like) – 2D array.
dsp_fact0 (int) – Downsampling factor along axis 0.
dsp_fact1 (int) – Downsampling factor along axis 1.
- Returns
array_like – 2D array. Downsampled image.
- algotom.prep.calculation.find_center_vo(sino_180, start=None, stop=None, step=0.25, radius=4, ratio=0.5, dsp=True, ncore=None, hor_drop=None, ver_drop=None)[source]¶
Find the center-of-rotation using the method described in Ref. [1].
- Parameters
sino_180 (array_like) – 2D array. 180-degree sinogram.
start (float) – Starting point for searching CoR. Use the value of (width/2 - width/16) if None.
stop (float) – Ending point for searching CoR. Use the value of (width/2 + width/16) if None.
step (float) – Sub-pixel accuracy of estimated CoR.
radius (float) – Searching range with the sub-pixel step.
ratio (float) – Ratio between the sample and the width of the sinogram.
dsp (bool) – Enable/disable downsampling.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
hor_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
ver_drop (int or None, optional) – Refer the method of “make_inverse_double_wedge_mask”
- Returns
float – Center-of-rotation.
References
- algotom.prep.calculation.calculate_curvature(list_metric)[source]¶
Calculate the curvature of a fitted curve going through the minimum value of a metric list.
- Parameters
list_metric (array_like) – 1D array. List of metrics.
- Returns
curvature (float) – Quadratic coefficient of the parabola fitting.
min_pos (float) – Position of the minimum value with sub-pixel accuracy.
- algotom.prep.calculation.correlation_metric(mat1, mat2)[source]¶
Calculate the correlation metric. Smaller metric corresponds to better correlation.
- Parameters
mat1 (array_like)
mat2 (array_like)
- Returns
float – Correlation metric.
- algotom.prep.calculation.search_overlap(mat1, mat2, win_width, side, denoise=True, norm=False, use_overlap=False, ncore=None)[source]¶
Calculate the correlation metrics between a rectangular region, defined by the window width, on the utmost left/right side of image 2 and the same size region in image 1 where the region is slided across image 1.
- Parameters
mat1 (array_like) – 2D array. Projection image or sinogram image.
mat2 (array_like) – 2D array. Projection image or sinogram image.
win_width (int) – Width of the searching window.
side ({0, 1}) – Only two options: 0 or 1. It is used to indicate the overlap side respects to image 1. “0” corresponds to the left side. “1” corresponds to the right side.
denoise (bool, optional) – Apply the Gaussian filter if True.
norm (bool, optional) – Apply the normalization if True.
use_overlap (bool, optional) – Use the combination of images in the overlap area for calculating correlation coefficients if True.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
list_metric (array_like) – 1D array. List of the correlation metrics.
offset (int) – Initial position of the searching window where the position corresponds to the center of the window.
- algotom.prep.calculation.find_overlap(mat1, mat2, win_width, side=None, denoise=True, norm=False, use_overlap=False, ncore=None)[source]¶
Find the overlap area and overlap side between two images (Ref. [1]) where the overlap side referring to the first image.
- Parameters
mat1 (array_like) – 2D array. Projection image or sinogram image.
mat2 (array_like) – 2D array. Projection image or sinogram image.
win_width (int) – Width of the searching window.
side ({None, 0, 1}, optional) – Only there options: None, 0, or 1. “None” corresponding to fully automated determination. “0” corresponding to the left side. “1” corresponding to the right side.
denoise (bool, optional) – Apply the Gaussian filter if True.
norm (bool, optional) – Apply the normalization if True.
use_overlap (bool, optional) – Use the combination of images in the overlap area for calculating correlation coefficients if True.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
overlap (float) – Width of the overlap area between two images.
side (int) – Overlap side between two images.
overlap_position (float) – Position of the window in the first image giving the best correlation metric.
References
- algotom.prep.calculation.find_overlap_multiple(list_mat, win_width, side=None, denoise=True, norm=False, use_overlap=False, ncore=None)[source]¶
Find the overlap-areas and overlap-sides of a list of images where the overlap side referring to the previous image.
- Parameters
list_mat (list of array_like) – List of 2D array. Projection image or sinogram image.
win_width (int) – Width of the searching window.
side ({None, 0, 1}, optional) – Only there options: None, 0, or 1. “None” corresponding to fully automated determination. “0” corresponding to the left side. “1” corresponding to the right side.
denoise (bool, optional) – Apply the Gaussian filter if True.
norm (bool, optional) – Apply the normalization if True.
use_overlap (bool, optional) – Use the combination of images in the overlap area for calculating correlation coefficients if True.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
list_overlap (list of tuple of floats) – List of [overlap, side, overlap_position]. overlap : Width of the overlap area between two images. side : Overlap side between two images. overlap_position : Position of the window in the first image giving the best correlation metric.
- algotom.prep.calculation.find_center_360(sino_360, win_width, side=None, denoise=True, norm=False, use_overlap=False, ncore=None)[source]¶
Find the center-of-rotation (COR) in a 360-degree scan with offset COR use the method presented in Ref. [1].
- Parameters
sino_360 (array_like) – 2D array. 360-degree sinogram.
win_width (int) – Window width used for finding the overlap area.
side ({None, 0, 1}, optional) – Overlap size. Only there options: None, 0, or 1. “None” corresponding to fully automated determination. “0” corresponding to the left side. “1” corresponding to the right side.
denoise (bool, optional) – Apply the Gaussian filter if True.
norm (bool, optional) – Apply the normalization if True.
use_overlap (bool, optional) – Use the combination of images in the overlap area for calculating correlation coefficients if True.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
cor (float) – Center-of-rotation.
overlap (float) – Width of the overlap area between two halves of the sinogram.
side (int) – Overlap side between two halves of the sinogram.
overlap_position (float) – Position of the window in the first image giving the best correlation metric.
References
- algotom.prep.calculation.find_shift_based_phase_correlation(mat1, mat2, gradient=True)[source]¶
Find relative translation in x and y direction between images with haft-pixel accuracy (Ref. [1]).
- Parameters
mat1 (array_like) – 2D array. Projection image or sinogram image.
mat2 (array_like) – 2D array. Projection image or sinogram image.
gradient (bool, optional) – Use the complex gradient of the input image for calculation.
- Returns
ty (float) – Translation in y-direction.
tx (float) – Translation in x-direction.
References
- algotom.prep.calculation.find_center_based_phase_correlation(mat1, mat2, flip=True, gradient=True)[source]¶
Find the center-of-rotation (COR) using projection images at 0-degree and 180-degree.
- Parameters
mat1 (array_like) – 2D array. Projection image at 0-degree.
mat2 (array_like) – 2D array. Projection image at 180-degree.
flip (bool, optional) – Flip the 180-degree projection in the left-right direction if True.
gradient (bool, optional) – Use the complex gradient of the input image for calculation.
- Returns
cor (float) – Center-of-rotation.
- algotom.prep.calculation.find_center_projection(mat1, mat2, flip=True, win_width=None, chunk_height=None, start_row=None, denoise=True, norm=False, use_overlap=False, ncore=None)[source]¶
Find the center-of-rotation (COR) using projection images at 0-degree and 180-degree based on a method in Ref. [1].
- Parameters
mat1 (array_like) – 2D array. Projection image at 0-degree.
mat2 (array_like) – 2D array. Projection image at 180-degree.
flip (bool, optional) – Flip the 180-degree projection in the left-right direction if True.
win_width (int, optional) – Width of the searching window.
chunk_height (int or float, optional) – Height of the sub-area of projection images. If a float is given, it must be in the range of [0.0, 1.0].
start_row (int, optional) – Starting row used to extract the sub-area.
denoise (bool, optional) – Apply the Gaussian filter if True.
norm (bool, optional) – Apply the normalization if True.
use_overlap (bool, optional) – Use the combination of images in the overlap area for calculating correlation coefficients if True.
ncore (int or None) – Number of cpu-cores used for computing. Automatically selected if None.
- Returns
cor (float) – Center-of-rotation.
References
- algotom.prep.calculation.calculate_reconstructable_height(y_start, y_stop, pitch, scan_type)[source]¶
Calculate reconstructable height in a helical scan.
- Parameters
y_start (float) – Y-position of the stage at the beginning of the scan.
y_stop (float) – Y-position of the stage at the end of the scan.
pitch (float) – The distance which the y-stage is translated in one full rotation.
scan_type ({“180”, “360”}) – One of two options: “180” for generating a 180-degree sinogram or “360” for generating a 360-degree sinogram.
- Returns
y_s (float) – Starting point of the reconstructable height.
y_e (float) – End point of the reconstructable height.
- algotom.prep.calculation.calculate_maximum_index(y_start, y_stop, pitch, pixel_size, scan_type)[source]¶
Calculate the maximum index of a reconstructable slice in a helical scan.
- Parameters
y_start (float) – Y-position of the stage at the beginning of the scan.
y_stop (float) – Y-position of the stage at the end of the scan.
pitch (float) – The distance which the y-stage is translated in one full rotation.
pixel_size (float) – Pixel size. The unit must be the same as y-position.
scan_type ({“180”, “360”}) – One of two options: “180” for generating a 180-degree sinogram or “360” for generating a 360-degree sinogram.
- Returns
int – Maximum index of reconstructable slices.