corrct.alignment package

Submodules

corrct.alignment.centering module

Volume centering classes.

@author: Nicola VIGANÒ, CEA-IRIG and ESRF - The European Synchrotron, Grenoble, France

class corrct.alignment.centering.RecenterVolume(proj_geom: ProjectionGeometry, angles_rad: Union[ndarray[Any, dtype[ScalarType]], _SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], precision: int = 2)[source]

Bases: object

Volume recentering class.

recenter_as(shifts_vu: ndarray[Any, dtype[ScalarType]], volume: ndarray[Any, dtype[ScalarType]], reference: ndarray[Any, dtype[ScalarType]], method: str = 'com') ndarray[Any, dtype[ScalarType]][source]

Recenter with respect to a given volume.

Parameters
  • shifts_vu (NDArray) – Current VU shifts.

  • volume (NDArray) – The volume to shift.

  • reference (NDArray) – The reference volume.

  • method (str, optional) – The method to use out of “com” | “xc” (cross-correlation), by default “com”

Returns

The corrected VU shifts.

Return type

NDArray

Raises

ValueError – In case of wrong method.

recenter_to(shifts_vu: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]], volume: ndarray[Any, dtype[ScalarType]], com_ref_zyx: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]]) ndarray[Any, dtype[ScalarType]][source]

Recenter to a given center-of-mass (CoM).

Parameters
  • shifts_vu (ArrayLike | NDArray) – The current VU shifts.

  • volume (NDArray) – The volume to shift.

  • com_ref_zyx (ArrayLike | NDArray) – The destination CoM.

Returns

The corrected VU shifts.

Return type

NDArray

corrct.alignment.fitting module

Fitting routines.

Created on Tue May 17 12:11:58 2022

@author: Nicola VIGANÒ, Computational Imaging group, CWI, The Netherlands, and ESRF - The European Synchrotron, Grenoble, France

corrct.alignment.fitting.extract_peak_region_nd(cc: ndarray[Any, dtype[floating]], peak_radius: int = 1, cc_coords: dtype[+ ScalarType]]]]] = None) tuple[numpy.ndarray[typing.Any, numpy.dtype[+ScalarType]], typing.Optional[typing.Sequence[numpy.ndarray[typing.Any, numpy.dtype[+ScalarType]]]]][source]

Extract a region around the maximum value.

Parameters
  • cc (NDArrayFloat) – Correlation image.

  • peak_radius (int, optional) – The l_inf radius of the area to extract around the peak. The default is 1.

  • cc_coords (ArrayLike, optional) – The coordinates of cc. The default is None.

Returns

  • f_vals (NDArrayFloat) – The extracted function values.

  • f_coords (Tuple[NDArrayFloat]) – The coordinates of the extracted values.

corrct.alignment.fitting.extract_peak_regions_1d(cc: ndarray[Any, dtype[floating]], axis: int = - 1, peak_radius: int = 1, cc_coords: Optional[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]]] = None) tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], typing.Optional[numpy.ndarray[typing.Any, numpy.dtype[+ScalarType]]]][source]

Extract a region around the maximum value.

Parameters
  • cc (NDArrayFloat) – Correlation image.

  • axis (int, optional) – Find the max values along the specified direction. The default is -1.

  • peak_radius (int, optional) – The l_inf radius of the area to extract around the peak. The default is 1.

  • cc_coords (ArrayLike, optional) – The coordinates of cc along the selected axis. The default is None.

Returns

  • f_vals (NDArrayFloat) – The extracted function values.

  • fc_ax (NDArrayFloat) – The coordinates of the extracted values, along the selected axis.

corrct.alignment.fitting.fit_shifts_u_sad(data_wu: ndarray[Any, dtype[floating]], proj_wu: ndarray[Any, dtype[floating]], error_norm: int = 1, decimals: int = 1) ndarray[Any, dtype[floating]][source]

Find the U shifts between two sets of lines, by means of the sum-of-absolute-difference (SAD).

Parameters
  • data_wu (NDArrayFloat) – The reference data.

  • proj_wu (NDArrayFloat) – The other data.

  • error_norm (int, optional) – The error norm to use, by default 1

  • decimals (int, optional) – The precision of the result, by default 1

Returns

A list of one shift for each row.

Return type

NDArrayFloat

corrct.alignment.fitting.fit_shifts_vu_xc(data_vwu: ndarray[Any, dtype[floating]], proj_vwu: ndarray[Any, dtype[floating]], pad_u: bool = False, normalize_fourier: bool = True, use_rfft: bool = True, stack_axis: int = - 2, decimals: int = 2) ndarray[Any, dtype[floating]][source]

Find the VU shifts of the projected data, through cross-correlation.

Parameters
  • data_vwu (NDArrayFloat) – The collected projection data.

  • proj_vwu (NDArrayFloat) – The forward-projected images from the reconstruction.

  • pad_u (bool, optional) – Pad the u coordinate. The default is False.

  • normalize_fourier (bool, optional) – Whether to normalize the Fourier representation of the cross-correlation. The default is True.

  • use_rfft (bool, optional) – Whether to use the rfft transform in place of the complex fft transform. The default is True.

  • stack_axis (int, optional) – The axis along which the VU images are stacked. The default is -2.

  • decimals (int, optional) – Decimals for the truncation of the sub-pixel The default is 2.

Returns

The VU shifts.

Return type

NDArrayFloat

corrct.alignment.fitting.fit_shifts_zyx_xc(ref_vol_zyx: ndarray[Any, dtype[floating]], rec_vol_zyx: ndarray[Any, dtype[floating]], pad_zyx: bool = False, normalize_fourier: bool = True, use_rfft: bool = True, decimals: int = 2) ndarray[Any, dtype[floating]][source]

Find the ZYX shifts of the volume, through cross-correlation.

Parameters
  • ref_vol_zyx (NDArrayFloat) – The reference volume.

  • rec_vol_zyx (NDArrayFloat) – The reconstructed volume to register.

  • pad_zyx (bool, optional) – Pad the ZYX coordinates. The default is False.

  • normalize_fourier (bool, optional) – Whether to normalize the Fourier representation of the cross-correlation. The default is True.

  • use_rfft (bool, optional) – Whether to use the rfft transform in place of the complex fft transform. The default is True.

  • decimals (int, optional) – Decimals for the truncation of the sub-pixel The default is 2.

Returns

The ZYX shifts.

Return type

NDArrayFloat

corrct.alignment.fitting.fit_sinusoid(angles: ndarray[Any, dtype[floating]], values: ndarray[Any, dtype[floating]], fit_l1: bool = False) tuple[float, float, float][source]

Fits a sinusoid to the given values.

Parameters
  • angles (NDArrayFloat) – Angles where to evaluate the sinusoid.

  • values (NDArrayFloat) – Values of the sinusoid.

  • fit_l1 (bool, optional) – Whether to use l1 fit instead of the l2 fit, by default False

Returns

The amplitude, phase and bias of the sinusoid.

Return type

Tuple[float, float, float]

corrct.alignment.fitting.refine_max_position_1d(f_vals: ndarray[Any, dtype[floating]], f_x: Optional[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]]] = None, return_vertex_val: bool = False, decimals: int = 2) Union[ndarray[Any, dtype[floating]], tuple[numpy.ndarray[Any, numpy.dtype[numpy.floating]], numpy.ndarray[Any, numpy.dtype[numpy.floating]]]][source]

Compute the sub-pixel max position of the given function sampling.

Parameters
  • f_vals (NDArrayFloat) – Function values of the sampled points

  • fx (ArrayLike, optional) – Coordinates of the sampled points

  • return_vertex_val (boolean, option) – Enables returning the vertex values. Defaults to False.

Raises

ValueError – In case position and values do not have the same size, or in case the fitted maximum is outside the fitting region.

Returns

Estimated function max, according to the coordinates in fx.

Return type

float

corrct.alignment.fitting.refine_max_position_2d(f_vals: ndarray[Any, dtype[floating]], fy: Optional[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]]] = None, fx: Optional[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[ScalarType]]]] = None) ndarray[Any, dtype[ScalarType]][source]

Compute the sub-pixel max position of the given function sampling.

Parameters
  • f_vals (NDArrayFloat) – Function values of the sampled points

  • fy (ArrayLike, optional) – Vertical coordinates of the sampled points

  • fx (ArrayLike, optional) – Horizontal coordinates of the sampled points

Raises

ValueError – In case position and values do not have the same size, or in case the fitted maximum is outside the fitting region.

Returns

Estimated (vertical, horizontal) function max, according to the coordinates in fy and fx.

Return type

tuple(float, float)

corrct.alignment.fitting.sinusoid(x: Union[ndarray[Any, dtype[floating]], float], a: Union[ndarray[Any, dtype[floating]], float], p: Union[ndarray[Any, dtype[floating]], float], b: Union[ndarray[Any, dtype[floating]], float]) ndarray[Any, dtype[floating]][source]

Compute the values of a sine function.

Parameters
  • x (NDArrayFloat | float) – The independent variable.

  • a (NDArrayFloat | float) – The amplitude of the sine.

  • p (NDArrayFloat | float) – The phase of the sine.

  • b (NDArrayFloat | float) – The bias of the sine.

Returns

The computed values.

Return type

NDArrayFloat

corrct.alignment.shifts module

Detector shifts finding classes.

@author: Nicola VIGANÒ, Computational Imaging group, CWI, The Netherlands, and ESRF - The European Synchrotron, Grenoble, France

class corrct.alignment.shifts.DetectorShiftsBase(data_dvwu: ndarray[Any, dtype[floating]], rot_angle_rad: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[floating]]], *, data_format: str = 'dvwu', data_mask_dvwu: Optional[ndarray[Any, dtype[ScalarType]]] = None, borders_dvwu: dict = {'d': None, 'u': None, 'v': None, 'w': None}, max_shifts: Optional[Union[float, ndarray[Any, dtype[floating]]]] = None, precision_decimals: int = 2, verbose: bool = True)[source]

Bases: object

Compute the detector shifts for a given dataset.

angles_rad: ndarray[Any, dtype[floating]]
data_vwu: ndarray[Any, dtype[floating]]
class corrct.alignment.shifts.DetectorShiftsPRE(data_dvwu: ndarray[Any, dtype[floating]], rot_angle_rad: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[floating]]], *, data_format: str = 'dvwu', data_mask_dvwu: Optional[ndarray[Any, dtype[ScalarType]]] = None, borders_dvwu: dict = {'d': None, 'u': None, 'v': None, 'w': None}, max_shifts: Optional[Union[float, ndarray[Any, dtype[floating]]]] = None, precision_decimals: int = 2, verbose: bool = True)[source]

Bases: DetectorShiftsBase

Compute the pre-alignment detector shifts for a given dataset.

angles_rad: ndarray[Any, dtype[floating]]
data_vwu: ndarray[Any, dtype[floating]]
fit_u(fit_l1: bool = False, background: Optional[Union[float, ndarray[Any, dtype[ScalarType]]]] = None, method: str = 'com') tuple[numpy.ndarray[Any, numpy.dtype[numpy.floating]], float][source]

Compute the pre-alignment shifts for the horizontal dimension.

The pre-alignment shifts, and center-of-rotation (CoR) are computed by fitting a sinusoid to the centers of mass of each angle in the sinogram. The bias of the sinusoid corresponds to the CoR, while the deviations from the fitted curve correspond to the shifts.

Parameters
  • fit_l1 (bool, optional) – Computes the l1-min fit of the sinusoid, by default False.

  • background (float | NDArray | None, optional) – Removes the given background, by default None.

  • method (str, optional) – The method used for the identification of the fiducial marker position. Options are “com” (center-of-mass) | “max” (maximum value), by default “com”.

Returns

The shifts and the CoR.

Return type

Tuple[NDArrayFloat, float]

fit_v(use_derivative: bool = True, use_rfft: bool = True, normalize_fourier: bool = True) ndarray[Any, dtype[floating]][source]

Compute the pre-alignment vertical shifts of a 3D dataset.

The pre-alignment shifts are computed by cross-correlation of one projection against the others. The projections are integrated in the horizontal direction.

In the vertical direction, it is suggested to use some high pass filter. The default option is to use the derivates of the intensity profiles.

Parameters
  • use_derivative (bool, optional) – Whether to use the derivate of the vertical profile, by default True

  • use_rfft (bool, optional) – Whether to use the rfft transform for the cross-correlation, by default True

  • normalize_fourier (bool, optional) – Whether to normalize the cross-correlation in Fourier space, by default True

Returns

The vertical shifts.

Return type

NDArrayFloat

Raises

ValueError – If the dataset is 2D.

class corrct.alignment.shifts.DetectorShiftsXC(data_dvwu: ndarray[Any, dtype[floating]], rot_angle_rad: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]], ndarray[Any, dtype[floating]]], *, data_format: str = 'dvwu', data_mask_dvwu: Optional[ndarray[Any, dtype[ScalarType]]] = None, borders_dvwu: dict = {'d': None, 'u': None, 'v': None, 'w': None}, max_shifts: Optional[Union[float, ndarray[Any, dtype[floating]]]] = None, precision_decimals: int = 2, verbose: bool = True)[source]

Bases: DetectorShiftsBase

Compute the center-of-rotation for a given dataset, by cross correlation.

angles_rad: ndarray[Any, dtype[floating]]
data_vwu: ndarray[Any, dtype[floating]]
find_shifts_vu(data_dvwu: ndarray[Any, dtype[floating]], proj_dvwu: ndarray[Any, dtype[floating]], use_derivative: bool = False) ndarray[Any, dtype[floating]][source]

Find shifts between two images or sets of lines.

Parameters
  • data_dvwu (NDArrayFloat) – The reference data.

  • proj_dvwu (NDArrayFloat) – The other data.

  • use_derivative (bool, optional) – Whether to use derivatives over the horizontal (U) coordinate, by default False.

Returns

The shifts in vertical (optional) and horizontal coordinates ([V]U).

Return type

NDArrayFloat

fit_u_180() float[source]

Find the center-of-rotation, using the 0 and 180 degrees projections.

Returns

The center-of-rotation.

Return type

float

fit_u_360() float[source]

Find the center of rotation over a 360 degrees scan, by taking the average of the 0-180 over all pairs of angles.

Returns

The center-of-rotation.

Return type

float

fit_vu(fit_l1: bool = False) ndarray[Any, dtype[ScalarType]][source]

Compute the pre-alignment vertical and horizontal shifts, using cross-correlation.

Parameters

fit_l1 (bool, optional) – Computes the l1-min fit of the sinusoid, by default False.

Returns

Pre-alignment shifts in VU coordinates.

Return type

NDArray

fit_vu_accum_drifts(ref_data_dvwu: Optional[ndarray[Any, dtype[floating]]] = None) ndarray[Any, dtype[ScalarType]][source]

Fit static image drifts.

Parameters

ref_data_dvwu (Optional[NDArrayFloat], optional) – Reference image, by default None. If None, the first image in the data stack will be used.

Returns

The shifts of the image stack.

Return type

NDArray

Raises

ValueError – When the number of reference images is either too many or not enough.

Module contents

Package for data alignment.