toupy.registration package¶
Submodules¶
toupy.registration.registration module¶
registration_gd.py¶
Tomographic data alignment using gradient descent optimisation.
Vertical alignment : Adam optimiser on the L2 mass-fluctuation cost. Horizontal alignment: Adam optimiser on the L2 sinogram-consistency cost.
Both replace the original discrete line-search / parabolic-fit approach (_search_vshift_direction / _search_hshift_direction) with a proper gradient-based update rule. All other helper functions are unchanged.
- toupy.registration.registration.alignprojections_horizontal(sinogram, theta, shiftstack, **params)[source]¶
Horizontal alignment of projections by tomographic consistency.
Shifts are optimised with a Newton gradient-descent algorithm that minimises the L2 distance between the measured sinogram and a synthetic sinogram computed from a filtered-back-projection reconstruction. Supports a multi-stage frequency-cutoff schedule and an optional spatial multiresolution warm-start.
- Parameters:
sinogram (ndarray, shape (nr, nc)) – Measured sinogram (rows = detector pixels, columns = projections).
theta (ndarray, shape (nc,)) – Projection angles in radians.
shiftstack (ndarray, shape (2, n)) – Current shift estimates
[vertical_shifts, horizontal_shifts]. Only the horizontal row (index 1) is updated.**params –
Algorithm parameters. Required keys:
- maxitint
Maximum number of outer iterations per stage.
- pixtolfloat
Convergence tolerance in pixels.
- freqcutofffloat
Fraction of Nyquist used as the sinogram low-pass cut-off.
- shiftmethstr
Interpolation method (
'linear','fourier','spline').- circlebool
Apply a circular mask to the reconstruction.
- cliplowfloat or None
Lower clipping threshold for the reconstruction values.
- cliphighfloat or None
Upper clipping threshold for the reconstruction values.
- derivativesbool
Whether projections are phase-gradient images.
- calc_derivativesbool
Whether to compute derivatives on the synthetic sinogram.
- freqcutoff_schedulelist of float, optional
Sequence of
freqcutoffvalues (coarsest first). Each stage warm-starts the next. Default is[params['freqcutoff']](single stage).- multiresolutionbool, optional
Spatial multiresolution warm-start at stage 0. Default
False.- mr_factorint, optional
Down-sampling factor for the spatial warm-start. Default
2.- n_coarse_iterint, optional
Number of coarse iterations in the warm-start. Default
5.- rtolfloat, optional
Relative improvement threshold for early stopping. Default
0.- silentbool, optional
If
True, suppress all matplotlib output (required for sub-process workers). DefaultFalse.
- Returns:
shiftstack – Updated shift array with optimised horizontal shifts in row 1.
- Return type:
ndarray, shape (2, n)
- toupy.registration.registration.alignprojections_vertical(input_stack, shiftstack, **params)[source]¶
Vertical alignment of projections using the mass-fluctuation approach.
Shifts are optimised with a Newton gradient-descent algorithm that minimises the L2 distance between each projection’s vertical mass-fluctuation signal and the stack mean.
- Parameters:
input_stack (ndarray, shape (n, nr, nc)) – Stack of projection images to be aligned.
shiftstack (ndarray, shape (2, n)) – Initial shift estimates
[vertical_shifts, horizontal_shifts]. Modified in-place and also returned.**params –
Algorithm parameters. Required keys:
- maxitint
Maximum number of outer iterations.
- pixtolfloat
Convergence tolerance in pixels.
- deltaxint
Horizontal margin (pixels) to exclude from the ROI.
- limsylist of int or None
Explicit
[row_start, row_end]row limits.Noneuses the full image height.- shiftmethstr
Interpolation method (
'linear','fourier','spline').- polyorderint
Polynomial order for baseline removal in the fluctuation signal.
- alignxbool, optional
If
True, estimate horizontal shifts from the centre-of-mass before the vertical loop. DefaultFalse.- subpixelbool, optional
Ignored (always runs at sub-pixel precision). Default
True.
- Returns:
shiftstack (ndarray, shape (2, n)) – Optimised shift array.
output_stack (ndarray, shape (n, nr, nc)) – Vertically aligned projection stack.
- toupy.registration.registration.center_of_mass_stack(input_stack, lims, shiftstack, shift_method='fourier')[source]¶
Compute the centre-of-mass position for each projection in the stack.
- Parameters:
input_stack (ndarray, shape (n, nr, nc)) – Stack of projection images.
lims (tuple of array_like) –
(limrow, limcol)— row and column index arrays defining the region of interest.shiftstack (ndarray, shape (2, n)) – Current shift estimates used to pre-align each projection before computing the centre of mass.
shift_method (str, optional) – Interpolation method for the shift operation. Default
'fourier'.
- Returns:
Array
[centerx, centery]wherecenterx[i]andcentery[i]are the horizontal and vertical centre-of-mass offsets (in pixels) for projectioni.- Return type:
ndarray, shape (2, n)
- toupy.registration.registration.compute_aligned_horizontal(input_stack, shiftstack, shift_method='linear', **kwargs)[source]¶
Horizontal-only alignment (copy variant).
**kwargsabsorbs any extra keys when called ascompute_aligned_horizontal(..., **params).shift_methoddefaults toparams["shiftmeth"]if that key is present andshift_methodwas not supplied explicitly.
- toupy.registration.registration.compute_aligned_sino(input_sino, shiftslice, shift_method='linear')[source]¶
Compute the aligned sinogram given the correction for object positions.
- Parameters:
input_sino (array_like) – Sinogram to be shifted.
shiftslice (array_like) – Per-projection horizontal shifts (1, n).
shift_method (str) – ‘linear’, ‘fourier’, or ‘spline’.
- Returns:
output_sino – Aligned sinogram.
- Return type:
array_like
- toupy.registration.registration.compute_aligned_stack(input_stack, shiftstack, shift_method='linear')[source]¶
Compute the aligned stack given the correction for object positions.
- Parameters:
input_stack (array_like) – Stack of images to be shifted.
shiftstack (array_like) – Array of object motion corrections (2, n).
shift_method (str) – ‘linear’, ‘fourier’, or ‘spline’.
- Returns:
output_stack – Aligned image stack.
- Return type:
array_like
- toupy.registration.registration.estimate_rot_axis(input_array, theta, **params)[source]¶
Interactively estimate the rotation-axis offset.
Opens a GUI figure showing the reconstructed slice and sinogram for the current offset (
params["rot_axis_offset"]). Enter a new integer value in the text box and click Update to recompute. Click Confirm when satisfied.- Parameters:
input_array (ndarray, shape (n, nr, nc)) – Stack of derivative projections.
theta (ndarray) – Projection angles (degrees).
**params – Must contain:
slicenum,rot_axis_offset,colormap,cliplow,cliphigh,sinolow,sinohigh,filtertype,freqcutoff,circle,algorithm,derivatives,calc_derivatives.
- Returns:
Terminal mode –
rot_axis_offset(int) — the confirmed offset value.Jupyter mode (two-cell workflow) –
picker(_RotAxisPicker) — accesspicker.rot_axis_offsetin the next cell after clicking Confirm.
- toupy.registration.registration.oneslicefordisplay(sinogram, theta, **params)[source]¶
Reconstruct and display one tomographic slice.
Pass updated values via
params(e.g.params["freqcutoff"] = 0.5) and re-call to change reconstruction settings.
- toupy.registration.registration.refine_horizontalalignment(input_stack, theta, shiftstack, **params)[source]¶
Run one refinement pass of horizontal alignment.
Refinement uses the current
freqcutoffonly —freqcutoff_scheduleandmultiresolutionare ignored because the input shifts are already a good warm-start. Useparams["rtol"](e.g. 1e-3) to stop early when the improvement per iteration becomes negligible.To run multiple refinement passes or change parameters between passes, update
paramsin your script/notebook and call this function again.- Parameters:
input_stack (ndarray, shape (n, nr, nc)) – Stack of derivative projections.
theta (ndarray) – Projection angles.
shiftstack (ndarray, shape (2, n)) – Current shift array (modified in-place and returned).
**params – Alignment parameters forwarded to
alignprojections_horizontal().
- Returns:
shiftstack (ndarray)
params (dict)
- toupy.registration.registration.register_2Darrays(image1, image2, subpixel=False)[source]¶
Image registration using phase cross-correlations.
- Parameters:
image1 (array_like) – Reference image.
image2 (array_like) – Image to be shifted relative to image1.
subpixel (bool, optional) – If
False(default) use pixel-precision registration. IfTrueuse sub-pixel precision (upsample factor = 100).
- Returns:
shift (list of floats) – [row_shift, col_shift].
diffphase (float) – Phase difference between the two images.
offset_image2 (array_like) – Shifted image2 aligned to image1.
- toupy.registration.registration.tomoconsistency_multiple(input_stack, theta, shiftstack, **params)[source]¶
Tomographic consistency alignment over multiple sinogram slices.
Each slice is aligned independently (using the same
shiftstackwarm-start) and the resulting per-slice horizontal shifts are averaged to produce a robust final estimate. Slices can be processed in parallel viaProcessPoolExecutor.- Parameters:
input_stack (ndarray, shape (n, nr, nc)) – Full 3-D projection stack.
theta (ndarray, shape (n,)) – Projection angles in radians.
shiftstack (ndarray, shape (2, n)) – Current shift estimates used as a warm-start for every slice. The horizontal row (index 1) is updated with the averaged result.
**params –
Algorithm parameters forwarded to
alignprojections_horizontal(). Extra keys:- slicenumint
Central slice index.
- n_slices_tcint, optional
Number of slices to process (centred on
slicenum). Default10.- n_workers_tcint, optional
Number of parallel worker processes. Default
max(1, cpu_count // 2). Pass1for sequential.
- Returns:
shiftstack – Updated shift array with the averaged horizontal shifts in row 1 (or the original shifts if the user declines the result).
- Return type:
ndarray, shape (2, n)
- toupy.registration.registration.vertical_fluctuations(input_stack, lims, shiftstack, shift_method='fourier', polyorder=2)[source]¶
Compute the vertical mass-fluctuation signal for each projection.
The signal is obtained by integrating each (shifted) projection over the horizontal axis within the region of interest, then subtracting a polynomial baseline fit.
- Parameters:
input_stack (ndarray, shape (n, nr, nc)) – Stack of projection images.
lims (tuple of array_like) –
(rows, cols)— row and column index arrays for the ROI.shiftstack (ndarray, shape (2, n)) – Current vertical and horizontal shift estimates.
shift_method (str, optional) – Interpolation method for the shift operation. Default
'fourier'.polyorder (int, optional) – Polynomial order used for baseline removal. Default
2.
- Returns:
vert_fluct – Array of vertical fluctuation signals, one row per projection.
- Return type:
ndarray, shape (n, n_rows_roi)
- toupy.registration.registration.vertical_shift(input_array, lims, vstep, maxshift, shift_method='linear', polyorder=2)[source]¶
Compute the vertical mass-fluctuation signal for a single projection shifted by a given amount.
- Parameters:
input_array (ndarray, shape (nr, nc)) – Single projection image.
lims (tuple of array_like) –
(rows, cols)— row and column index arrays for the ROI.vstep (float) – Vertical shift to apply (pixels).
maxshift (int) – Safety margin (pixels) around the ROI to absorb border effects.
shift_method (str, optional) – Interpolation method for the shift operation. Default
'linear'.polyorder (int, optional) – Polynomial order used for baseline removal. Default
2.
- Returns:
shift_calc – Vertical fluctuation signal after shifting and polynomial baseline subtraction.
- Return type:
ndarray, shape (n_rows_roi,)
toupy.registration.shift module¶
- class toupy.registration.shift.ShiftFunc(**params)[source]¶
Bases:
VariablesCollection of shift functions for 1-D and 2-D arrays.
The desired method is selected at construction time via the
shiftmethparameter and exposed through__call__so that the same interface works regardless of the underlying algorithm.- Parameters:
**params –
Must contain:
- shiftmethstr
Shift algorithm to use. One of
'linear'(bilinear interpolation),'fourier'(FFT-based sub-pixel shift via pyFFTW), or'spline'(spline interpolation viascipy.ndimage.interpolation.shift()).
- __call__(*args)[source]¶
Apply the selected shift method to an array.
- Parameters:
*args –
- args[0]array_like
Input 1-D or 2-D array.
- args[1]float or tuple of float
Shift amplitude in pixels. For 1-D arrays pass a scalar; for 2-D arrays pass
(row_shift, col_shift).- args[2]str, optional
Padding mode (overrides the instance default).
- args[3]bool, optional
Trueto return a complex array;False(default) to return a real array.
- Returns:
Shifted array with the same shape as
args[0]. Returns the input unchanged when the shift is zero.- Return type:
ndarray
- shift_fft(input_array, shift)[source]¶
Performs pixel and subpixel shift (with wraping) using pyFFTW.
Since FFTW has efficient functions for array sizes which can be decompose in prime factor, the input_array is padded to the next fast size given by pyFFTW.next_fast_len. The padding is done in mode = ‘reflect’ by default to reduce border artifacts.
- Parameters:
- Returns:
output_array – Shifted image
- Return type:
array_like
- shift_linear(input_array, shift)[source]¶
Shifts an image with wrap around and bilinear interpolation
- Parameters:
- Returns:
output_array – Shifted image
- Return type:
array_like