tomopy.recon.rotation
¶
Module for functions related to finding axis of rotation.
Functions:
|
Find rotation axis location. |
|
Find rotation axis location using Nghia Vo's method. |
|
Find rotation axis location by finding the offset between the first projection and a mirrored projection 180 degrees apart using phase correlation in Fourier space. |
|
Save images reconstructed with a range of rotation centers. |
- tomopy.recon.rotation.find_center(tomo, theta, ind=None, init=None, tol=0.5, mask=True, ratio=1.0, sinogram_order=False)[source]¶
Find rotation axis location.
The function exploits systematic artifacts in reconstructed images due to shifts in the rotation center. It uses image entropy as the error metric and ‘’Nelder-Mead’’ routine (of the scipy optimization module) as the optimizer [B12].
- Parameters
tomo (ndarray) – 3D tomographic data.
theta (array) – Projection angles in radian.
ind (int, optional) – Index of the slice to be used for reconstruction.
init (float) – Initial guess for the center.
tol (scalar) – Desired sub-pixel accuracy.
mask (bool, optional) – If
True
, apply a circular mask to the reconstructed image to limit the analysis into a circular region.ratio (float, optional) – The ratio of the radius of the circular mask to the edge of the reconstructed image.
sinogram_order (bool, optional) – Determins whether data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False, theta first axis).
- Returns
float – Rotation axis location.
- tomopy.recon.rotation.find_center_pc(proj1, proj2, tol=0.5, rotc_guess=None)[source]¶
Find rotation axis location by finding the offset between the first projection and a mirrored projection 180 degrees apart using phase correlation in Fourier space. The
phase_cross_correlation
function uses cross-correlation in Fourier space, optionally employing an upsampled matrix-multiplication DFT to achieve arbitrary subpixel precision. [B16].- Parameters
proj1 (ndarray) – 2D projection data.
proj2 (ndarray) – 2D projection data.
tol (scalar, optional) – Subpixel accuracy
rotc_guess (float, optional) – Initual guess value for the rotation center
- Returns
float – Rotation axis location.
- tomopy.recon.rotation.find_center_vo(tomo, ind=None, smin=-50, smax=50, srad=6, step=0.25, ratio=0.5, drop=20, ncore=None)[source]¶
Find rotation axis location using Nghia Vo’s method. [B25].
- Parameters
tomo (ndarray) – 3D tomographic data or a 2D sinogram.
ind (int, optional) – Index of the slice to be used for reconstruction.
smin, smax (int, optional) – Coarse search radius. Reference to the horizontal center of the sinogram.
srad (float, optional) – Fine search radius.
step (float, optional) – Step of fine searching.
ratio (float, optional) – The ratio between the FOV of the camera and the size of object. It’s used to generate the mask.
drop (int, optional) – Drop lines around vertical center of the mask.
ncore (int, optional) – Number of cores that will be assigned to jobs.
- Returns
float – Rotation axis location.
- tomopy.recon.rotation.mask_empty_slice(tomo, threshold=0.25)[source]¶
Generate a mask to indicate whether current slice contains sample
At APS 1ID, some of the projection images contains large empty area above the sample, resulting in empty layers.
- Parameters
tomo (ndarray) – 3D tomographic data.
threshold (float, optional) – determine whether a layer is considered to be empty
- Returns
nparray – a mask indicate the emptyness of each layer
- tomopy.recon.rotation.write_center(tomo, theta, dpath='tmp/center', cen_range=None, ind=None, mask=False, ratio=1.0, sinogram_order=False, algorithm='gridrec', filter_name='parzen')[source]¶
Save images reconstructed with a range of rotation centers.
Helps finding the rotation center manually by visual inspection of images reconstructed with a set of different centers.The output images are put into a specified folder and are named by the center position corresponding to the image.
- Parameters
tomo (ndarray) – 3D tomographic data.
theta (array) – Projection angles in radian.
dpath (str, optional) – Folder name to save output images.
cen_range (list, optional) – [start, end, step] Range of center values.
ind (int, optional) – Index of the slice to be used for reconstruction.
mask (bool, optional) – If
True
, apply a circular mask to the reconstructed image to limit the analysis into a circular region.ratio (float, optional) – The ratio of the radius of the circular mask to the edge of the reconstructed image.
sinogram_order (bool, optional) – Determins whether data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False, theta first axis).
algorithm ({str, function}) – One of the following string values.
- ‘art’
Algebraic reconstruction technique [B2].
- ‘bart’
Block algebraic reconstruction technique.
- ‘fbp’
Filtered back-projection algorithm.
- ‘gridrec’
- ‘mlem’
Maximum-likelihood expectation maximization algorithm [B3].
- ‘osem’
Ordered-subset expectation maximization algorithm [B18].
- ‘ospml_hybrid’
Ordered-subset penalized maximum likelihood algorithm with weighted linear and quadratic penalties.
- ‘ospml_quad’
Ordered-subset penalized maximum likelihood algorithm with quadratic penalties.
- ‘pml_hybrid’
Penalized maximum likelihood algorithm with weighted linear and quadratic penalties [B19].
- ‘pml_quad’
Penalized maximum likelihood algorithm with quadratic penalty.
- ‘sirt’
Simultaneous algebraic reconstruction technique.
- ‘tv’
Total Variation reconstruction technique [B8].
- ‘grad’
Gradient descent method with a constant step size
- ‘tikh’
Tikhonov regularization with identity Tikhonov matrix.
- filter_namestr, optional
Name of the filter for analytic reconstruction.
- ‘none’
No filter.
- ‘shepp’
Shepp-Logan filter (default).
- ‘cosine’
Cosine filter.
- ‘hann’
Cosine filter.
- ‘hamming’
Hamming filter.
- ‘ramlak’
Ram-Lak filter.
- ‘parzen’
Parzen filter.
- ‘butterworth’
Butterworth filter.
- ‘custom’
A numpy array of size next_power_of_2(num_detector_columns)/2 specifying a custom filter in Fourier domain. The first element of the filter should be the zero-frequency component.
- ‘custom2d’
A numpy array of size num_projections*next_power_of_2(num_detector_columns)/2 specifying a custom angle-dependent filter in Fourier domain. The first element of each filter should be the zero-frequency component.