tomopy.prep.alignment
¶
Functions:
|
Aligns the projection image stack using the sequential re-projection algorithm [B9]. |
|
Aligns the projection image stack using the joint re-projection algorithm [B9]. |
|
Simulates jitter in projection images. |
|
Adds Gaussian noise with zero mean and a given standard deviation as a ratio of the maximum value in data. |
|
Blurs the edge of the projection images. |
|
Shift projections images for a given set of shift values in horizontal and vertical directions. |
|
Linearly scales the projection images in the range between -1 and 1. |
|
Tilt object at a given angle from the rotation axis. |
|
Apply distortion correction to projections using the polynomial model. |
|
Generate an unwarped sinogram of a 3D tomographic data using the polynomial model. |
-
tomopy.prep.alignment.
add_jitter
(prj, low=0, high=1)[source]¶ Simulates jitter in projection images. The jitter is simulated by drawing random samples from a uniform distribution over the half-open interval [low, high).
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
low (float, optional) – Lower boundary of the output interval. All values generated will be greater than or equal to low. The default value is 0.
high (float) – Upper boundary of the output interval. All values generated will be less than high. The default value is 1.0.
- Returns
ndarray – 3D stack of projection images with jitter.
-
tomopy.prep.alignment.
add_noise
(prj, ratio=0.05)[source]¶ Adds Gaussian noise with zero mean and a given standard deviation as a ratio of the maximum value in data.
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
ratio (float, optional) – Ratio of the standard deviation of the Gaussian noise distribution to the maximum value in data.
- Returns
ndarray – 3D stack of projection images with added Gaussian noise.
-
tomopy.prep.alignment.
align_joint
(prj, ang, fdir='.', iters=10, pad=0, 0, blur=True, center=None, algorithm='sirt', upsample_factor=10, rin=0.5, rout=0.8, save=False, debug=True)[source]¶ Aligns the projection image stack using the joint re-projection algorithm [B9].
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
ang (ndarray) – Projection angles in radians as an array.
iters (scalar, optional) – Number of iterations of the algorithm.
pad (list-like, optional) – Padding for projection images in x and y-axes.
blur (bool, optional) – Blurs the edge of the image before registration.
center (array, optional) – Location of rotation axis.
algorithm ({str, function}) – One of the following string values.
upsample_factor (integer, optional) – The upsampling factor. Registration accuracy is inversely propotional to upsample_factor.
rin (scalar, optional) – The inner radius of blur function. Pixels inside rin is set to one.
rout (scalar, optional) – The outer radius of blur function. Pixels outside rout is set to zero.
save (bool, optional) – Saves projections and corresponding reconstruction for each algorithm iteration.
debug (book, optional) – Provides debugging info such as iterations and error.
- Returns
ndarray – 3D stack of projection images with jitter.
ndarray – Error array for each iteration.
-
tomopy.prep.alignment.
align_seq
(prj, ang, fdir='.', iters=10, pad=0, 0, blur=True, center=None, algorithm='sirt', upsample_factor=10, rin=0.5, rout=0.8, save=False, debug=True)[source]¶ Aligns the projection image stack using the sequential re-projection algorithm [B9].
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
ang (ndarray) – Projection angles in radians as an array.
iters (scalar, optional) – Number of iterations of the algorithm.
pad (list-like, optional) – Padding for projection images in x and y-axes.
blur (bool, optional) – Blurs the edge of the image before registration.
center (array, optional) – Location of rotation axis.
algorithm ({str, function}) – One of the following string values.
upsample_factor (integer, optional) – The upsampling factor. Registration accuracy is inversely propotional to upsample_factor.
rin (scalar, optional) – The inner radius of blur function. Pixels inside rin is set to one.
rout (scalar, optional) – The outer radius of blur function. Pixels outside rout is set to zero.
save (bool, optional) – Saves projections and corresponding reconstruction for each algorithm iteration.
debug (book, optional) – Provides debugging info such as iterations and error.
- Returns
ndarray – 3D stack of projection images with jitter.
ndarray – Error array for each iteration.
-
tomopy.prep.alignment.
blur_edges
(prj, low=0, high=0.8)[source]¶ Blurs the edge of the projection images.
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
low (scalar, optional) – Min ratio of the blurring frame to the image size.
high (scalar, optional) – Max ratio of the blurring frame to the image size.
- Returns
ndarray – Edge-blurred 3D stack of projection images.
-
tomopy.prep.alignment.
calc_slit_box_aps_1id
(slit_box_corners, inclip=1, 10, 1, 10)[source]¶ Calculate the clip box based on given slip corners.
- Parameters
slit_box_corners (np.ndarray) – Four corners of the slit box as a 4x2 matrix
inclip (tuple, optional) – Extra inclipping to avoid clipping artifacts
- Returns
Tuple – Cliping indices as a tuple of four (clipFromTop, clipToBottom, clipFromLeft, clipToRight)
-
tomopy.prep.alignment.
distortion_correction_proj
(tomo, xcenter, ycenter, list_fact, ncore=None, nchunk=None)[source]¶ Apply distortion correction to projections using the polynomial model. Coefficients are calculated using Vounwarp package.:cite:Vo:15
- Parameters
tomo (ndarray) – 3D tomographic data.
xcenter (float) – Center of distortion in x-direction. From the left of the image.
ycenter (float) – Center of distortion in y-direction. From the top of the image.
list_fact (list of floats) – Polynomial coefficients of the backward model.
ncore (int, optional) – Number of cores that will be assigned to jobs.
nchunk (int, optional) – Chunk size for each core.
- Returns
ndarray – Corrected 3D tomographic data.
-
tomopy.prep.alignment.
distortion_correction_sino
(tomo, ind, xcenter, ycenter, list_fact)[source]¶ Generate an unwarped sinogram of a 3D tomographic data using the polynomial model. Coefficients are calculated using Vounwarp package [B26]
- Parameters
tomo (ndarray) – 3D tomographic data.
ind (int) – Index of the unwarped sinogram.
xcenter (float) – Center of distortion in x-direction. From the left of the image.
ycenter (float) – Center of distortion in y-direction. From the top of the image.
list_fact (list of floats) – Polynomial coefficients of the backward model.
- Returns
2D array – Corrected sinogram.
-
tomopy.prep.alignment.
find_slits_corners_aps_1id
(img, method='quadrant+', medfilt2_kernel_size=3, medfilt_kernel_size=23)[source]¶ Automatically locate the slit box location by its four corners.
NOTE: The four slits that form a binding box is the current setup at aps_1id, which reduce the illuminated region on the detector. Since the slits are stationary, they can serve as a reference to check detector drifting during the scan. Technically, the four slits should be used to find the transformation matrix (not necessarily affine) to correct the image. However, since we are dealing with 2D images with very little distortion, affine transformation matrices were used for approximation. Therefore the “four corners” are used instead of all four slits.
- Parameters
img (np.ndarray) – 2D images
method (str, [‘simple’, ‘quadrant’, ‘quadrant+’], optional) –
- method for auto detecting slit corners
- simple :: assume a rectange slit box, fast but less accurate
(1 pixel precision)
- quadrant :: subdivide the image into four quandrant, then use
an explicit method to find the corner (1 pixel precision)
- quadrant+ :: similar to quadrant, but use curve_fit (gauss1d) to
find the corner (0.1 pixel precision)
medfilt2_kernel_size (int, optional) – 2D median filter kernel size for noise reduction
medfilt_kernel_size (int, optional) – 1D median filter kernel size for noise reduction
- Returns
tuple – autodetected slit corners (counter-clockwise order) (upperLeft, lowerLeft, lowerRight, upperRight)
-
tomopy.prep.alignment.
load_distortion_coefs
(file_path)[source]¶ Load distortion coefficients from a text file. Order of the infor in the text file: xcenter ycenter factor_0 factor_1 factor_2 ..
- Parameters
file_path (Path to the file.)
- Returns
Tuple of (xcenter, ycenter, list_fact).
-
tomopy.prep.alignment.
remove_slits_aps_1id
(imgstacks, slit_box_corners, inclip=1, 10, 1, 10)[source]¶ Remove the slits from still images
- Parameters
imgstacks (np.ndarray) – tomopy images stacks (axis_0 is the oemga direction)
slit_box_corners (np.ndarray) – four corners of the slit box
inclip (tuple, optional) – Extra inclipping to avoid clipping artifacts
- Returns
np.ndarray – tomopy images stacks without regions outside slits
-
tomopy.prep.alignment.
scale
(prj)[source]¶ Linearly scales the projection images in the range between -1 and 1.
- Parameters
prj (ndarray) – 3D stack of projection images. The first dimension is projection axis, second and third dimensions are the x- and y-axes of the projection image, respectively.
- Returns
ndarray – Scaled 3D stack of projection images.
-
tomopy.prep.alignment.
shift_images
(prj, sx, sy)[source]¶ Shift projections images for a given set of shift values in horizontal and vertical directions.
-
tomopy.prep.alignment.
tilt
(obj, rad=0, phi=0)[source]¶ Tilt object at a given angle from the rotation axis.
Warning
Not implemented yet.
- Parameters
obj (ndarray) – 3D discrete object.
rad (scalar, optional) – Radius in polar cordinates to define tilt angle. The value is between 0 and 1, where 0 means no tilt and 1 means a tilt of 90 degrees. The tilt angle can be obtained by arcsin(rad).
phi (scalar, optional) – Angle in degrees to define tilt direction from the rotation axis. 0 degree means rotation in sagittal plane and 90 degree means rotation in coronal plane.
- Returns
ndarray – Tilted 3D object.