nabu.estimation
source package nabu.estimation
Classes
-
AlignmentBase — Alignment basic functions.
-
CenterOfRotation — Alignment basic functions.
-
CenterOfRotationSlidingWindow — Alignment basic functions.
-
CenterOfRotationGrowingWindow — Alignment basic functions.
-
CenterOfRotationAdaptiveSearch — 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 )
-
CameraFocus — Alignment basic functions.
-
CameraTilt — Alignment basic functions.
-
DetectorTranslationAlongBeam — Alignment basic functions.
Functions
-
estimate_flat_distortion — Estimate the wavefront distortion on a flat image, from another image.
source class AlignmentBase(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
refine_max_position_2d — Computes the sub-pixel max position of the given function sampling.
-
refine_max_position_1d — Computes the sub-pixel max position of the given function sampling.
-
extract_peak_region_2d — Extracts a region around the maximum value.
-
extract_peak_regions_1d — Extracts a region around the maximum value.
-
close_plot_window — Close a plot window. Applicable only if the class was instantiated with verbose=True.
-
close_last_plot_windows — Close the last "n" plot windows. Applicable only if the class was instanciated with verbose=True.
source staticmethod AlignmentBase.refine_max_position_2d(f_vals: np.ndarray, fy=None, fx=None)
Computes the sub-pixel max position of the given function sampling.
Parameters
-
f_vals : numpy.ndarray — Function values of the sampled points
-
fy : numpy.ndarray, optional — Vertical coordinates of the sampled points
-
fx : numpy.ndarray, 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
-
tuple(float, float) — Estimated (vertical, horizontal) function max, according to the coordinates in fy and fx.
source staticmethod AlignmentBase.refine_max_position_1d(f_vals, fx=None, return_vertex_val=False, return_all_coeffs=False)
Computes the sub-pixel max position of the given function sampling.
Parameters
-
f_vals : numpy.ndarray — Function values of the sampled points
-
fx : numpy.ndarray, 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
-
float — Estimated function max, according to the coordinates in fx.
source staticmethod AlignmentBase.extract_peak_region_2d(cc, peak_radius=1, cc_vs=None, cc_hs=None)
Extracts a region around the maximum value.
Parameters
-
cc : numpy.ndarray — Correlation image.
-
peak_radius : int, optional — The l_inf radius of the area to extract around the peak. The default is 1.
-
cc_vs : numpy.ndarray, optional — The vertical coordinates of
cc
. The default is None. -
cc_hs : numpy.ndarray, optional — The horizontal coordinates of
cc
. The default is None.
Returns
-
f_vals : numpy.ndarray — The extracted function values.
-
fv : numpy.ndarray — The vertical coordinates of the extracted values.
-
fh : numpy.ndarray — The horizontal coordinates of the extracted values.
source staticmethod AlignmentBase.extract_peak_regions_1d(cc, axis=-1, peak_radius=1, cc_coords=None)
Extracts a region around the maximum value.
Parameters
-
cc : numpy.ndarray — 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 : numpy.ndarray, optional — The coordinates of
cc
along the selected axis. The default is None.
Returns
-
f_vals : numpy.ndarray — The extracted function values.
-
fc_ax : numpy.ndarray — The coordinates of the extracted values, along the selected axis.
Raises
-
ValueError
source method AlignmentBase.close_plot_window(n, errors='raise')
Close a plot window. Applicable only if the class was instantiated with verbose=True.
Parameters
-
n : int — Figure number to close
-
errors : str, optional — What to do with errors. It can be either "raise", "log" or "ignore".
Raises
-
ValueError
source method AlignmentBase.close_last_plot_windows(n=1)
Close the last "n" plot windows. Applicable only if the class was instanciated with verbose=True.
Parameters
-
n : int, optional — Integer indicating how many plot windows should be closed.
source class CenterOfRotation(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_shift — Find the Center of Rotation (CoR), given two images.
source method CenterOfRotation.find_shift(img_1: np.ndarray, img_2: np.ndarray, side=None, 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, return_relative_to_middle=None)
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
-
float — Estimated center of rotation position from the center of the RoI in pixels.
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))
source class CenterOfRotationSlidingWindow(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_shift — Semi-automatically find the Center of Rotation (CoR), given two images or sinograms. Suitable for half-aquisition scan.
source method CenterOfRotationSlidingWindow.find_shift(img_1: np.ndarray, img_2: np.ndarray, side='center', window_width=None, roi_yxhw=None, median_filt_shape=None, peak_fit_radius=1, high_pass=None, low_pass=None, return_validity=False, return_relative_to_middle=None)
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.
Parameters and usage is the same as CenterOfRotation, except for the following two parameters.
Parameters
-
side : string or float, optional — Expected region of the CoR. Allowed values: 'left', 'center' or 'right'. Default is 'center'
-
window_width : int, optional — Width of window that will slide on the other image / part of the sinogram. Default is None.
Raises
-
ValueError
source class CenterOfRotationGrowingWindow(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_shift — Automatically find the Center of Rotation (CoR), given two images or sinograms. Suitable for half-aquisition scan.
source method CenterOfRotationGrowingWindow.find_shift(img_1: np.ndarray, img_2: np.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, return_relative_to_middle=None)
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.
Usage and parameters are the same as CenterOfRotationSlidingWindow, except for the following parameter.
Parameters
-
min_window_width : int, optional — Minimum window width that covers the common region of the two images / sinograms. Default is 11.
Raises
-
ValueError
source class CenterOfRotationAdaptiveSearch(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_shift — Find the Center of Rotation (CoR), given two images.
source method CenterOfRotationAdaptiveSearch.find_shift(img_1: np.ndarray, img_2: np.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, return_relative_to_middle=None)
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.
Usage and parameters are the same as CenterOfRotation, except for the following parameters.
Parameters
-
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.
Raises
-
ValueError
source class SinoCor(img_1, img_2, logger=None)
This class has 2 methods
- overlap. Find a rough estimate of COR
- accurate. Try to refine COR to 1/10 pixel
Methods
source staticmethod SinoCor.schift(mat, val)
source method SinoCor.overlap(side='right', window_width=None)
Compute COR by minimizing difference of circulating ROI
- side: preliminary knowledge if the COR is on right or left
- window_width: width of ROI that will slide on the other part of the sinogram by default, 20% of the width of the detector.
Raises
-
ValueError
source method SinoCor.accurate(neighborhood=7, shift_value=0.1)
refine the calculation around COR integer pre-calculated value The search will be executed in the defined neighborhood
Parameters
-
neighborhood : int — Parameter for accurate calculation in the vicinity of the rough estimate. It must be an odd number. 0.1 pixels float shifts will be performed over this number of pixel
source estimate_flat_distortion(flat, image, tile_size=100, interpolation_kind='linear', padding_mode='edge', correction_spike_threshold=None, logger=None)
Estimate the wavefront distortion on a flat image, from another image.
Parameters
-
flat : np.array — The flat-field image to be corrected
-
image : np.ndarray — The image to correlate the flat against.
-
tile_size : int — The wavefront corrections are calculated by correlating the image to the flat, region by region. The regions are tiles of size tile_size
-
interpolation_kind : "linear" or "cubic" — The interpolation method used for interpolation
-
padding_mode : string — Padding mode. Must be valid for np.pad when wavefront correction is applied, the corrections are first found for the tiles, which gives the shift at the center of each tiled. Then, to interpolate the corrections, at the positions f every pixel, on must add also the border of the extremal tiles. This is done by padding with a width of 1, and using the mode given 'padding_mode'.
-
correction_spike_threshold : float, optional — By default it is None and no spike correction is performed on the shifts grid which is found by correlation. If set to a float, a spike removal will be applied using such threshold
Returns
-
coordinates : np.ndarray — An array having dimensions (flat.shape[0], flat.shape[1], 2) where each coordinates[i,j] contains the coordinates of the position in the image "flat" which correlates to the pixel (i,j) in the image "im2".
source class CameraFocus(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_distance — Find the focal distance of the camera system.
-
find_scintillator_tilt — Finds the scintillator tilt and focal distance of the camera system.
source method CameraFocus.find_distance(img_stack: np.ndarray, img_pos: np.array, metric='std', roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None)
Find the focal distance of the camera system.
This routine computes the motor position that corresponds to having the scintillator on the focal plain of the camera system.
Parameters
-
img_stack : numpy.ndarray — A stack of images at different distances.
-
img_pos : numpy.ndarray — Position of the images along the translation axis
-
metric : string, optional — The property, whose maximize occurs at the focal position. Defaults to 'std' (standard deviation). All options are: 'std' | 'grad' | 'psd'
-
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
.
Returns
-
focus_pos : float — Estimated position of the focal plane of the camera system.
-
focus_ind : float — Image index of the estimated position of the focal plane of the camera system (starting from 1!).
Examples
Given the focal stack associated to multiple positions of the camera
focus motor called img_stack
, and the associated positions img_pos
,
the following code computes the highest focus position:
focus_calc = alignment.CameraFocus()
focus_pos, focus_ind = focus_calc.find_distance(img_stack, img_pos)
where focus_pos
is the corresponding motor position, and focus_ind
is the associated image position (starting from 1).
source method CameraFocus.find_scintillator_tilt(img_stack: np.ndarray, img_pos: np.array, regions_number=4, metric='std', roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None)
Finds the scintillator tilt and focal distance of the camera system.
This routine computes the mounting tilt of the scintillator and the motor position that corresponds to having the scintillator on the focal plain of the camera system.
The input is supposed to be a stack of square images, whose sizes are
multiples of the regions_number
parameter. If images with a different
size are passed, this function will crop the images. This also generates
a warning. To suppress the warning, it is suggested to specify a ROI
that satisfies those criteria (see examples).
The computed tilts tilt_vh
are in unit-length per pixel-size. To
obtain the tilts it is necessary to divide by the pixel-size:
tilt_vh_deg = np.rad2deg(np.arctan(tilt_vh / pixel_size))
The correction to be applied is:
tilt_corr_vh_deg = - np.rad2deg(np.arctan(tilt_vh / pixel_size))
The legacy octave macros computed the approximation of these values in radians:
tilt_corr_vh_rad = - tilt_vh / pixel_size
Note that pixel_size
should be in the same unit scale as img_pos
.
Parameters
-
img_stack : numpy.ndarray — A stack of images at different distances.
-
img_pos : numpy.ndarray — Position of the images along the translation axis
-
regions_number : int, optional — The number of regions to subdivide the image into, along each direction. Defaults to 4.
-
metric : string, optional — The property, whose maximize occurs at the focal position. Defaults to 'std' (standard deviation). All options are: 'std' | 'grad' | 'psd'
-
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 -> auto-suggest correct size.
-
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
.
Returns
-
focus_pos : float — Estimated position of the focal plane of the camera system.
-
focus_ind : float — Image index of the estimated position of the focal plane of the camera system (starting from 1!).
-
tilts_vh : tuple(float, float) — Estimated scintillator tilts in the vertical and horizontal direction respectively per unit-length per pixel-size.
Examples
Given the focal stack associated to multiple positions of the camera
focus motor called img_stack
, and the associated positions img_pos
,
the following code computes the highest focus position:
focus_calc = alignment.CameraFocus()
focus_pos, focus_ind, tilts_vh = focus_calc.find_scintillator_tilt(img_stack, img_pos)
tilt_corr_vh_deg = - np.rad2deg(np.arctan(tilt_vh / pixel_size))
or to keep compatibility with the old octave macros:
tilt_corr_vh_rad = - tilt_vh / pixel_size
For non square images, or images with sizes that are not multiples of
the regions_number
parameter, and no ROI is being passed, this function
will try to crop the image stack to the correct size.
If you want to remove the warning message, it is suggested to set a ROI
like the following:
regions_number = 4
img_roi = (np.array(img_stack.shape[1:]) // regions_number) * regions_number
img_roi = np.fmin(img_roi, img_roi.min())
focus_calc = alignment.CameraFocus()
focus_pos, focus_ind, tilts_vh = focus_calc.find_scintillator_tilt(
img_stack, img_pos, roi_yxhw=img_roi, regions_number=regions_number)
source class CameraTilt(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
compute_angle — Find the camera tilt, given two opposite images.
source method CameraTilt.compute_angle(img_1: np.ndarray, img_2: np.ndarray, method='1d-correlation', roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None)
Find the camera tilt, given two opposite images.
This method finds the tilt between the camera pixel columns and the rotation axis, by performing a 1-dimensional correlation between two opposite images.
The output of this function, allows to compute motor movements for aligning the camera tilt.
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).
-
method : str — Tilt angle computation method. Default is "1d-correlation" (traditional). All options are: - "1d-correlation": fastest, but works best for small tilts - "fft-polar": slower, but works well on all ranges of tilts
-
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
Raises
-
ValueError — In case images are not 2-dimensional or have different sizes.
Returns
-
cor_offset_pix : float — Estimated center of rotation position from the center of the RoI in pixels.
-
tilt_deg : float — Estimated camera tilt angle in degrees.
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, :, :])
tilt_calc = CameraTilt()
cor_offset, camera_tilt = tilt_calc.compute_angle(radio1, radio2)
Or for noisy images:
cor_offset, camera_tilt = tilt_calc.compute_angle(radio1, radio2, median_filt_shape=(3, 3))
source class DetectorTranslationAlongBeam(vert_fft_width=False, horz_fft_width=False, verbose=False, logger=None, data_type=np.float32, extra_options=None)
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.<span class="mkapi-tooltip" title="numpy._core.float32">float32</span>
— Computation data type.
Methods
-
find_shift — Find the vertical and horizontal shifts for translations of the detector along the beam direction.
source method DetectorTranslationAlongBeam.find_shift(img_stack: np.ndarray, img_pos: np.array, roi_yxhw=None, median_filt_shape=None, padding_mode=None, peak_fit_radius=1, high_pass=None, low_pass=None, return_shifts=False, use_adjacent_imgs=False)
Find the vertical and horizontal shifts for translations of the detector along the beam direction.
These shifts are in pixels-per-unit-translation, and they are due to the misalignment of the translation stage, with respect to the beam propagation direction.
To compute the vertical and horizontal tilt angles from the obtained shift_pix
:
tilt_deg = np.rad2deg(np.arctan(shift_pix * pixel_size))
where pixel_size
and and the input parameter img_pos
have to be
expressed in the same units.
Parameters
-
img_stack : numpy.ndarray — A stack of images (usually 4) at different distances
-
img_pos : numpy.ndarray — Position of the images along the translation axis
-
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_shifts : boolean, optional — Adds a third returned argument, containing the pixel shifts of each image with respect to the first one in the stack. Defaults to False.
-
use_adjacent_imgs : boolean, optional — Compute correlation between adjacent images. It can be used when dealing with large shifts, to avoid overflowing the shift. This option allows to replicate the behavior of the reference function
alignxc.m
However, it is detrimental to shift fitting accuracy. Defaults to False.
Returns
-
coeff_v : float — Estimated vertical shift in pixel per unit-distance of the detector translation.
-
coeff_h : float — Estimated horizontal shift in pixel per unit-distance of the detector translation.
-
shifts_vh : list, optional — The pixel shifts of each image with respect to the first image in the stack. Activated if return_shifts is True.
Examples
The following example creates a stack of shifted images, and retrieves the computed shift. Here we use a high-pass filter, due to the presence of some low-frequency noise component.
import numpy as np
import scipy as sp
import scipy.ndimage
from nabu.preproc.alignment import DetectorTranslationAlongBeam
tr_calc = DetectorTranslationAlongBeam()
stack = np.zeros([4, 512, 512])
# Add low frequency spurious component
for i in range(4):
stack[i, 200 - i * 10, 200 - i * 10] = 1
stack = sp.ndimage.filters.gaussian_filter(stack, [0, 10, 10.0]) * 100
# Add the feature
x, y = np.meshgrid(np.arange(stack.shape[-1]), np.arange(stack.shape[-2]))
for i in range(4):
xc = x - (250 + i * 1.234)
yc = y - (250 + i * 1.234 * 2)
stack[i] += np.exp(-(xc * xc + yc * yc) * 0.5)
# Image translation along the beam
img_pos = np.arange(4)
# Find the shifts from the features
shifts_v, shifts_h = tr_calc.find_shift(stack, img_pos, high_pass=1.0)
print(shifts_v, shifts_h)
( -2.47 , -1.236 )
and the following commands convert the shifts in angular tilts:
tilt_v_deg = np.rad2deg(np.arctan(shifts_v * pixel_size))
tilt_h_deg = np.rad2deg(np.arctan(shifts_h * pixel_size))
To enable the legacy behavior of alignxc.m
(correlation between adjacent images):
shifts_v, shifts_h = tr_calc.find_shift(stack, img_pos, use_adjacent_imgs=True)
To plot the correlation shifts and the fitted straight lines for both directions:
tr_calc = DetectorTranslationAlongBeam(verbose=True)
shifts_v, shifts_h = tr_calc.find_shift(stack, img_pos)