Modules

flexDATA includes three modules: geometry, data and display.

A wide range of acquisition geometries can be defined using the geometry module. For I/O operations and some utilities for out-of-memory arrays are implemented in the data module. Basic display utilities can be found in display module.

flexdata.geometry module

This module contains acquisition geometry classes: circular, linear, helical. The circular class corresponds to the simplest case of circular orbit cone-beam CT with minimal number of parameters. Additional parameters can be use to define a non-conventional geometry. For instance: offsets and rotations (‘det_roll’, ‘det_tan’, ‘det_ort’, etc.), axis tilts (‘axs_roll’, ‘axs_pitch’), volume transformations (‘vol_tra’, ‘vol_rot’), recostruction resolution and anisotropic sampling (‘img_pixel’, ‘det_sample’, ‘vol_sample’).

flexdata.geometry.astra_projection_geom(geom, data_shape, index=None)[source]

Initialize ASTRA projection geometry.

Args:
geom : geometry class data_shape: [detector_count_z, theta_count, detector_count_x] index : if provided - sequence of the rotation angles
flexdata.geometry.astra_volume_geom(geom, vol_shape, slice_first=None, slice_last=None)[source]

Initialize ASTRA volume geometry.

class flexdata.geometry.basic(src2obj=None, det2obj=None, det_pixel=None, img_pixel=None, unit='mm')[source]

Bases: object

Base geometry class. Needs only SDD, ODD and pixel size to initialize. Use ‘parameters’ to store parameters and ‘description’ to store relevant metadata.

astra_projection_geom(data_shape, index=None)[source]

Get ASTRA projection geometry.

Args:
data_shape: [detector_count_z, theta_count, detector_count_x] index : if provided - sequence of the rotation angles
Returns:
geometry : ASTRA cone-beam geometry.
astra_volume_geom(vol_shape, slice_first=None, slice_last=None)[source]

Initialize ASTRA volume geometry.

Args:
vol_shape : volume array shape slice_first: first slice of an ROI to update slice_last : last slice of an ROI to update
copy()[source]

Copy me.

det2obj

Detector-to-object distance.

det_sample

Pixel shape.

detector_bounds(proj_shape)[source]

Get the boundaries of the detector at the start of the scan in length units.

detector_centre()[source]

Get the centre coordinate of the first position of the detector.

detector_size(proj_shape)[source]

Get the size of detector in length units.

from_astra_cone_vec(vectors)[source]

Take vectors from ASTRA ‘cone_vec’ geometry. This will override any other parameters.

from_dictionary(dictionary)[source]

Use dictionary records to initialize this geometry.

from_matrix(R, T)[source]

Rotates and translates the reconstruction volume.

Args:
R (3x3 array): rotation matrix T (1x3 array): translation vector
get_detector_orbit(proj_count=None, index=None)[source]

Get the detector orbit. In the base class it is a circular orbit.

Args:
proj_count: number of projections index : index of the projection subset
Returns:
src_pos : array of the source positions. Can be generated by circular_orbit(…), for instance.
get_source_orbit(proj_count=None, index=None)[source]

Get the source orbit. In the base class it is a circular orbit.

Args:
proj_count: number of projections index : index of the projection subset
Returns:
src_pos : array of the source positions. Can be generated by circular_orbit(…), for instance.
get_vectors(proj_count, index=None)[source]

Get source, detector and detector orientation vectors.

Args:
angle_count : number of rotation angles index : index of angles that should be used
magnification

Magnification.

pixel

Pixel size (mm).

src2det

Source-to-detector distance.

src2obj

Source-to-object distance.

to_dictionary()[source]

Return a dictionary describing this geometry.

vol_sample

Voxel shape.

volume_bounds(vol_shape)[source]

Return volume bounds:

volume_size(vol_shape)[source]

Return volume size in length units.

volume_xyz(shape, offset=[0.0, 0.0, 0.0])[source]

Coordinate grid in units of the geometry.

Args:
shape : volume shape offset: offset in length units
voxel

Voxel size (mm).

class flexdata.geometry.circular(src2obj=None, det2obj=None, det_pixel=None, img_pixel=None, ang_range=(0, 360), unit='mm')[source]

Bases: flexdata.geometry.basic

Circular orbit geometry class. Includes additional parameters such as detector and source shifts and rotations.

get_detector_orbit(proj_count=None, index=None)[source]

Get the detector orbit.

Args:
proj_count: number of projections index : index of the projection subset
Returns:
det_pos : array of the detector positions. det_tan : array of detector tangential directional vector det_rad : array of detector radial directional vector det_orth: array of detector orthogonal directional vector
get_source_orbit(proj_count=None, index=None)[source]

Get the source orbit.

Args:
proj_count: number of projections index : index of the projection subset
Returns:
src_pos : array of the source positions.
get_thetas(proj_count=None, index=None)[source]

Get rotation angles. Either returns an equidistant array or self.thetas if that is defined.

Args:
proj_count : number of angles in equidistant case. Not used if ‘thetas’ are defined explicitly in parameters.
Returns:
thetas : angle array in degrees
flexdata.geometry.circular_orbit(radius, thetas, roll=0, pitch=0, yaw=0, origin=[0, 0, 0], tan_shift=0, index=None)[source]

Generate a circular orbit vector.

Args:
radius : orbit radius angle_count: number of rotation angles roll, pitch: define orientation of the rotation axis yaw : initial angular position origin : xyz vector of the orbit centre tan_shift : tangential shift from the default position (scalar or array) index : index of the subset of total rotation angles
Returns:
position : position vector tangent : tangent direction radius : radal direction orthogonal : orthogonal direction
flexdata.geometry.detector_bounds(geom, proj_shape)[source]

Get the boundaries of the detector in length units.

flexdata.geometry.detector_size(geom, proj_shape)[source]

Get the size of detector in length units.

flexdata.geometry.get_vectors(geom, angle_count, index=None)[source]

Get source, detector and detector orientation vectors.

Args:
geom : geometry class angle_count : number of rotation angles index : index of angles that should be used
class flexdata.geometry.helical(src2obj=None, det2obj=None, det_pixel=None, img_pixel=None, axis_range=(0, 100), ang_range=(0, 720), unit='mm')[source]

Bases: flexdata.geometry.circular

Helical orbit geometry class. Similar to the ‘circular’ class with additional parameter of helix

get_detector_orbit(proj_count=None, index=None)[source]

Get the detector orbit. In the base class it is a circular orbit.

get_source_orbit(proj_count=None, index=None)[source]

Get the source orbit. In the base class it is a circular orbit.

class flexdata.geometry.linear(src2obj=None, det2obj=None, det_pixel=None, img_pixel=None, src_hrz_rng=(0, 1), src_vrt_rng=(0, 1), det_hrz_rng=(1, 0), det_vrt_rng=(1, 0), unit='mm')[source]

Bases: flexdata.geometry.basic

A simple linear orbit geometry class.

get_detector_orbit(proj_count=None, index=None)[source]

Get the detector orbit. In the base class it is a circular orbit.

get_source_orbit(proj_count=None, index=None)[source]

Get the source orbit. In the base class it is a circular orbit.

flexdata.geometry.linear_orbit(hrz_rng, rad_rng, vrt_rng, proj_count, index=None)[source]

Generate a linear orbit vector.

Args:
hrz_rng: horizontal range of motion rad_rng: radial range of motion vrt_rng: vertical range of motion
Returns:
position : position vector tangent : tangent direction radius : radal direction orthogonal : orthogonal direction
flexdata.geometry.tiles_shape(shape, geometry_list)[source]

Compute the size of the stiched dataset.

Args:
shape: shape of a single projection stack. geometry_list: list of geometries.
flexdata.geometry.volume_bounds(geom, proj_shape)[source]

A very simplified version of volume bounds…

flexdata.geometry.volume_shape(geom, proj_shape)[source]

Based on physical volume bnounds compute shape in pixels:

flexdata.data module

This module contains some input / output routines for stacks of images and parsers for translating scanner log files into geometry definitions.

Most of the basic image formats are supported through imageio module. Raw binaries and matlab binary files can be loaded.

Utility functions to hande big arrays of data. All routines support memmap arrays. However, some operations will need enough memory for at least one copy of the data for intermediate results. This can be improved through better use of memmaps.

flexdata.data.add_dim(array_1, array_2, dim=None)[source]

Add two arrays with arbitrary dimensions. We assume that one or two dimensions match.

flexdata.data.anyslice(array, index, dim)[source]

Slice an array along an arbitrary dimension.

flexdata.data.autocorrelation(array, axes=(0, 1, 2))[source]

Compute autocorrelation.

flexdata.data.bin(array, dim=None, geometry=None)[source]

Simple binning of the data along the chosen direction.

flexdata.data.cast2shape(array, shape)[source]

Make the array to conform with the given shape.

flexdata.data.cast2type(array, dtype, bounds=None)[source]

Cast from float to int or float to float rescaling values if needed.

flexdata.data.convolve_filter(array, filtr)[source]

Apply a filter defined in Fourier space (a CTF) via convolution.

Args:
array : data array (implicit) filtr : a filter defined in Fourier space (1D - 3D)
flexdata.data.convolve_kernel(array, kernel)[source]

Apply a kernel defined in real space (center in the center of the array) via convolution.

Args:
array (ndarray) : data array (implicit) kernel(ndarray) : real space kernel (1D - 3D)
flexdata.data.crop(array, dim, width, geometry=None)[source]

Crop an array along the given dimension. Provide geometry if cropping the projection data, it will update the detector center.

flexdata.data.deconvolve_filter(array, filtr, epsilon)[source]

Inverse convolution with Tikhonov regularization.

Args:
array (ndarray) : data array (implicit) filtr (ndarray) : Fourier space filter (1D - 3D) epsilon : regularization parameter axes : list of axes to apply deconvolution to.
flexdata.data.deconvolve_kernel(array, kernel, epsilon, axes=(0, 2))[source]

Inverse convolution with Tikhonov regularization.

Args:
array (ndarray) : data array (implicit) filtr (ndarray) : Fourier space filter (1D - 3D) epsilon : regularization parameter axes : list of axes to apply deconvolution to.
flexdata.data.divergence(array, axes=[0, 1, 2])[source]

Compute the divergence of an array.

Args:
axes : list of axes where the divergence is applied.
Returns:
ndarray: divergence of the input array
flexdata.data.file_to_dictionary(file_path, separator=':', translation=None)[source]

Read a text file and return a dictionary with records.

Args:
file_path (str): file to read separator (str): separator between the keys and values translation (dict): dictionary for translating initial keys to a new naming
flexdata.data.flipdim(array, transpose=[1, 0, 2], updown=True)[source]

Convert a given numpy array (sorted: index, hor, vert) to ASTRA-compatible projections stack

flexdata.data.free_disk(path)[source]

Return amount of free memory on disk in GB.

Args:
percent (bool): percentage of the total or in GB.
flexdata.data.free_memory(percent=False)[source]

Return amount of free RAM memory in GB.

Args:
percent (bool): percentage of the total or in GB.
flexdata.data.get_files_sorted(path, name)[source]

Sort file entries using the natural (human) sorting

flexdata.data.get_folders_sorted(path)[source]

Get all paths from a path with a star (using glob)

flexdata.data.gradient(array, axes=[0, 1, 2])[source]

Compute the gradient of an array.

Args:
axes : list of axes to apply gradient to.
Returns:
ndarray: shape = (3, k, l, m) where k,l,m - dimensions of the input array.
class flexdata.data.logger[source]

Bases: object

A class for logging and printing messages.

static error(message)[source]

Raise an error.

file = ''
static print(message)[source]

Simply prints and saves a message.

static title(message)[source]

Print something important.

static warning(message)[source]

Raise a warning.

flexdata.data.medipix2astra(array)[source]

Convert a given numpy array (sorted: index, hor, vert) to ASTRA-compatible projections stack

class flexdata.data.memmap[source]

Bases: numpy.core.memmap.memmap

Standard memmaps don’t seem to reliably delete files that are created on disk. This fixes it…

delete()[source]
flexdata.data.mult_dim(array_1, array_2, dim=None)[source]

Multiply a 3D array by a 1D or a 2D vector along one of the dimensions.

flexdata.data.pad(array, dim, width, mode='edge', geometry=None)[source]

Pad an array along a given dimension. numpy.pad seems to be very memory hungry! Don’t use it for large arrays.

Args:
array: input array dim : dim to apply the ramp width: width of the ramp mode :’linear’ - creates linear decay of intensity; ‘edge’ - smears data in a costant manner; ‘zero’ - sets values to zeroes. geometry: geometry record to update (updates detector offset).
flexdata.data.ramp(array, dim, width, mode='linear')[source]

Create ramps at the ends of the array (without changing its size).

Args:
array: input array dim : dim to apply the ramp width: width of the ramp mode :’linear’ - creates linear decay of intensity; ‘edge’ - smears data in a costant manner; ‘zero’ - sets values to zeroes.
flexdata.data.raw2astra(array)[source]

Convert a given numpy array (sorted: index, hor, vert) to ASTRA-compatible projections stack

flexdata.data.read_flexray(path, sample=1, skip=1, memmap=None, proj_number=None)[source]

Read projecition data for the FLex-Ray scaner. Read, dark-, flat-field images and scan parameters.

Args:
path (str): path to flexray data. skip (int): read every ## image sample (int): keep every ## x ## pixel memmap (str): output a memmap array using the given path proj_number (int): force projection number (treat lesser numbers as missing)
Returns:
proj (numpy.array): projections stack flat (numpy.array): reference flat field images dark (numpy.array): dark field images geom (geometry) : description of the geometry, physical settings and comments
flexdata.data.read_flexraylog(path, sample=1)[source]

Read the log file of FLexRay scanner and return dictionaries with parameters of the scan.

Args:
path (str): path to the files location sample (int): subsampling of the input data
Returns:
geometry : circular geometry class
flexdata.data.read_flexraymeta(path, sample=1)[source]

Read the metafile produced by Flexray scripting.

Args:
path (str): path to the files location sample (int): subsampling of the input data
Returns:
geometry : circular geometry class
flexdata.data.read_geometry(path, sample=1)[source]

Read a native meta file.

Args:
path (str): path to the file location sample (int): subsampling of the input data
Returns:
geometry : circular geometry class
flexdata.data.read_image(file, sample=1, shape=None, format=None, dtype=None)[source]

Read a single image. Use sampling and roi parameters to reduce the array size. Use shape, format and dtype to force file reading settings.

flexdata.data.read_stack(path, name, skip=1, sample=1, shape=None, dtype=None, format=None, transpose=[1, 0, 2], updown=True, memmap=None, success=None)[source]

Read stack of files and return a numpy array.

Args:
path (str): path to the files location name (str): common part of the files name skip (int): read every so many files sample (int): sampling factor in x/y direction shape (array): shape of the files. Use it when the format is ‘raw’. dtype (str or numpy.dtype): data type of the files format (str): file format (‘tif’, ‘raw’, etc) flipdim (bool): apply dimension switch for ASTRA compatibility memmap (str): if provided, return a disk mapped array to save RAM success(array): map of the files that could be loaded (equals 0 in case of a read failure)
Returns:
numpy.array : 3D array with the first dimension representing the image index
flexdata.data.read_toml(file_path)[source]

Read a toml file.

Args:
file_path (str): read file form that location
flexdata.data.rewrite_memmap(old_array, new_array)[source]

Reshaping memmaps is tough. We will recreate one instead hoping that this will not overflow our RAM… This is a dirty qick fix! Try to use resize instead!

flexdata.data.shape_alike(array_1, array_2)[source]
Make sure two arrays have the same shape by padding either array_1 or array_2:
Returns: array1, array2 - reshaped.
flexdata.data.stack_shape(path, name, skip=1, sample=1, shape=None, dtype=None, format=None)[source]

Determine the shape of stack on disk.

Args:
path (str): path to the files location name (str): common part of the files name skip (int): read every so many files sample (int): sampling factor in x/y direction shape (array): shape of the files. Use it when the format is ‘raw’. dtype (str or numpy.dtype): data type of the files format (str): file format (‘tif’, ‘raw’, etc)
flexdata.data.write_astra(filename, data_shape, geom)[source]

Write an astra-readable projection geometry vector.

flexdata.data.write_image(filename, image, compress=0)[source]

Write a single image. Use compression if needed (0-9).

flexdata.data.write_stack(path, name, data, dim=1, skip=1, dtype=None, zip=False, format='tiff', updown=False)[source]

Write an image stack.

Args:
path (str): destination path name (str): first part of the files name data (numpy.array): data to write dim (int): dimension along which array is separated into images skip (int): how many images to skip in between dtype (type): forse this data type compress (str): use None, ‘zip’ or ‘jp2’. format (str): file extension (‘raw’, ‘tiff’, ‘jp2’, etc)
flexdata.data.write_toml(filename, record)[source]

Write a toml file.

Args:
filename (str): location to write the file to record (dict, geometry): geomety class record or an arbitrary dictionary

flexdata.display module

This module contains a few simple routines for displaying data: * 2D displays like: slice, projection, etc. * Interactive slicer: pyqt_graph * Other displays: mesh, color_project

flexdata.display.color_project(array, dim=1, sample=2, bounds=[0.01, 0.1], title=None, cmap='nipy_spectral', file=None)[source]

Create a pseudo color projection of a 3D volume.

flexdata.display.max_projection(array, dim=0, bounds=None, title=None, cmap='gray', file=None)[source]

Projection of maximum values.

flexdata.display.mesh(stl_mesh)[source]

Display an stl mesh. Use flexCompute.generate_stl(volume) to generate mesh.

flexdata.display.min_projection(array, dim=0, title=None, cmap='gray', file=None)[source]

Projection of minimum values.

flexdata.display.plot(x, y=None, semilogy=False, title=None, legend=None)[source]

A standard 2D plot.

flexdata.display.plot2d(x, y=None, semilogy=False, title=None, legend=None)[source]

A standard 2D plot.

flexdata.display.plot3d(x, y, z, connected=False, title=None)[source]

Plot a 3D line or a scatter plot.

flexdata.display.projection(array, dim=1, bounds=None, title=None, cmap='gray', file=None)[source]

A simple projection of the volume along one of the dimensions.

flexdata.display.pyqt_graph(array, dim=0, title=None)[source]

Create a PYQT window to display a 3D arrayset.

flexdata.display.slice(array, index=None, dim=0, bounds=None, title=None, cmap='gray', file=None, cbar=True)[source]