Preprocessing Algorithms

leap_preprocessing_algorithms.gain_correction(leapct, g, air_scan, dark_scan, calibration_scans=None, ROI=None, badPixelMap=None, flux_response=None)

Performs gain correction

This function processes raw radiographs by subtracting off the dark current and correcting for the pixel-to-pixel gain variations which reduces ring artifacts

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (C contiguous float32 numpy array or torch tensor) – radiograph data

  • air_scan (C contiguous float32 numpy array or torch tensor) – air scan radiograph data; if not given, assumes that the input data is transmission data

  • dark_scan (C contiguous float32 numpy array or torch tensor) – dark scan radiograph data; if not given assumes the inputs have already been dark subtracted

  • calibrartion_scans (C contiguous 3D float32 numpy array or torch tensor) – calibration scan data

  • ROI (4-element integer array) – specifies a bounding box ([first row, last row, first column, last column]) for which to estimate a mean value for flux correction

  • badPixelMap (C contiguous float32 numpy array or torch tensor) – 2D bad pixel map (numRows x numCols) where a value of 1.0 marks a pixel as bad

  • flux_response (C contiguous 1D float32 numpy array or torch tensor) – transfer function of dark subtracted raw data

leap_preprocessing_algorithms.makeAttenuationRadiographs(leapct, g, air_scan=None, dark_scan=None, ROI=None, isAttenuationData=False)

Converts data to attenuation radiographs (flat fielding and negative log)

\[\begin{split}\begin{eqnarray} transmission\_data &=& (raw - dark\_scan) / (air\_scan - dark\_scan) \\ attenuation\_data &=& -log(transmission\_data) \end{eqnarray}\end{split}\]

This function assumes the input data is never attenuation data. See argument descriptions below for how the input data is treated.

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (C contiguous float32 numpy array or torch tensor) – radiograph data

  • air_scan (C contiguous float32 numpy array or torch tensor) – air scan radiograph data; if not given, assumes that the input data is transmission data

  • dark_scan (C contiguous float32 numpy array or torch tensor) – dark scan radiograph data; if not given assumes the inputs have already been dark subtracted

  • ROI (4-element integer array) – specifies a bounding box ([first row, last row, first column, last column]) for which to estimate a mean value for flux correction

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.badPixelCorrection(leapct, g, air_scan=None, dark_scan=None, badPixelMap=None, windowSize=3, isAttenuationData=True)

Removes bad pixels from CT projections

LEAP CT geometry parameters must be set prior to running this function and can be applied to any CT geometry type. This algorithm processes each projection independently and removes bad pixels specified by the user using a median filter

If no bad pixel map is provided, this routine will estimate it from the average of all projections.

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation or transmission projection data

  • badPixelMap (C contiguous float32 numpy array or torch tensor) – 2D bad pixel map (numRows x numCols) where a value of 1.0 marks a pixel as bad

  • windowSize (int) – the window size; can be 3, 5, or 7

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.outlierCorrection(leapct, g, threshold=0.03, windowSize=3, isAttenuationData=True)

Removes outliers (zingers) from CT projections

No LEAP parameters need to be set for this function to work and can be applied to any CT geometry type. This algorithm processes each projection independently and removes outliers by a thresholded median filter

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation or transmission projection data

  • threshold (float) – A pixel will be replaced by the median of its neighbors if |g - median(g)|/median(g) > threshold

  • windowSize (int) – The window size of the median filter applied is windowSize x windowSize

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.outlierCorrection_highEnergy(leapct, g, isAttenuationData=True)

Removes outliers (zingers) from CT projections

No LEAP parameters need to be set for this function to work and can be applied to any CT geometry type. This algorithm processes each projection independently and removes outliers by a series of three thresholded median filters

This outlier correction algorithm should most be used for MV CT or neutron CT where outliers effect a larger neighborhood of pixels. This function uses a three-stage thresholded median filter to remove outliers.

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation or transmission projection data

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.detectorDeblur_FourierDeconv(leapct, g, H, WienerParam=0.0, isAttenuationData=True)

Removes detector blur by fourier deconvolution

Parameters:
  • g (contiguous float32 numpy array or torch tensor) – attenuation or transmission projection data

  • H (2D contiguous float32 numpy array or torch tensor) – Magnitude of the frequency response of blurring psf, DC is at [0,0]

  • WienerParam (float) – Parameter for Wiener deconvolution, number should be between 0.0 and 1.0

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.detectorDeblur_RichardsonLucy(leapct, g, H, numIter=10, isAttenuationData=True)

Removes detector blur by Richardson-Lucy iterative deconvolution

Richardson-Lucy iterative deconvolution is developed for Poisson-distributed data and inherently preserve the non-negativity of the input. It uses the following update step, where t = transmission data

\[\begin{eqnarray} t_{n+1} &=& t_n H^T\left[ \frac{t_0}{Ht_n} \right] \end{eqnarray}\]
Parameters:
  • g (contiguous float32 numpy array or torch tensor) – attenuation or transmission projection data

  • H (2D contiguous float32 numpy array or torch tensor) – Magnitude of the frequency response of blurring psf, DC is at [0,0]

  • numIter (int) – Number of iterations

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.ringRemoval_fast(leapct, g, delta=0.01, beta=1000.0, numIter=30, maxChange=0.05)

Removes detector pixel-to-pixel gain variations that cause ring artifacts in reconstructed images

This algorithm estimates the rings by first averaging all projections. Then denoises this signal by minimizing the TV norm. Finally the gain correction map to correct the data is estimated by the difference of the TV-smoothed and averaged projection data (these are all 2D signals). This is summarized by the math equations below.

\[\begin{split}\begin{eqnarray} \overline{g} &:=& \frac{1}{numAngles}\sum_{angles} g \\ gain\_correction &:=& argmin_x \; \left[\frac{1}{2} \|x - \overline{g}\|^2 + \beta TV(x) \right] - \overline{g} \\ g\_corrected &:=& g + gain\_correction \end{eqnarray}\end{split}\]

Assumes the input data is in attenuation space. No LEAP parameters need to be set for this function to work and can be applied to any CT geometry type. This algorithm is effective at removing ring artifacts and runs fast, but can sometimes create new rings of tangents of sharp transitions.

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation projection data

  • delta (float) – The delta parameter of the Total Variation Functional

  • numIter (int) – Number of iterations

  • maxChange (float) – An upper limit on the maximum difference that can be applied to a detector pixels

  • beta (float) – The strength of the regularization

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.ringRemoval_median(leapct, g, threshold=0.0, windowSize=5, numIter=1)

Removes detector pixel-to-pixel gain variations that cause ring artifacts in reconstructed images

\[\begin{split}\begin{eqnarray} gain\_correction &:=& \frac{1}{numAngles}\sum_{angles} [MedianFilter2D(g) - g] \\ g\_corrected &:=& g + gain\_correction \end{eqnarray}\end{split}\]

Assumes the input data is in attenuation space. No LEAP parameters need to be set for this function to work and can be applied to any CT geometry type. This algorithm is effective at removing ring artifacts without creating new ring artifacts.

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation projection data

  • threshold (float) – The threshold for the thresholded median filter, where a pixel will be replaced by the median of its neighbors if |g - median(g)|/median(g) > threshold

  • windowSize (int) – The window size of the median filter applied is windowSize x windowSize

  • numIter (int) – Number of iterations

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.ringRemoval(leapct, g, delta=0.01, beta=10.0, numIter=30, maxChange=0.05)

Removes detector pixel-to-pixel gain variations that cause ring artifacts in reconstructed images

This algorithm estimates the gain correction necessary to remove ring artifacts by solving denoising the data by minimizing the TV norm with the gradient step determined by averaging the gradients over all angles.

Assumes the input data is in attenuation space. No LEAP parameters need to be set for this function to work and can be applied to any CT geometry type. This algorithm is effective at removing ring artifacts without creating new ring artifacts, but is more computationally expensive than ringRemoval_fast.

This algorithm works by minimizing the standard TV denoising cost function with gradient descent except we replace the gradient by averaging it over all projections. This ensures that the projections are smoothed by a function that is constant over projection angles. The cost function is given by

\[\begin{eqnarray} \frac{1}{2} \|x - g\|^2 + \beta TV(x) \end{eqnarray}\]
Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation projection data

  • delta (float) – The delta parameter of the Total Variation Functional

  • beta (float) – The strength of the regularization

  • numIter (int) – Number of iterations

Returns:

True if successful, False otherwise

leap_preprocessing_algorithms.transmission_shift(leapct, g, shift, isAttenuationData=True)

Subtracts constant from transmission data which is a simple method for scatter correction

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation projection data

  • shift (float) – the amount to subtract from the transmission data

  • isAttenuationData (bool) – True if g is attenuation data, False otherwise

Returns:

True is successful, False otherwise

leap_preprocessing_algorithms.parameter_sweep(leapct, g, values, param='centerCol', iz=None, algorithmName='FBP')

Performs single-slice reconstructions of several values of a given parameter

The CT geometry parameters and the CT volume parameters must be set prior to running this function.

The parameters to sweep are all standard LEAP CT geometry parameter names, except ‘tilt’ which is only available for cone- and modular-beam data. (note that the data g is not rotated, just the model of the detector orientation which is better because no interpolation is necessary).

Parameters:
  • leapct (tomographicModels object) – This is just needed to access LEAP algorithms

  • g (contiguous float32 numpy array or torch tensor) – attenuation projection data

  • values (list of floats) – the values to reconstruct with

  • param (string) – the name of the parameter to sweep; can be ‘centerCol’, ‘centerRow’, ‘tau’, ‘sod’, ‘sdd’, ‘tilt’, ‘vertical_shift’, ‘horizontal_shift’

  • iz (integer) – the z-slice index to perform the reconstruction; if not given, uses the central slice

  • algorithmName (string) – the name of the algorithm to use for reconstruction; can be ‘FBP’ or ‘inconsistencyReconstruction’

Returns:

stack of single-slice reconstructions (i.e., 3D numpy array or torch tensor) for all parameter values