nabu.estimation.cor module

class nabu.estimation.cor.CenterOfRotation(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: AlignmentBase

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

find_shift(img_1: ndarray, img_2: ndarray, shift_axis: int = -1, roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None, return_validity=False, cor_options=None)[source]

Find the Center of Rotation (CoR), given two images.

This method finds the half-shift between two opposite images, by means of correlation computed in Fourier space.

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • shift_axis (int) – Axis along which we want the shift to be computed. Default is -1 (horizontal).

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • peak_fit_radius (int, optional) – Radius size around the max correlation pixel, for sub-pixel fitting. Minimum and default value is 1.

  • low_pass (float or sequence of two floats) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) – High-pass filter properties, as described in nabu.misc.fourier_filters

  • return_validity (a boolean, defaults to false) – if set to True adds a second return value which may have three string values. These values are “unknown”, “sound”, “questionable”. It will be “uknown” if the validation method is not implemented and it will be “sound” or “questionable” if it is implemented.

Raises:

ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotation()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3))
class nabu.estimation.cor.CenterOfRotationSlidingWindow(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: CenterOfRotation

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

find_shift(img_1: ndarray, img_2: ndarray, side, window_width=None, roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None, return_validity=False, cor_options=None)[source]

Semi-automatically find the Center of Rotation (CoR), given two images or sinograms. Suitable for half-aquisition scan.

This method finds the half-shift between two opposite images, by minimizing difference over a moving window.

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • side (string) – Expected region of the CoR. Allowed values: ‘left’, ‘center’ or ‘right’.

  • window_width (int, optional) – Width of window that will slide on the other image / part of the sinogram. Default is None.

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • peak_fit_radius (int, optional) – Radius size around the max correlation pixel, for sub-pixel fitting. Minimum and default value is 1.

  • low_pass (float or sequence of two floats) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) – High-pass filter properties, as described in nabu.misc.fourier_filters

  • return_validity (a boolean, defaults to false) – if set to True adds a second return value which may have three string values. These values are “unknown”, “sound”, “questionable”. It will be “uknown” if the validation method is not implemented and it will be “sound” or “questionable” if it is implemented.

Raises:

ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotationSlidingWindow()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3))
class nabu.estimation.cor.CenterOfRotationGrowingWindow(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: CenterOfRotation

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

find_shift(img_1: ndarray, img_2: ndarray, side='all', min_window_width=11, roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None, return_validity=False)[source]

Automatically find the Center of Rotation (CoR), given two images or sinograms. Suitable for half-aquisition scan.

This method finds the half-shift between two opposite images, by minimizing difference over a moving window.

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • side (string, optional) – Expected region of the CoR. Allowed values: ‘left’, ‘center’, ‘right’, or ‘all’. Default is ‘all’.

  • min_window_width (int, optional) – Minimum window width that covers the common region of the two images / sinograms. Default is 11.

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • peak_fit_radius (int, optional) – Radius size around the max correlation pixel, for sub-pixel fitting. Minimum and default value is 1.

  • low_pass (float or sequence of two floats) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) – High-pass filter properties, as described in nabu.misc.fourier_filters

  • return_validity (a boolean, defaults to false) – if set to True adds a second return value which may have three string values. These values are “unknown”, “sound”, “questionable”. It will be “uknown” if the validation method is not implemented and it will be “sound” or “questionable” if it is implemented.

Raises:

ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotationGrowingWindow()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3))
class nabu.estimation.cor.CenterOfRotationAdaptiveSearch(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: CenterOfRotation

This adaptive method works by applying a gaussian which highlights, by apodisation, a region which can possibly contain the good center of rotation. The whole image is spanned during several applications of the apodisation. At each application the apodisation function, which is a gaussian, is moved to a new guess position. The lenght of the step, by which the gaussian is moved, and its sigma are obtained by multiplying the shortest distance from the left or right border with a self.step_fraction and self.sigma_fraction factors which ensure global overlapping. for each step a region around the CoR of each image is selected, and the regions of the two images are compared to calculate a cost function. The value of the cost function, at its minimum is used to select the best step at which the CoR is taken as final result. The option filtered_cost= True (default) triggers the filtering (according to low_pass and high_pass) of the two images which are used for he cost function. ( Note: the low_pass and high_pass options are used, if given, also without the filtered_cost option, by being passed to the base class CenterOfRotation )

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

sigma_fraction = 0.25
step_fraction = 0.16666666666666666
find_shift(img_1: ndarray, img_2: ndarray, roi_yxhw=None, median_filt_shape=None, padding_mode=None, high_pass=None, low_pass=None, margins=None, filtered_cost=True, return_validity=False)[source]

Find the Center of Rotation (CoR), given two images.

This method finds the half-shift between two opposite images, by means of correlation computed in Fourier space. A global search is done on on the detector span (minus a margin) without assuming centered scan conditions.

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • low_pass (float or sequence of two floats.) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) – High-pass filter properties, as described in nabu.misc.fourier_filters.

  • margins (None or a couple of floats or ints) – if margins is None or in the form of (margin1,margin2) the search is done between margin1 and dim_x-1-margin2. If left to None then by default (margin1,margin2) = ( 10, 10 ).

  • filtered_cost (boolean.) – True by default. It triggers the use of filtered images in the calculation of the cost function.

  • return_validity (a boolean, defaults to false) – if set to True adds a second return value which may have three string values. These values are “unknown”, “sound”, “questionable”. It will be “uknown” if the validation method is not implemented and it will be “sound” or “questionable” if it is implemented.

Raises:

ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotationAdaptiveSearch()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3), high_pass=20, low_pass=1   )
class nabu.estimation.cor.CenterOfRotationFourierAngles(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: CenterOfRotation

This CoR estimation algo is proposed by V. Valls (BCU). It is based on the Fourier transform of the columns on the sinogram. It requires an initial guesss of the CoR wich is retrieved from dataset_info.dataset_scanner.estimated_cor_from_motor. It is assumed in mm and pixel size in um. Options are (for the moment) hard-coded in the SinoCORFinder.cor_finder.extra_options dict.

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

gaussian(p, x)[source]
tukey(p, x)[source]
sinlet(p, x)[source]
find_shift(img_1, img_2, angles, side, roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None)[source]

Find the Center of Rotation (CoR), given two images.

This method finds the half-shift between two opposite images, by means of correlation computed in Fourier space.

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • shift_axis (int) – Axis along which we want the shift to be computed. Default is -1 (horizontal).

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • peak_fit_radius (int, optional) – Radius size around the max correlation pixel, for sub-pixel fitting. Minimum and default value is 1.

  • low_pass (float or sequence of two floats) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) – High-pass filter properties, as described in nabu.misc.fourier_filters

  • return_validity (a boolean, defaults to false) – if set to True adds a second return value which may have three string values. These values are “unknown”, “sound”, “questionable”. It will be “uknown” if the validation method is not implemented and it will be “sound” or “questionable” if it is implemented.

Raises:

ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotation()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3))
class nabu.estimation.cor.CenterOfRotationOctaveAccurate(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=<class 'numpy.float32'>, extra_options=None, cor_options=None)[source]

Bases: AlignmentBase

This is a Python implementation of Octave/fastomo3/accurate COR estimator. The Octave ‘accurate’ function is renamed local_correlation. The Nabu standard find_shift has the same API as the other COR estimators (sliding, growing…)

The class inherits directly from AlignmentBase.

Alignment basic functions.

Parameters:
  • vert_fft_width (boolean, optional) –

    If True, restrict the vertical size to a power of 2:

    >>> new_v_dim = 2 ** math.floor(math.log2(v_dim))
    

  • horz_fft_width (boolean, optional) –

    If True, restrict the horizontal size to a power of 2:

    >>> new_h_dim = 2 ** math.floor(math.log2(h_dim))
    

  • verbose (boolean, optional) – When True it will produce verbose output, including plots.

  • data_type (numpy.float32) – Computation data type.

find_shift(img_1: ndarray, img_2: ndarray, side: str, roi_yxhw=None, median_filt_shape=None, padding_mode=None, low_pass=0.01, high_pass=None)[source]

Automatically finds the Center of Rotation (CoR), given two images (projections/radiographs). Suitable for half-aquisition scan.

This method finds the half-shift between two opposite images, by minimizing the variance of small ROI around a global COR estimate (obtained by maximizing Fourier-space computed global correlations).

The output of this function, allows to compute motor movements for aligning the sample rotation axis. Given the following values:

  • L1: distance from source to motor

  • L2: distance from source to detector

  • ps: physical pixel size

  • v: output of this function

displacement of motor = (L1 / L2 * ps) * v

Parameters:
  • img_1 (numpy.ndarray) – First image

  • img_2 (numpy.ndarray) – Second image, it needs to have been flipped already (e.g. using numpy.fliplr).

  • side (string) – Expected region of the CoR. Must be ‘center’ in that case.

  • roi_yxhw ((2, ) or (4, ) numpy.ndarray, tuple, or array, optional) – 4 elements vector containing: vertical and horizontal coordinates of first pixel, plus height and width of the Region of Interest (RoI). Or a 2 elements vector containing: plus height and width of the centered Region of Interest (RoI). Default is None -> deactivated. The ROI will be used for the global estimate.

  • median_filt_shape ((2, ) numpy.ndarray, tuple, or array, optional) – Shape of the median filter window. Default is None -> deactivated.

  • padding_mode (str in numpy.pad's mode list, optional) –

    Padding mode, which determines the type of convolution. If None or ‘wrap’ are passed, this resorts to the traditional circular convolution. If ‘edge’ or ‘constant’ are passed, it results in a linear convolution. Default is the circular convolution. All options are:

    None | ‘constant’ | ‘edge’ | ‘linear_ramp’ | ‘maximum’ | ‘mean’ | ‘median’ | ‘minimum’ | ‘reflect’ | ‘symmetric’ |’wrap’

  • low_pass (float or sequence of two floats) – Low-pass filter properties, as described in nabu.misc.fourier_filters

  • high_pass (float or sequence of two floats) –

    High-pass filter properties, as described in nabu.misc.fourier_filters

    Raises

  • ------

  • ValueError – In case images are not 2-dimensional or have different sizes.

Returns:

Estimated center of rotation position from the center of the RoI in pixels.

Return type:

float

Examples

The following code computes the center of rotation position for two given images in a tomography scan, where the second image is taken at 180 degrees from the first.

>>> radio1 = data[0, :, :]
... radio2 = np.fliplr(data[1, :, :])
... CoR_calc = CenterOfRotationOctaveAccurate()
... cor_position = CoR_calc.find_shift(radio1, radio2)

Or for noisy images:

>>> cor_position = CoR_calc.find_shift(radio1, radio2, median_filt_shape=(3, 3))