
TomoPy is an open-source Python package for tomographic data processing and image reconstruction.
Features¶
- Image reconstruction algorithms for tomography.
- Various filters, ring removal algorithms, phase retrieval algorithms.
- Forward projection operator for absorption and wave propagation.
Contribute¶
- Issue Tracker: https://github.com/tomopy/tomopy/issues
- Documentation: https://github.com/tomopy/tomopy/tree/master/doc
- Source Code: https://github.com/tomopy/tomopy/tree/master/tomopy
- Tests: https://github.com/tomopy/tomopy/tree/master/test
Table of Contents¶
About¶
Tomographic reconstruction creates three-dimensional views of an object by combining two-dimensional images taken from multiple directions, for example, this is how a CAT (computer-aided tomography) scanner generates 3D views of the heart or brain.
Data collection can be rapid, but the required computations are massive and often the beamline staff can be overwhelmed by data that are collected far faster than corrections and reconstruction can be performed [C15]. Further, many common experimental perturbations can degrade the quality of tomographs, unless corrections are applied.
To address the needs for image correction and tomographic reconstruction in an instrument independent manner, the TomoPy code was developed [A1], which is a parallelizable high performance reconstruction code.
Install directions¶
This section covers the basics of how to download and install TomoPy.
Contents:
Supported Environments¶
TomoPy is tested, built, and distributed for python 2.7 3.5 3.6 on Linux/macOS and python 3.5 3.6 on Windows 10.
Installing from Conda (Recommended)¶
If you only want to run TomoPy, not develop it, then you should install through a package manager. Conda, our supported package manager, can install TomoPy and its dependencies for you.
First, you must have Conda installed, then open a terminal or a command prompt window and run:
$ conda install -c conda-forge tomopy
This will install TomoPy and all the dependencies from the conda-forge channel.
Updating the installation¶
TomoPy is an active project, so we suggest you update your installation frequently. To update the installation run:
$ conda update -c conda-forge tomopy
For some more information about using Conda, please refer to the docs.
Installing from source with Conda¶
Sometimes an adventurous user may want to get the source code, which is always more up-to-date than the one provided by Conda (with more bugs of course!).
For this you need to get the source from the TomoPy repository on GitHub. Download the source to your local computer using git by opening a terminal and running:
$ git clone https://github.com/tomopy/tomopy.git
in the folder where you want the source code. This will create a folder called tomopy which contains a copy of the source code.
Installing dependencies¶
You will need to install all the dependencies listed in
requirements.txt
or meta.yaml
files. For example, requirements can be
installed using Conda by running:
$ conda install --file requirements.txt
After navigating to inside the tomopy directory, you can install TomoPy by building/compiling the shared libraries and running the install script:
$ python build.py
$ pip install .
Common issues¶
No issues with the current build system have been reported.
Importing TomoPy¶
When importing, it is best to import TomoPy before importing numpy. See this thread for details.
Tomographic data files¶
For reading tomography files formatted in different ways, please go check the DXchange package. There are various examples and demonstration scripts about how to use the package for loading your datasets.
The package can be installed by simply running the following in a terminal:
conda install -c conda-forge dxchange
For a repository of experimental and simulated data sets please check TomoBank [C6].
Development¶
This section explains the basics for developers who wish to contribute to the TomoPy project.
Contents:
Cloning the repository¶
The project is maintained on GitHub, which is a version control and a collaboration platform for software developers. To start first register on GitHub and fork the TomoPy repository by clicking the Fork button in the header of the TomoPy repository:

This creates a copy of the project in your personal GitHub space. The next thing you want to do is to clone it to your local machine. You can do this by clicking the Clone in Desktop button in the bottom of the right hand side bar:

This will launch the GitHub desktop application (available for both Mac and Win) and ask you where you want to save it. Select a location in your computer and feel comfortable with making modifications in the code.
Running the Tests¶
Tomopy has a suite of python unit tests that live in the /test
directory,
where they follow the same tree structure as the packages under /tomopy
.
These are automatically run by TravisCI when you make a pull request
(See below for how to do that) and you can run them manually using pytest,
or whichever python test runner you prefer. To make it easier to run tests
on the changes you make to the code, it is recommended that you install TomoPy
in development mode. (python setup.py develop)
The pytest test runner, is available through pip or anaconda.
To run the tests open a terminal, navigate to your project folder,
then run py.test
.
To run sections of tests, pass py.test a
directory or filepath, as in py.test test/test_recon
or
py.test test/test_recon/test_rotation.py
.
When writing tests, at minimum we try to check all function returns with synthetic data, together with some dimension, type, etc. Writing tests is highly encouraged!
Coding conventions¶
We try to keep our code consistent and readable. So, please keep in mind the following style and syntax guidance before you start coding.
First of all the code should be well documented, easy to understand, and integrate well into the rest of the project. For example, when you are writing a new function always describe the purpose and the parameters:
def my_awesome_func(a, b):
"""
Adds two numbers.
Parameters
----------
a : scalar (float)
First number to add
b : scalar (float)
Second number to add
Returns
-------
output : scalar (float)
Added value
"""
return a+b
Package versioning¶
We follow the X.Y.Z (Major.Minor.Patch) semantic for package versioning. The version should be updated before each pull request accordingly. The patch number is incremented for minor changes and bug fixes which do not change the software’s API. The minor version is incremented for releases which add new, but backward-compatible, API features, and the major version is incremented for API changes which are not backward-compatible. For example, software which relies on version 2.1.5 of an API is compatible with version 2.2.3, but not necessarily with 3.2.4.
Committing changes¶
After making some changes in the code, you may want to take a snapshot of the edits you made. That’s when you make a commit. To do this, launch the GitHub desktop application and it should provide you all the changes in your code since your last commit. Write a brief Summary and Description about the changes you made and click the Commit button:

You can continue to make changes, add modules, write your own functions, and take more Commit snapshots of your code writing process.
Contributing back¶
Once you feel that the functionality you added would benefit the community, then you should consider contributing back to the TomoPy project. You will need to push your local commits to GitHub, then go to your online GitHub repository of the project and click on the green button to compare, review, and create a pull request.

After clicking on this button, you are presented with a review page where you can get a high-level overview of what exactly has changed between your forked branch and the original TomoPy repository. When you’re ready to submit your pull request, click Create pull request:

Clicking on Create pull request sends you to a discussion page, where you can enter a title and optional description. It’s important to provide as much useful information and a rationale for why you’re making this Pull Request in the first place.
When you’re ready typing out your heartfelt argument, click on Send pull request. You’re done!
Release Notes¶
TomoPy 1.0.0 Release Notes¶
New features¶
- FFTW implementation is now adopted. All functions that rely on FFTs such as gridrec, phase retrieval, stripe removal, etc. are now using the FFTW implementation through PyFFTW.
sinogram_order
is added torecon
as an additional argument. It determines whether data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False). Default is False, but we plan to make it True in the upcoming release.- Reconstruction algorithms only copies data if necessary.
- Updated library to support new mproc and recon functions. The data is now passed in sinogram order to recon functions. Also updated tests.
ncores
andnchunks
are now independent.- Setting
nchunks
to zero removes the dimension. That allows for the functions work on 2D data rather than 3D data. - Sliced data are used so that each process only receives the data it needs. No more
istart
andiend
variables for setting up indices in parallel processes. - Functions will reuse sharedmem arrays if they can.
Deprecated features¶
- All data I/O related functions are deprecated. They are available through DXchange package.
- Removed fft.h and fft.c, they are now completely replaced with FFTW.
Backward incompatible changes¶
emission
argument is removed fromrecon
. After this change the tomographic image reconstruction algorithms always assume data to be normalized.
Contributors¶
- Arthur Glowacki (@aglowacki)
- Daniel Pelt (@dmpelt)
- Dake Feng (@dakefeng)
- Doga Gursoy (@dgursoy)
- Francesco De Carlo (@decarlof)
- Lin Jiao (@yxqd)
- Luis Barroso-Luque (@lbluque)
- Michael Sutherland (@michael-sutherland)
- Rafael Vescovi (@ravescovi)
- Thomas Caswell (@tacaswell)
- Pete R. Jemian (@prjemian)
- Wei Xu (@celiafish)
API reference¶
This section contains the API reference and usage information for TomoPy.
TomoPy Modules:
tomopy.misc.corr
¶
Module for data correction and masking functions.
Functions:
adjust_range (arr[, dmin, dmax]) |
Change dynamic range of values in an array. |
circ_mask (arr, axis[, ratio, val, ncore]) |
Apply circular mask to a 3D array. |
gaussian_filter (arr[, sigma, order, axis, ncore]) |
Apply Gaussian filter to 3D array along specified axis. |
median_filter (arr[, size, axis, ncore]) |
Apply median filter to 3D array along specified axis. |
median_filter_cuda (arr[, size, axis]) |
Apply median filter to 3D array along 0 axis with GPU support. |
sobel_filter (arr[, axis, ncore]) |
Apply Sobel filter to 3D array along specified axis. |
remove_nan (arr[, val, ncore]) |
Replace NaN values in array with a given value. |
remove_neg (arr[, val, ncore]) |
Replace negative values in array with a given value. |
remove_outlier (arr, dif[, size, axis, …]) |
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_outlier_cuda (arr, dif[, size, axis]) |
Remove high intensity bright spots from a 3D array along axis 0 dimension using GPU. |
remove_ring (rec[, center_x, center_y, …]) |
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.
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.
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.
-
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=u'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) – minimum angle in degrees (int) to be considered ring artifact
- 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.morph
¶
Module for data size morphing functions.
Functions:
downsample (arr[, level, axis]) |
Downsample along specified axis of a 3D array. |
upsample (arr[, level, axis]) |
Upsample along specified axis of a 3D array. |
pad (arr, axis[, npad, mode, ncore]) |
Pad an array along specified axis. |
sino_360_t0_180 (data[, overlap, rotation]) |
Converts 0-360 degrees sinogram to a 0-180 sinogram. |
trim_sinogram (data, center, x, y, diameter) |
Provide sinogram corresponding to a circular region of interest by trimming the complete sinogram of a compact object. |
-
tomopy.misc.morph.
downsample
(arr, level=1, axis=2)[source]¶ Downsample along specified axis of a 3D array.
Parameters: - arr (ndarray) – 3D input array.
- level (int, optional) – Downsampling level in powers of two.
- axis (int, optional) – Axis along which downsampling will be performed.
Returns: ndarray – Downsampled 3D array in float32.
-
tomopy.misc.morph.
upsample
(arr, level=1, axis=2)[source]¶ Upsample along specified axis of a 3D array.
Parameters: - arr (ndarray) – 3D input array.
- level (int, optional) – Downsampling level in powers of two.
- axis (int, optional) – Axis along which upsampling will be performed.
Returns: ndarray – Upsampled 3D array in float32.
-
tomopy.misc.morph.
pad
(arr, axis, npad=None, mode=u'constant', ncore=None, **kwargs)[source]¶ Pad an array along specified axis.
Parameters: arr (ndarray) – Input array.
npad (int, optional) – New dimension after padding.
axis (int, optional) – Axis along which padding will be performed.
mode (str or function) – One of the following string values or a user supplied function.
- ‘constant’
Pads with a constant value.
- ‘edge’
Pads with the edge values of array.
constant_values (float, optional) – Used in ‘constant’. Pad value
ncore (int, optional) – Number of cores that will be assigned to jobs.
Returns: ndarray – Padded 3D array.
-
tomopy.misc.morph.
sino_360_to_180
(data, overlap=0, rotation=u'left')[source]¶ Converts 0-360 degrees sinogram to a 0-180 sinogram.
If the number of projections in the input data is odd, the last projection will be discarded.
Parameters: - data (ndarray) – Input 3D data.
- overlap (scalar, optional) – Overlapping number of pixels.
- rotation (string, optional) – Left if rotation center is close to the left of the field-of-view, right otherwise.
Returns: ndarray – Output 3D data.
-
tomopy.misc.morph.
sino_360_t0_180
(data, overlap=0, rotation=u'left')¶ Converts 0-360 degrees sinogram to a 0-180 sinogram.
If the number of projections in the input data is odd, the last projection will be discarded.
Parameters: - data (ndarray) – Input 3D data.
- overlap (scalar, optional) – Overlapping number of pixels.
- rotation (string, optional) – Left if rotation center is close to the left of the field-of-view, right otherwise.
Returns: ndarray – Output 3D data.
-
tomopy.misc.morph.
trim_sinogram
(data, center, x, y, diameter)[source]¶ Provide sinogram corresponding to a circular region of interest by trimming the complete sinogram of a compact object.
Parameters: - data (ndarray) – Input 3D data.
- center (float) – Rotation center location.
- x, y (int, int) – x and y coordinates in pixels (image center is (0, 0))
- diameter (float) – Diameter of the circle of the region of interest.
Returns: ndarray – Output 3D data.
tomopy.misc.phantom
¶
Module for generating synthetic phantoms.
Functions:
baboon ([size, dtype]) |
Load test baboon image array. |
barbara ([size, dtype]) |
Load test Barbara image array. |
cameraman ([size, dtype]) |
Load test cameraman image array. |
checkerboard ([size, dtype]) |
Load test checkerboard image array. |
lena ([size, dtype]) |
Load test Lena image array. |
peppers ([size, dtype]) |
Load test peppers image array. |
phantom (size, params[, dtype]) |
Generate a cube of given size using a list of ellipsoid parameters. |
shepp2d ([size, dtype]) |
Load test Shepp-Logan image array. |
shepp3d ([size, dtype]) |
Load 3D Shepp-Logan image array. |
-
tomopy.misc.phantom.
baboon
(size=512, dtype=u'float32')[source]¶ Load test baboon image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
barbara
(size=512, dtype=u'float32')[source]¶ Load test Barbara image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
cameraman
(size=512, dtype=u'float32')[source]¶ Load test cameraman image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
checkerboard
(size=512, dtype=u'float32')[source]¶ Load test checkerboard image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
lena
(size=512, dtype=u'float32')[source]¶ Load test Lena image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
peppers
(size=512, dtype=u'float32')[source]¶ Load test peppers image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
shepp2d
(size=512, dtype=u'float32')[source]¶ Load test Shepp-Logan image array.
Parameters: - size (int or tuple of int, optional) – Size of the output image.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
shepp3d
(size=128, dtype=u'float32')[source]¶ Load 3D Shepp-Logan image array.
Parameters: - size (int or tuple, optional) – Size of the 3D data.
- dtype (str, optional) – The desired data-type for the array.
Returns: ndarray – Output 3D test image.
-
tomopy.misc.phantom.
phantom
(size, params, dtype=u'float32')[source]¶ Generate a cube of given size using a list of ellipsoid parameters.
Parameters: - size (tuple of int) – Size of the output cube.
- params (list of dict) – List of dictionaries with the parameters defining the ellipsoids to include in the cube.
- dtype (str, optional) – Data type of the output ndarray.
Returns: ndarray – 3D object filled with the specified ellipsoids.
tomopy.prep.alignment
¶
Functions:
align_seq (prj, ang[, fdir, iters, pad, …]) |
Aligns the projection image stack using the sequential re-projection algorithm [C8]. |
align_joint (prj, ang[, fdir, iters, pad, …]) |
Aligns the projection image stack using the joint re-projection algorithm [C8]. |
add_jitter (prj[, low, high]) |
Simulates jitter in projection images. |
add_noise (prj[, ratio]) |
Adds Gaussian noise with zero mean and a given standard deviation as a ratio of the maximum value in data. |
blur_edges (prj[, low, high]) |
Blurs the edge of the projection images. |
shift_images (prj, sx, sy) |
Shift projections images for a given set of shift values in horizontal and vertical directions. |
scale (prj) |
Linearly scales the projection images in the range between -1 and 1. |
tilt (obj[, rad, phi]) |
Tilt object at a given angle from the rotation axis. |
-
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 [C8].
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_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 [C8].
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.
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.
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.
-
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.
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.normalize
¶
Module for data normalization.
Functions:
minus_log (arr[, ncore, out]) |
Computation of the minus log of a given array. |
normalize (arr, flat, dark[, cutoff, ncore, out]) |
Normalize raw projection data using the flat and dark field projections. |
normalize_bg (tomo[, air, ncore, nchunk]) |
Normalize 3D tomgraphy data based on background intensity. |
normalize_nf (tomo, flats, dark, flat_loc[, …]) |
Normalize raw 3D projection data with flats taken more than once during tomography. |
normalize_roi (arr[, roi, ncore]) |
Normalize raw projection data using an average of a selected window on projection images. |
-
tomopy.prep.normalize.
minus_log
(arr, ncore=None, out=None)[source]¶ Computation of the minus log of a given array.
Parameters: - arr (ndarray) – 3D stack of projections.
- 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 – Minus-log of the input data.
-
tomopy.prep.normalize.
normalize
(arr, flat, dark, cutoff=None, ncore=None, out=None)[source]¶ Normalize raw projection data using the flat and dark field projections.
Parameters: - arr (ndarray) – 3D stack of projections.
- flat (ndarray) – 3D flat field data.
- dark (ndarray) – 3D dark field data.
- cutoff (float, optional) – Permitted maximum vaue for the normalized data.
- 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 – Normalized 3D tomographic data.
-
tomopy.prep.normalize.
normalize_bg
(tomo, air=1, ncore=None, nchunk=None)[source]¶ Normalize 3D tomgraphy data based on background intensity.
Weight sinogram such that the left and right image boundaries (i.e., typically the air region around the object) are set to one and all intermediate values are scaled linearly.
Parameters: - tomo (ndarray) – 3D tomographic data.
- air (int, optional) – Number of pixels at each boundary to calculate the scaling factor.
- 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.normalize.
normalize_roi
(arr, roi=[0, 0, 10, 10], ncore=None)[source]¶ Normalize raw projection data using an average of a selected window on projection images.
Parameters: - arr (ndarray) – 3D tomographic data.
- roi (list of int, optional) – [top-left, top-right, bottom-left, bottom-right] pixel coordinates.
- ncore (int, optional) – Number of cores that will be assigned to jobs.
Returns: ndarray – Normalized 3D tomographic data.
-
tomopy.prep.normalize.
normalize_nf
(tomo, flats, dark, flat_loc, cutoff=None, ncore=None, out=None)[source]¶ Normalize raw 3D projection data with flats taken more than once during tomography. Normalization for each projection is done with the mean of the nearest set of flat fields (nearest flat fields).
Parameters: - tomo (ndarray) – 3D tomographic data.
- flats (ndarray) – 3D flat field data.
- dark (ndarray) – 3D dark field data.
- flat_loc (list of int) – Indices of flat field data within tomography
- 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 – Normalized 3D tomographic data.
tomopy.prep.phase
¶
Module for phase retrieval.
Functions:
retrieve_phase (tomo[, pixel_size, dist, …]) |
Perform single-step phase retrieval from phase-contrast measurements [C9]. |
-
tomopy.prep.phase.
retrieve_phase
(tomo, pixel_size=0.0001, dist=50, energy=20, alpha=0.001, pad=True, ncore=None, nchunk=None)[source]¶ Perform single-step phase retrieval from phase-contrast measurements [C9].
Parameters: - tomo (ndarray) – 3D tomographic data.
- pixel_size (float, optional) – Detector pixel size in cm.
- dist (float, optional) – Propagation distance of the wavefront in cm.
- energy (float, optional) – Energy of incident wave in keV.
- alpha (float, optional) – Regularization parameter.
- pad (bool, optional) – If True, extend the size of the projections by padding with zeros.
- ncore (int, optional) – Number of cores that will be assigned to jobs.
- nchunk (int, optional) – Chunk size for each core.
Returns: ndarray – Approximated 3D tomographic phase data.
tomopy.prep.stripe
¶
Module for pre-processing tasks.
Functions:
remove_stripe_fw (tomo[, level, wname, …]) |
Remove horizontal stripes from sinogram using the Fourier-Wavelet (FW) based method [C4]. |
remove_stripe_ti (tomo[, nblock, alpha, …]) |
Remove horizontal stripes from sinogram using Titarenko’s approach [C11]. |
remove_stripe_sf (tomo[, size, ncore, nchunk]) |
Normalize raw projection data using a smoothing filter approach. |
-
tomopy.prep.stripe.
remove_stripe_fw
(tomo, level=None, wname=u'db5', sigma=2, pad=True, ncore=None, nchunk=None)[source]¶ Remove horizontal stripes from sinogram using the Fourier-Wavelet (FW) based method [C4].
Parameters: - tomo (ndarray) – 3D tomographic data.
- level (int, optional) – Number of discrete wavelet transform levels.
- wname (str, optional) – Type of the wavelet filter. ‘haar’, ‘db5’, sym5’, etc.
- sigma (float, optional) – Damping parameter in Fourier space.
- pad (bool, optional) – If True, extend the size of the sinogram by padding with zeros.
- 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.stripe.
remove_stripe_ti
(tomo, nblock=0, alpha=1.5, ncore=None, nchunk=None)[source]¶ Remove horizontal stripes from sinogram using Titarenko’s approach [C11].
Parameters: - tomo (ndarray) – 3D tomographic data.
- nblock (int, optional) – Number of blocks.
- alpha (int, optional) – Damping factor.
- 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.stripe.
remove_stripe_sf
(tomo, size=5, ncore=None, nchunk=None)[source]¶ Normalize raw projection data using a smoothing filter approach.
Parameters: - tomo (ndarray) – 3D tomographic data.
- size (int, optional) – Size of the smoothing filter.
- 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.recon.algorithm
¶
Module for reconstruction algorithms.
Functions:
recon (tomo, theta[, center, sinogram_order, …]) |
Reconstruct object from projection data. |
-
tomopy.recon.algorithm.
recon
(tomo, theta, center=None, sinogram_order=False, algorithm=None, init_recon=None, ncore=None, nchunk=None, **kwargs)[source]¶ Reconstruct object from projection data.
Parameters: tomo (ndarray) – 3D tomographic data.
theta (array) – Projection angles in radian.
center (array, optional) – Location of rotation axis.
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 [C2].
- ‘bart’
Block algebraic reconstruction technique.
- ‘fbp’
Filtered back-projection algorithm.
- ‘gridrec’
- ‘mlem’
Maximum-likelihood expectation maximization algorithm [C3].
- ‘osem’
Ordered-subset expectation maximization algorithm [C16].
- ‘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 [C17].
- ‘pml_quad’
Penalized maximum likelihood algorithm with quadratic penalty.
- ‘sirt’
Simultaneous algebraic reconstruction technique.
- ‘tv’
Total Variation reconstruction technique [C7].
- ‘grad’
Gradient descent method with a constant step size
num_gridx, num_gridy (int, optional) – Number of pixels along x- and y-axes in the reconstruction grid.
filter_name (str, 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.
filter_par (list, optional) – Filter parameters as a list.
num_iter (int, optional) – Number of algorithm iterations performed.
num_block (int, optional) – Number of data blocks for intermediate updating the object.
ind_block (array of int, optional) – Order of projections to be used for updating.
reg_par (float, optional) – Regularization parameter for smoothing.
init_recon (ndarray, optional) – Initial guess of the reconstruction.
ncore (int, optional) – Number of cores that will be assigned to jobs.
nchunk (int, optional) – Chunk size for each core.
Returns: ndarray – Reconstructed 3D object.
Warning
Filtering is not implemented for fbp.
Example
>>> import tomopy >>> obj = tomopy.shepp3d() # Generate an object. >>> ang = tomopy.angles(180) # Generate uniformly spaced tilt angles. >>> sim = tomopy.project(obj, ang) # Calculate projections. >>> rec = tomopy.recon(sim, ang, algorithm='art') # Reconstruct object. >>> >>> # Show 64th slice of the reconstructed object. >>> import pylab >>> pylab.imshow(rec[64], cmap='gray') >>> pylab.show()
Example using the ASTRA toolbox for recontruction
For more information, see http://sourceforge.net/p/astra-toolbox/wiki/Home/ and https://github.com/astra-toolbox/astra-toolbox. To install the ASTRA toolbox with conda, use:
conda install -c https://conda.binstar.org/astra-toolbox astra-toolbox
>>> import tomopy >>> obj = tomopy.shepp3d() # Generate an object. >>> ang = tomopy.angles(180) # Generate uniformly spaced tilt angles. >>> sim = tomopy.project(obj, ang) # Calculate projections. >>> >>> # Reconstruct object: >>> rec = tomopy.recon(sim, ang, algorithm=tomopy.astra, >>> options={'method':'SART', 'num_iter':10*180, >>> 'proj_type':'linear', >>> 'extra_options':{'MinConstraint':0}}) >>> >>> # Show 64th slice of the reconstructed object. >>> import pylab >>> pylab.imshow(rec[64], cmap='gray') >>> pylab.show()
tomopy.recon.rotation
¶
Module for functions related to finding axis of rotation.
Functions:
find_center (tomo, theta[, ind, init, tol, …]) |
Find rotation axis location. |
find_center_vo (tomo[, ind, smin, smax, …]) |
Find rotation axis location using Nghia Vo’s method. |
find_center_pc (proj1, proj2[, tol]) |
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. |
write_center (tomo, theta[, dpath, …]) |
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 [C10].
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_vo
(tomo, ind=None, smin=-50, smax=50, srad=6, step=0.5, ratio=0.5, drop=20)[source]¶ Find rotation axis location using Nghia Vo’s method. [C22].
Parameters: - tomo (ndarray) – 3D tomographic data.
- 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.
Returns: float – Rotation axis location.
-
tomopy.recon.rotation.
find_center_pc
(proj1, proj2, tol=0.5)[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
register_translation
function uses cross-correlation in Fourier space, optionally employing an upsampled matrix-multiplication DFT to achieve arbitrary subpixel precision. [C14].Parameters: - proj1 (ndarray) – 2D projection data.
- proj2 (ndarray) – 2D projection data.
- tol (scalar, optional) – Subpixel accuracy
Returns: float – Rotation axis location.
-
tomopy.recon.rotation.
write_center
(tomo, theta, dpath=u'tmp/center', cen_range=None, ind=None, mask=False, ratio=1.0, sinogram_order=False, algorithm=u'gridrec', filter_name=u'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 [C2].
- ‘bart’
Block algebraic reconstruction technique.
- ‘fbp’
Filtered back-projection algorithm.
- ‘gridrec’
- ‘mlem’
Maximum-likelihood expectation maximization algorithm [C3].
- ‘osem’
Ordered-subset expectation maximization algorithm [C16].
- ‘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 [C17].
- ‘pml_quad’
Penalized maximum likelihood algorithm with quadratic penalty.
- ‘sirt’
Simultaneous algebraic reconstruction technique.
- ‘tv’
Total Variation reconstruction technique [C7].
- ‘grad’
Gradient descent method with a constant step size
filter_name (str, 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.
tomopy.sim.project
¶
Module for simulation of x-rays.
Functions:
angles (nang[, ang1, ang2]) |
Return uniformly distributed projection angles in radian. |
project (obj, theta[, center, emission, pad, …]) |
Project x-rays through a given 3D object. |
fan_to_para (tomo, dist, geom) |
Convert fan-beam data to parallel-beam data. |
para_to_fan (tomo, dist, geom) |
Convert parallel-beam data to fan-beam data. |
add_gaussian (tomo[, mean, std]) |
Add Gaussian noise. |
add_poisson (tomo) |
Add Poisson noise. |
add_salt_pepper (tomo[, prob, val]) |
Add salt and pepper noise. |
add_focal_spot_blur (tomo, spotsize) |
Add focal spot blur. |
-
tomopy.sim.project.
angles
(nang, ang1=0.0, ang2=180.0)[source]¶ Return uniformly distributed projection angles in radian.
Parameters: - nang (int, optional) – Number of projections.
- ang1 (float, optional) – First projection angle in degrees.
- ang2 (float, optional) – Last projection angle in degrees.
Returns: array – Projection angles
-
tomopy.sim.project.
project
(obj, theta, center=None, emission=True, pad=True, sinogram_order=False, ncore=None, nchunk=None)[source]¶ Project x-rays through a given 3D object.
Parameters: - obj (ndarray) – Voxelized 3D object.
- theta (array) – Projection angles in radian.
- center (array, optional) – Location of rotation axis.
- emission (bool, optional) – Determines whether output data is emission or transmission type.
- pad (bool, optional) – Determines if the projection image width will be padded or not. If True, then the diagonal length of the object cross-section will be used for the output size of the projection image width.
- sinogram_order (bool, optional) – Determines whether output data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False, theta first axis).
- ncore (int, optional) – Number of cores that will be assigned to jobs.
- nchunk (int, optional) – Chunk size for each core.
Returns: ndarray – 3D tomographic data.
-
tomopy.sim.project.
project2
(objx, objy, theta, center=None, emission=True, pad=True, sinogram_order=False, axis=0, ncore=None, nchunk=None)[source]¶ Project x-rays through a given 3D object.
Parameters: - objx (ndarray) – (x, y) components of vector of a voxelized 3D object.
- theta (array) – Projection angles in radian.
- center (array, optional) – Location of rotation axis.
- emission (bool, optional) – Determines whether output data is emission or transmission type.
- pad (bool, optional) – Determines if the projection image width will be padded or not. If True, then the diagonal length of the object cross-section will be used for the output size of the projection image width.
- sinogram_order (bool, optional) – Determines whether output data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False, theta first axis).
- ncore (int, optional) – Number of cores that will be assigned to jobs.
- nchunk (int, optional) – Chunk size for each core.
Returns: ndarray – 3D tomographic data.
-
tomopy.sim.project.
project3
(objx, objy, objz, theta, center=None, emission=True, pad=True, sinogram_order=False, axis=0, ncore=None, nchunk=None)[source]¶ Project x-rays through a given 3D object.
Parameters: - objx (ndarray) – (x, y) components of vector of a voxelized 3D object.
- theta (array) – Projection angles in radian.
- center (array, optional) – Location of rotation axis.
- emission (bool, optional) – Determines whether output data is emission or transmission type.
- pad (bool, optional) – Determines if the projection image width will be padded or not. If True, then the diagonal length of the object cross-section will be used for the output size of the projection image width.
- sinogram_order (bool, optional) – Determines whether output data is a stack of sinograms (True, y-axis first axis) or a stack of radiographs (False, theta first axis).
- ncore (int, optional) – Number of cores that will be assigned to jobs.
- nchunk (int, optional) – Chunk size for each core.
Returns: ndarray – 3D tomographic data.
-
tomopy.sim.project.
fan_to_para
(tomo, dist, geom)[source]¶ Convert fan-beam data to parallel-beam data.
Warning
Not implemented yet.
Parameters: - tomo (ndarray) – 3D tomographic data.
- dist (float) – Distance from fan-beam vertex to rotation center.
- geom (str) – Fan beam geometry. ‘arc’ or ‘line’.
Returns: ndarray – Transformed 3D tomographic data.
-
tomopy.sim.project.
para_to_fan
(tomo, dist, geom)[source]¶ Convert parallel-beam data to fan-beam data.
Warning
Not implemented yet.
Parameters: - tomo (ndarray) – 3D tomographic data.
- dist (float) – Distance from fan-beam vertex to rotation center.
- geom (str) – Fan beam geometry. ‘arc’ or ‘line’.
Returns: ndarray – Transformed 3D tomographic data.
-
tomopy.sim.project.
add_gaussian
(tomo, mean=0, std=None)[source]¶ Add Gaussian noise.
Parameters: - tomo (ndarray) – 3D tomographic data.
- mean (float, optional) – Mean of the Gaussian distribution.
- std (float, optional) – Standard deviation of the Gaussian distribution.
Returns: ndarray – 3D tomographic data after Gaussian noise added.
-
tomopy.sim.project.
add_poisson
(tomo)[source]¶ Add Poisson noise.
Parameters: tomo (ndarray) – 3D tomographic data. Returns: ndarray – 3D tomographic data after Poisson noise added.
-
tomopy.sim.project.
add_salt_pepper
(tomo, prob=0.01, val=None)[source]¶ Add salt and pepper noise.
Parameters: - tomo (ndarray) – 3D tomographic data.
- prob (float, optional) – Independent probability that each element of a pixel might be corrupted by the salt and pepper type noise.
- val (float, optional) – Value to be assigned to the corrupted pixels.
Returns: ndarray – 3D tomographic data after salt and pepper noise added.
tomopy.sim.propagate
¶
Module for simulation of x-rays.
Functions:
propagate_tie (mu, delta, pixel_size, dist) |
Propagate emitting x-ray wave based on Transport of Intensity. |
-
tomopy.sim.propagate.
calc_intensity
(probe, proj, shift_x=None, shift_y=None, mode=u'near')[source]¶ Calculate far field intensity.
Parameters: - probe (ndarray) – Rectangular x-ray source kernel.
- proj (ndarray) – Object plane intensity image.
- shift_x, shift_y (int, optional) – Shift amount of probe along x and y axes.
- mode (str, optional) – Specify the regime. ‘near’ or ‘far’
Returns: ndarray – Individual raster scanned far field images as 3D array.
-
tomopy.sim.propagate.
propagate_tie
(mu, delta, pixel_size, dist)[source]¶ Propagate emitting x-ray wave based on Transport of Intensity.
Parameters: - mu (ndarray, optional) – 3D tomographic data for attenuation.
- delta (ndarray) – 3D tomographic data for refractive index.
- pixel_size (float) – Detector pixel size in cm.
- dist (float) – Propagation distance of the wavefront in cm.
Returns: ndarray – 3D propagated tomographic intensity.
-
tomopy.sim.propagate.
probe_gauss
(nx, ny, fwhm=None, center=None, max_int=1)[source]¶ Simulate incident x-ray beam (probe) as a square Gaussian kernel.
Parameters: - nx, ny (int) – Grid size along x and y axes.
- fwhm (float, optional) – Effective radius of the source.
- center (array_like, optional) – x and y coordinates of the center of the gaussian function.
- max_int (int) – Maximum x-ray intensity.
Returns: ndarray – 2D source intensity distribution.
Examples¶
This section contains Jupyter Notebooks and Python scripts examples for various tomoPy functions.
To run these examples in a notebooks install Jupyter or run the python scripts from here
Gridrec¶
Here is an example on how to use the gridrec [C5]
reconstruction algorithm with TomoPy
[A1].
You can download the python scritp here
or the Jupyter notebook here
%pylab inline
Install TomoPy then:
import tomopy
Tomographic data input in TomoPy is supported by DXchange.
import dxchange
matplotlib provide plotting of the result in this notebook. Paraview or other tools are available for more sophisticated 3D rendering.
import matplotlib.pyplot as plt
Set the path to the micro-CT data to reconstruct.
fname = '../../tomopy/data/tooth.h5'
Select the sinogram range to reconstruct.
start = 0
end = 2
tooth.h5
data set file format
follows the APS beamline 2-BM and 32-ID
data-exchange file format definition.
Major synchrotron file format examples
are available at DXchange.
proj, flat, dark, theta = dxchange.read_aps_32id(fname, sino=(start, end))
Plot the sinogram:
plt.imshow(proj[:, 0, :], cmap='Greys_r')
plt.show()

If the angular information is not avaialable from the raw data you need to set the data collection angles. In this case theta is set as equally spaced between 0-180 degrees.
if (theta is None):
theta = tomopy.angles(proj.shape[0])
else:
pass
Perform the flat-field correction of raw data:
proj = tomopy.normalize(proj, flat, dark)
Tomopy provides various methods ([C10], [C22], [C14]) to find the rotation center.
rot_center = tomopy.find_center(proj, theta, init=290, ind=0, tol=0.5)
tomopy.rotation:Trying center: [ 290.]
tomopy.rotation:Trying center: [ 304.5]
tomopy.rotation:Trying center: [ 275.5]
tomopy.rotation:Trying center: [ 282.75]
tomopy.rotation:Trying center: [ 297.25]
tomopy.rotation:Trying center: [ 304.5]
tomopy.rotation:Trying center: [ 304.5]
tomopy.rotation:Trying center: [ 293.625]
tomopy.rotation:Trying center: [ 290.]
tomopy.rotation:Trying center: [ 295.4375]
tomopy.rotation:Trying center: [ 291.8125]
tomopy.rotation:Trying center: [ 294.53125]
tomopy.rotation:Trying center: [ 295.4375]
tomopy.rotation:Trying center: [ 294.078125]
Calculate
proj = tomopy.minus_log(proj)
Reconstruction using Gridrec algorithm.Tomopy provides various reconstruction methods including the one part of the ASTRA toolbox.
recon = tomopy.recon(proj, theta, center=rot_center, algorithm='gridrec')
Mask each reconstructed slice with a circle.
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Using the ASTRA toolbox through TomoPy¶
Here is an example on how to use the ASTRA toolbox through its integration with TomoPy, as published in [A2].
%pylab inline
Install the ASTRA toolbox and TomoPy then:
import tomopy
DXchange is installed with tomopy to provide support for tomographic data loading. Various data format from all major synchrotron facilities are supported.
import dxchange
matplotlib provide plotting of the result in this notebook. Paraview or other tools are available for more sophisticated 3D rendering.
import matplotlib.pyplot as plt
Set the path to the micro-CT data to reconstruct.
fname = '../../tomopy/data/tooth.h5'
Select the sinogram range to reconstruct.
start = 0
end = 2
This data set file format follows the APS beamline 2-BM and 32-ID definition. Other file format readers are available at DXchange.
proj, flat, dark, theta = dxchange.read_aps_32id(fname, sino=(start, end))
Plot the sinogram:
plt.imshow(proj[:, 0, :], cmap='Greys_r')
plt.show()

If the angular information is not avaialable from the raw data you need to set the data collection angles. In this case theta is set as equally spaced between 0-180 degrees.
if (theta is None):
theta = tomopy.angles(proj.shape[0])
else:
pass
Perform the flat-field correction of raw data:
proj = tomopy.normalize(proj, flat, dark)
Tomopy provides various methods to find rotation center.
rot_center = tomopy.find_center(proj, theta, init=290, ind=0, tol=0.5)
Calculate
proj = tomopy.minus_log(proj)
Reconstruction with TomoPy¶
Reconstruction can be performed using either TomoPy’s algorithms, or the algorithms of the ASTRA toolbox.
To compare, we first show how to reconstruct an image using TomoPy’s Gridrec algorithm:
recon = tomopy.recon(proj, theta, center=rot_center, algorithm='gridrec')
Mask each reconstructed slice with a circle.
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction with the ASTRA toolbox¶
To reconstruct the image with the ASTRA toolbox instead of TomoPy,
change the algorithm
keyword to tomopy.astra
, and specify the
projection kernel to use (proj_type
) and which ASTRA algorithm to
reconstruct with (method
) in the options
keyword.
More information about the projection kernels and algorithms that are supported by the ASTRA toolbox can be found in the documentation: projection kernels and algorithms. Note that only the 2D (i.e. slice-based) algorithms are supported when reconstructing through TomoPy.
For example, to use a line-based CPU kernel and the FBP method, use:
options = {'proj_type':'linear', 'method':'FBP'}
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.astra, options=options)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

If you have a CUDA-capable NVIDIA GPU, reconstruction times can be greatly reduced by using GPU-based algorithms of the ASTRA toolbox, especially for iterative reconstruction methods.
To use the GPU, change the proj_type
option to 'cuda'
, and use
CUDA-specific algorithms (e.g. 'FBP_CUDA'
for FBP):
options = {'proj_type':'cuda', 'method':'FBP_CUDA'}
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.astra, options=options)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Many algorithms of the ASTRA toolbox support additional options, which
can be found in the
documentation.
These options can be specified using the extra_options
keyword.
For example, to use the GPU-based iterative SIRT method with a nonnegativity constraint, use:
extra_options ={'MinConstraint':0}
options = {'proj_type':'cuda', 'method':'SIRT_CUDA', 'num_iter':200, 'extra_options':extra_options}
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.astra, options=options)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction with UFO¶
UFO is a general-purpose image processing framework developed at the Karlsruhe Institute of Technology and uses OpenCL to execute processing tasks on multiple accelerator devices such as NVIDIA and AMD GPUs, AMD and Intel CPUs as well as Intel Xeon Phi cards.
Here is an example on how to use TomoPy with UFO and its accompanying reconstruction algorithms.
Install
TomoPy,
ufo-core and
ufo-filters. Make sure to
install the Python Numpy interfaces in the python
subdirectory of ufo-core.
DXchange is installed with tomopy to provide support for tomographic data loading. Various data format from all major synchrotron facilities are supported.
import dxchange
matplotlib allows us to plot the result in this notebook.
import matplotlib.pyplot as plt
Set the path to the micro-CT dataset and the sinogram range to reconstruct.
fname = 'tooth.h5'
start, end = (0, 2)
This dataset file format follows the APS beamline 2-BM and 32-ID definition. Other file format readers are available at DXchange.
proj, flat, dark, theta = dxchange.read_aps_32id(fname, sino=(start, end))
Plot the sinogram:
plt.imshow(proj[:, 0, :], cmap='Greys_r')
plt.show()

If the angular information is not available from the raw data you need to set the data collection angles. In this case theta is set as equally spaced between 0-180 degrees.
if (theta is None):
theta = tomopy.angles(proj.shape[0])
else:
pass
Perform the flat-field correction of raw data:
proj = tomopy.normalize(proj, flat, dark)
Tomopy provides various methods to find rotation center.
center = tomopy.find_center(proj, theta, init=290, ind=0, tol=0.5)
Calculate
proj = tomopy.minus_log(proj)
Now, reconstruct using UFO’s filtered backprojection algorithm. Note,
that we must set ncore
to 1 in order to let UFO do the
multi-threading. If left to the default value or set to a value other
than 1 will crash the reconstruction.
recon = tomopy.recon(proj, theta, center=center, algorithm=ufo_fbp, ncore=1)
Mask each reconstructed slice with a circle.
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Vector Reconstruction¶
The vector reconstruction algorithm can be used for instance, to reconstruct the magnetization vector field inside a magnetic sample.
Here is an example on how to use the vector reconstruction algorithm [B2] [A4] with TomoPy[A1].
From a reconstructed 3D object to its projections and back
In order to test the algorithm, the projections of a reconstructed oject can be computed, and from these projections we can come back to the reconstructed model object. Finally we will compare the results of the vector field reconstruction against the initial object.
All datasets used in this tutorial are available in tomoBank.
First, let’s make the necessary imports
import dxchange
import tomopy
import numpy as np
import matplotlib.pyplot as plt
import time
Let’s load the object: the three components of the magnetization vector all throughout the object. The object will be padded in order to have a cubic object. Afterwards it will be downsampled to make faster computations.
obx = dxchange.read_tiff('M4R1_mx.tif').astype('float32')
oby = dxchange.read_tiff('M4R1_my.tif').astype('float32')
obz = dxchange.read_tiff('M4R1_mz.tif').astype('float32')
npad = ((182, 182), (64, 64), (0, 0))
obx = np.pad(obx, npad, mode='constant', constant_values=0)
oby = np.pad(oby, npad, mode='constant', constant_values=0)
obz = np.pad(obz, npad, mode='constant', constant_values=0)
obx = tomopy.downsample(obx, level=2, axis=0)
obx = tomopy.downsample(obx, level=2, axis=1)
obx = tomopy.downsample(obx, level=2, axis=2)
oby = tomopy.downsample(oby, level=2, axis=0)
oby = tomopy.downsample(oby, level=2, axis=1)
oby = tomopy.downsample(oby, level=2, axis=2)
obz = tomopy.downsample(obz, level=2, axis=0)
obz = tomopy.downsample(obz, level=2, axis=1)
obz = tomopy.downsample(obz, level=2, axis=2)
Define the projection angles: 31 angles, from 90 to 270 degrees:
ang = tomopy.angles(31, 90, 270)
And calculate the projections of the object taking rotation axes around the three perpendicular cartesian axes:
prj1 = tomopy.project3(obx, oby, obz, ang, axis=0, pad=False)
prj2 = tomopy.project3(obx, oby, obz, ang, axis=1, pad=False)
prj3 = tomopy.project3(obx, oby, obz, ang, axis=2, pad=False)
The three coordinates of a given projection can be visualized as follows:
fig = plt.figure(figsize=(15, 8))
ax1 = fig.add_subplot(1, 3, 1)
ax1.imshow(obx[52,:,:])
ax2 = fig.add_subplot(1, 3, 2)
ax2.imshow(oby[52,:,:])
ax3 = fig.add_subplot(1, 3, 3)
ax3.imshow(obz[52,:,:])

Finally we will reconstruct the vector field components, taking as input the projections that we have calculated thanks to the first 3D initial object. The number of iterations can be adapted to have a faster but more imprecise reconstruction, or to have a more precise reconstruction.
rec1, rec2, rec3 = tomopy.vector3(prj1, prj2, prj3, ang, ang, ang, axis1=0, axis2=1, axis3=2, num_iter=100)
dxchange.write_tiff(rec1)
dxchange.write_tiff(rec2)
dxchange.write_tiff(rec3)
Comparison of results against input object
In this section, we compare the results of the vector field components obtained thanks to the tomopy reconstruction, against the vector field components of the object given as input:
Comparison of the first magnetization vector component against the input data object (for a given slice).
fig = plt.figure(figsize=(9, 7))
ax1 = fig.add_subplot(1, 2, 1)
ax1.imshow(obx[52,:,:])
ax2 = fig.add_subplot(1, 2, 2)
ax2.imshow(rec1[52,:,:])

Comparison of the second magnetization vector component against the input data object (for a given slice):
fig = plt.figure(figsize=(9, 7))
ax1 = fig.add_subplot(1, 2, 1)
ax1.imshow(oby[52,:,:])
ax2 = fig.add_subplot(1, 2, 2)
ax2.imshow(rec2[52,:,:])

Comparison of the third magnetization vector component against the input data object (for a given slice):
fig = plt.figure(figsize=(9, 7))
ax1 = fig.add_subplot(1, 2, 1)
ax1.imshow(obz[52,:,:])
ax2 = fig.add_subplot(1, 2, 2)
ax2.imshow(rec3[52,:,:])

Other examples
Three jupyter notebooks with examples as well as with some mathematical concepts related to the vector reconstruction, can be found in the tomopy/doc/demo folder:
Examples using vector3: input data projections from 3 orthogonal tilt angles:
- vectorrec_1.ipynb
- vectorrec_disk.ipynb
Example using vector2: input data projections from 2 orthogonal tilt angles:
- vector_heterostructure.ipynb
The Vector Reconstruction examples html slides can be build by applying (from the doc/demo folder) the following commands:
jupyter-nbconvert --to slides --post serve vectorrec_1.ipynb
jupyter-nbconvert --to slides --post serve vectorrec_disk.ipynb
jupyter-nbconvert --to slides --post serve vector_heterostructure.ipynb
LPrec¶
Here is an example on how to use the log-polar based method (https://github.com/math-vrn/lprec) for reconstruction in Tomopy
%pylab inline
Install lprec from github, then
import tomopy
DXchange is installed with tomopy to provide support for tomographic data loading. Various data format from all major synchrotron facilities are supported.
import dxchange
matplotlib provide plotting of the result in this notebook. Paraview or other tools are available for more sophisticated 3D rendering.
import matplotlib.pyplot as plt
Set the path to the micro-CT data to reconstruct.
fname = '../../tomopy/data/tooth.h5'
Select the sinogram range to reconstruct.
start = 0
end = 2
This data set file format follows the APS beamline 2-BM and 32-ID definition. Other file format readers are available at DXchange.
proj, flat, dark = dxchange.read_aps_32id(fname, sino=(start, end))
Plot the sinogram:
plt.imshow(proj[:, 0, :], cmap='Greys_r')
plt.show()

If the angular information is not avaialable from the raw data you need to set the data collection angles. In this case theta is set as equally spaced between 0-180 degrees.
theta = tomopy.angles(proj.shape[0])
Perform the flat-field correction of raw data:
proj = tomopy.normalize(proj, flat, dark)
Select the rotation center manually
rot_center = 296
Calculate
proj = tomopy.minus_log(proj)
Reconstruction using FBP method with the log-polar coordinates
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.lprec, lpmethod='fbp', filter_name='parzen')
Mask each reconstructed slice with a circle.
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction using the gradient descent method with the log-polar coordinates
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.lprec, lpmethod='grad', ncore=1, num_iter=64, reg_par=-1)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction using the conjugate gradient method with the log-polar coordinates
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.lprec, lpmethod='cg', ncore=1, num_iter=16, reg_par=-1)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction using the TV method with the log-polar coordinates
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.lprec, lpmethod='tv', ncore=1, num_iter=256, reg_par=1e-3)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Reconstruction using the MLEM method with the log-polar coordinates
recon = tomopy.recon(proj, theta, center=rot_center, algorithm=tomopy.lprec, lpmethod='em', ncore=1, num_iter=64, reg_par=0.05)
recon = tomopy.circ_mask(recon, axis=0, ratio=0.95)
plt.imshow(recon[0, :,:], cmap='Greys_r')
plt.show()

Frequently asked questions¶
Here’s a list of questions.
Questions
- How can I report bugs?
- Are there any video tutorials?
- Can I run this on a HPC cluster?
- Are there any segmentation routines?
- Are there any tools for aligning projections?
- What is ASTRA toolbox?
- Why TomoPy and ASTRA were integrated?
- Which platforms are supported?
- Does ASTRA support all GPUs?
- What is UFO?
How can I report bugs?¶
The easiest way to report bugs or get help is to open an issue on GitHub. Simply go to the project GitHub page, click on Issues in the right menu tab and submit your report or question.
Are there any video tutorials?¶
We currently do not have specific plans in this direction, but we agree that it would be very helpful.
Can I run this on a HPC cluster?¶
In their default installation packages, TomoPy and the ASTRA toolbox are limited to running on a multicore single machine. The ASTRA toolbox, and TomoPy through the presented ASTRA integration, are able to use multiple GPUs that are installed in a single machine. Both toolboxes can be run on a HPC cluster through parallelization using MPI, but since installation and running on a HPC cluster is often cluster specific, the default installation packages do not include these capabilities.
As such, the integrated packages that is presented in the manuscript currently does not support running on a HPC cluster. Note that the ASTRA toolbox provides a separate MPI enabled package for use on a HPC cluster. We refer to [C23] for more details about TomoPy’s planned HPC implementation. It is a MapReduce type MPI implementation layer, which was succesfully used on many clusters, i.e. Stampede, Cori, Mira. There are plans to allow user access to TomoPy on a HPC cluster (e.g. through a client or webportal), but these projects will take some time before they are being matured for user’s use.
Are there any segmentation routines?¶
Some data processing operations can be applied after reconstruction. Examples of these type of operations are image based ring removal methods, and gaussian filtering or median filtering the reconstructed image. Typicaly, these methods are called “postprocessing algorithms, since they occur after the reconstruction.
The package does not include segmentation algorithms, since we are currently focused on tomography, while we feel that segmentation are more part of the application specific data analysis that occurs after tomographic processing. An important exception is when segmentation steps are used as part of the tomographic reconstruction algorithm, such as in the DART algorithm.
Are there any tools for aligning projections?¶
Yes we have. Please check the Examples section for details.
What is ASTRA toolbox?¶
The ASTRA toolbox provides highly efficient tomographic reconstruction methods by implementing them on graphic processing units (GPUs). It includes advanced iterative methods and allows for very flexible scanning geometries. The ASTRA toolbox also includes building blocks which can be used to develop new reconstruction methods, allowing for easy and efficient implementation and modification of advanced reconstruction methods. However, the toolbox is only focused on reconstruction, and does not include pre-processing or post-processing methods that are typically required for correctly processing synchrotron data. Furthermore, no routines to read data from disk are provided by the toolbox.
Why TomoPy and ASTRA were integrated?¶
The TomoPy toolbox is specifically designed to be easy to use and deploy at a synchrotron facility beamline. It supports reading many common synchrotron data formats from disk [C13], and includes several other processing algorithms commonly used for synchrotron data. TomoPy also includes several reconstruction algorithms, which can be run on multi-core workstations and large-scale computing facilities. The algorithms in TomoPy are all CPU-based, however, which can make them prohibitively slow in the case of iterative methods, which are often required for advanced tomographic experiments.
By integrating the ASTRA toolbox in the TomoPy framework, the optimized GPU-based reconstruction methods become easily available for synchrotron beamline users, and users of the ASTRA toolbox can more easily read data and use TomoPy’s other functionality for data filtering and cleaning.
Which platforms are supported?¶
TomoPy supports Linux and Mac OS X, and the ASTRA toolbox supports Linux and Windows. As such, the combined package currently supports only Linux, but we are working on supporting more operating systems.
Does ASTRA support all GPUs?¶
The GPU algorithms are all implemented used nVidia CUDA. As a result, only nVidia CUDA enabled video cards can be used to run them.
What is UFO?¶
UFO is a general purpose image processing framework, optimized for heterogeneous compute systems and streams of data. Arbitrary data processing tasks are plugged together to form larger processing pipelines. These pipelines are then mapped to the hardware resources available at run-time, i.e. both multiple GPUs and CPUs.
One specific use case that has been integrated into the TomoPy is fast reconstruction using the filtered backprojection and direct Fourier inversion methods although others for pre- and post-processing might be added in the future.
Credits¶
We kindly request that you cite the following article(s) [A1] if you use TomoPy (and also cite [A2] if you use ASTRA or [A3] if you use UFO). For vector reconstructions please additionally cite [A4].
[A1] | Gürsoy D, De Carlo F, Xiao X, and Jacobsen C. Tomopy: a framework for the analysis of synchrotron tomographic data. Journal of Synchrotron Radiation, 21(5):1188–1193, 2014. |
[A2] | Pelt D, Gürsoy D, Palenstijn WJ, Sijbers J, De Carlo F, and Batenburg KJ. Integration of tomopy and the astra toolbox for advanced processing and reconstruction of tomographic synchrotron data. Journal of Synchrotron Radiation, 23(3):842–849, 2016. |
[A3] | Vogelgesang M, Chilingaryan S, Rolo T dos Santos, and Kopmann A. Ufo: a scalable gpu-based image processing framework for on-line monitoring. In Proceedings of The 14th IEEE Conference on High Performance Computing and Communication & The 9th IEEE International Conference on Embedded Software and Systems (HPCC-ICESS), 824–829. 6 2012. |
[A4] | Hierro-Rodriguez A, Gürsoy D, Phatak C, Quiros C, Sorrentino A, Alvarez-Prado LM, Velez M, Martin JI, Alameda JM, Pereiro E, and Ferrer S. 3d reconstruction of magnetization from dichroic soft x-ray transmission tomography. Journal of Synchrotron Radiation, 2018. |
Applications¶
[B1] | Patterson BM, Cordes NL, Henderson K, Williams JJ, Stannard T, Singh SS, Ovejero AR, Xiao X, Robinson M, and Chawla N. In situ x-ray synchrotron tomographic imaging during the compression of hyper-elastic polymeric materials. Journal of Materials Science, 51(1):171–187, 2016. |
[B2] | Phatak C and Gürsoy D. Iterative reconstruction of magnetic induction using lorentz transmission electron tomography. Ultramicroscopy, 150:54–64, 2015. |
[B3] | Gürsoy D, Biçer T, Lanzirotti A, Newville MG, and De Carlo F. Hyperspectral image reconstruction for x-ray fluorescence tomography. Optics Express, 23(7):9014–9023, 2015. |
[B4] | Gürsoy D, Biçer T, Almer JD, Kettimuthu R, Stock SR, and De Carlo F. Maximum a posteriori estimation of crystallographic phases in x-ray diffraction tomography. Philosophical Transactions A, 2015. |
[B5] | Duke DJ, Swantek AB, Sovis N, Tilocco FZ, Powell CF, AL Kastengren, Gürsoy D, and Biçer T. Time-resolved x-ray tomography of gasoline direct injection sprays. SAE International Journal of Engines, 2015. |
[B6] | Kamke FA, McKinley PE, Ching DJ, Zauner M, and Xiao X. Micro x-ray computed tomography of adhesive bonds in wood. Wood and Fiber Science, 2016. |
[B7] | Birkbak ME, Leemreize H, Frohlich S, Stock SR, and Birkedal H. Diffraction scattering computed tomography: a window into the structures of complex nanomaterials. Nanoscale, 2015. |
[B8] | Miller SM, Xiao X, and Faber KT. Freeze-cast alumina pore networks: effects of freezing conditions and dispersion medium. Journal of the European Ceramic Society, 35(13):3595––3605, 2015. |
[B9] | Roncal WG, Dyer EL, Gürsoy D, Kording K, and Kasthuri N. From sample to knowledge: towards an integrated approach for neuroscience discovery. arXiv, 2016. |
References¶
[C1] | Bergamaschi A, Medjoubi K, Messaoudi C, Marco S, and Somogyi A. Mmx-i: data-processing software for multimodal x-ray imaging and tomography. Journal od Synchrotron Radiation, 23:783–794, 2016. |
[C2] | Kak AC and Slaney M. Principles of computerized tomographic imaging. Volume 33. SIAM, 1988. |
[C3] | Dempster AP, Laird NM, and Rubin DB. Maximum likelihood from incomplete data via the em algorithm. journal of the Royal Statistical Society, Series B, 39(1):1–38, 1977. |
[C4] | Münch B, Trtik P, Marone F, and Stampanoni M. Stripe and ring artifact removal with combined wavelet–fourier filtering. Optics Express, 17(10):8567–8591, 2009. |
[C5] | Dowd BA, Campbell GH, Marr RB, Nagarkar VV, Tipnis SV, Axe L, and Siddons DP. Developments in synchrotron x-ray computed microtomography at the national synchrotron light source. In Proc. SPIE, volume 3772, 224–236. 1999. |
[C6] | Francesco De Carlo, Doga Gursoy, Daniel Jackson Ching, Kees Joost Batenburg, Wolfgang Ludwig, Lucia Mancini, Federica Marone, Rajmund Mokso, Daniel M. Pelt, Jan Sijbers, and Mark Rivers. Tomobank: a tomographic data repository for computational x-ray science. Measurement Science and Technology, 2017. URL: https://doi.org/10.1088/1361-6501/aa9c19. |
[C7] | Antonin Chambolle and Thomas Pock. A first-order primal-dual algorithm for convex problems with applications to imaging. Journal of mathematical imaging and vision, 40(1):120–145, 2011. |
[C8] | Gursoy D, Hong YP, He K, Hujsak K, Yoo S, Chen S, Li Y, Ge M, Miller LM, Chu YS, De Andrade V, He K, Cossairt O, Katsaggelos AK, and Jacobsen C. Rapid alignment of nanotomography data using joint iterative reconstruction and reprojection. Scientific Reports, 2017. |
[C9] | Paganin D, Mayo SC, Gureyev TE, Miller PR, and Wilkins SW. Simultaneous phase and amplitude extraction from a single defocused image of a homogeneous object. Journal of Microscopy, 206(1):33–40, 2002. |
[C10] | Tilman D, Felix B, and Andreas S. Automated determination of the center of rotation in tomography data. Journal of the Optical Society of America A, 23(5):1048–1057, 2006. |
[C11] | Miqueles EX, Rinkel J, O’Dowd F, and Bermúdez JSV. Generalized titarenko’s algorithm for ring artefacts reduction. Journal of Synchrotron Radiation, 21(6):1333–1346, 2014. |
[C12] | Brun F, Pacile S, Accardo A, Kourousias G, Dreossi D, Mancini L, Tromba G, and Pugliese R. Enhanced and flexible software tools for x-ray computed tomography at the italian synchrotron radiation facility elettra. Fundamenta Informaticae, 141(2-3):233–243, 2015. |
[C13] | De Carlo F, Gursoy D, Marone F, Rivers M, Parkinson YD, Khan F, Schwarz N, Vine DJ, Vogt S, Gleber SC, Narayanan S, Newville M, Lanzirotti T, Sun Y, Hong YP, and Jacobsen C. Scientific data exchange: a schema for hdf5-based storage of raw and analyzed data. Journal of Synchrotron Radiation, 21(6):1224–1230, 2014. |
[C14] | Manuel GS, Thurman ST, and Fienup JR. Efficient subpixel image registration algorithms. Optics Letters, 33(2):156–158, 2008. |
[C15] | Toby HB, Gürsoy D, De Carlo F, Schwarz N, Sharma H, and Jacobsen CJ. Practices and standards for data and processing at the APS. Synchrotron Radiation News, 28(2):15–21, 2015. |
[C16] | Hudson HM and Larkin RS. Accelerated image reconstruction using ordered subsets of projection data. Medical Imaging, IEEE Transactions on, 13(4):601–609, 1994. |
[C17] | Chang J-H, Anderson JMM, and Votaw JT. Regularized image reconstruction algorithms for positron emission tomography. Medical Imaging, IEEE Transactions on, 23(9):1165–1175, 2004. |
[C18] | Mertens JCE and Chawla JJWN. A method for zinger artifact reduction in high-energy x-ray computed tomography. Nuclear Instruments and Methods in Physics Research Section A, 800:82–92, 2015. |
[C19] | Vogelgesang M, Rota L, Ardila Perez Luis E, Caselle M, Chilingaryan S, and Kopmann A. High-throughput data acquisition and processing for real-time x-ray imaging. In Proc. SPIE, volume 9967, 996715–996715–9. 2016. doi:10.1117/12.2237611. |
[C20] | De Jonge MD, Ryan CG, and Jacobsen C. X-ray nanoprobes and diffraction-limited storage rings: opportunities and challenges of fluorescence tomography of biological specimens. Journal of Synchrotron Radiation, 21(5):1031–1047, 2014. |
[C21] | Rivers ML. Tomorecon: high-speed tomography reconstruction on workstations using multi-threading. In Proc. SPIE, volume 8506, 85060U–85060U–13. 2012. |
[C22] | Vo N, Drakopoulos M, Atwood RC, and Reinhard C. Reliable method for calculating the center of rotation in parallel-beam tomography. Optics Express, 22(16):19078–19086, 2014. |
[C23] | Biçer T, Gürsoy D, Kettimuthu R, De Carlo F, Agrawal G, and Foster IT. Rapid tomographic image reconstruction via large-scale parallelization. In Lecture Notes in Computer Science, volume 9233, 289–302. 2015. |
[C24] | Bhimji W, Bard D, Roumanus M, Paul D, Ovsyannikov A, Friesen B, Bryson M, Correa, Lockwood GK, Tsulaia V, Byna S, Farrell S, Gürsoy D, Daley C, Beckner V, Van Straalen B, Wright NJ, Antypas K, and Prabhat M. Accelerating science with the nersc burts buffer early user program. In Cray User Group Conference. 2016. |
[C25] | Xu W and Feng D. Studying performance of a penalized maximum likelihood method for pet reconstruction on nvidia gpu and intel xeon phi coprocessor. Proc. 4th Intl Conf. Image Formation in X-ray Computed Tomography, pages 191–194, 2016. |