API Reference

Holography UDFs

class libertem_holo.udf.HoloReconstructUDF(*, out_shape: tuple[int, int], sb_position: tuple[float, float], aperture: ndarray, precision: bool = True)[source]

Reconstruct off-axis electron holograms using a Fourier-based method.

Running run_udf() on an instance of this class will reconstruct a complex electron wave. Use the wave key to access the raw data in the result.

See Usage for detailed application example

New in version 0.3.0.

Examples

>>> from libertem_holo.base.filters import butterworth_disk
>>> shape = tuple(dataset.shape.sig)
>>> # NOTE: in real use, one would use `HoloParams.from_hologram` to
>>> # directly estimate parameters from the data, instead of specifying manually
>>> sb_position = [2, 3]
>>> sb_size = 4.4
>>> aperture = butterworth_disk(shape=shape, radius=sb_size)
>>> holo_udf = HoloReconstructUDF(
...     out_shape=shape,
...     sb_position=sb_position,
...     aperture=np.fft.fftshift(aperture),
... )
>>> wave = ctx.run_udf(dataset=dataset, udf=holo_udf)['wave'].data

Base holography functions

Reconstruction

Core reconstruction functions of LiberTEM-holo.

This includes both FFT-based and phase-shifting approaches.

libertem_holo.base.reconstr.display_fft_image(image: ndarray, sb_position: tuple, slice_fft: ndarray, mask: ndarray = 1, detail: bool = True) None[source]

Display an fft image.

This function helps to show the steps of the reconstruction and to define the best length and width of line mask

Parameters:
  • image (array) – the input image

  • sb_position (tuple) – the sideband position (y, x), referred to the non-shifted FFT.

  • fft (slice) – contain minimum and maximum value in y and x axis.

  • mask (array) – contain the aperture and line mask

  • detail (bolean) –

libertem_holo.base.reconstr.estimate_omega(image: ndarray, sideband_position: tuple[float, float], *, flip: bool = False) tuple[float, float][source]

Estimates the frequency carrier of the hologram.

Parameters:
  • image (array) – Holographic data array

  • sideband_position (tuple) – The sideband position (y, x), referred to the non-shifted FFT.

  • flip – ??? TODO

Returns:

omega – frequency carrier in y and x axis

Return type:

tuple

libertem_holo.base.reconstr.get_phase(hologram: ~numpy.ndarray, params: ~libertem_holo.base.utils.HoloParams, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Reconstruct hologram using HoloParams and extract and unwrap phase.

libertem_holo.base.reconstr.phase_offset_correction(aligned_stack, wtype: ~typing.Literal['weighted'] | ~typing.Literal['unweighted'] = 'weighted', threshold: float = 1e-12, return_stack: bool = False, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) tuple[ndarray, ndarray, ndarray | None][source]

This part of the code is to correct for the phase drift in the holograms due to the biprism drift with time. Since we are dealing with the phase of the image which wraps from -pi to pi a simple scalar correction of the phase doesnt work as the phase wraps around the phase axis. So the options for a solution would be either iterative solution which is computationally heavy or an eigenvalue solution which is implemented here.

We start with a stack of holograms acquired using Holoworks or stack acquisition and this function returns a phase corrected complex image and optionally a stack of the phase corrected complex images before averaging.

Adapted from a function by oleh.melnyk: https://github.com/Ptychography-4-0/ptychography/blob/master/src/ptychography40/stitching/stitching.py for more details, see https://arxiv.org/pdf/2005.02032.pdf

Parameters:
  • aligned_stack – Array of shape (N, sy, sx) where N is the number of complex images; should have dtype complex64 or complex128.

  • threshold – Minimum absolute value to be considered in finding the phase match for stitching.

  • wtype – Selected type of weights to be used in the angular synchronization There are 2 possible choices: ‘unweighted’ or ‘weighted’, depending on usage or not of the weights for angular synchronization (explained above)

  • return_stack – Also returns a stack shaped like the input stack, where the phase offset correction is applied.

  • xp – Either numpy or cupy for GPU support

libertem_holo.base.reconstr.reconstruct_bf(frame: ~numpy.ndarray, aperture: ~numpy.ndarray, slice_fft: tuple[slice, slice], *, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Reconstruct a brightfield image from a hologram.

Please use libertem_holo.base.filter.central_line_filter() to filter out fresnel fringes as appropriate.

libertem_holo.base.reconstr.reconstruct_direct(stack: ndarray, omega: tuple[float, float]) ndarray[source]

Reconstruct a stack of phase shifted holograms.

This is using the direct reconstruction method.

Parameters:
  • stack (array_like) – Stack of holograms

  • number_of_images (int) – The number of images inside the stack

  • omega (tuple) – frequency carrier in y and x axis

Returns:

phase final – the reconstructed phase image

Return type:

2D array

libertem_holo.base.reconstr.reconstruct_direct_euler(image: ndarray, omega: tuple[float, float]) ndarray[source]

Reconstruct a stack of phase shifted holograms.

This is using the direct reconstruction method by Ru et.al. 1994 (euler form)

Parameters:
  • image (array) – Holographic data array

  • omega (tuple) – frequency carrier in y and x axis

Returns:

phase final – the reconstructed phase image

Return type:

2D array

libertem_holo.base.reconstr.reconstruct_double_resolution(frames: ~numpy.ndarray, sb_pos: tuple[float, float], aperture: ~numpy.ndarray, slice_fft: tuple[slice, slice], *, precision: bool = True, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Reconstruct a stack of phase shifted holography with double resolution method.

Parameters:
  • frames (array_like) – Two holograms taken at a phase offset of π; shape (2, height, width)

  • sb_pos – The sideband position, for example as returned by estimate_sideband_position

  • aperture – A numpy or cupy array containing the aperture to apply in fourier space

  • slice_fft – A slice to crop to the selected output shape in fourier space, as returned by get_slice_fft.

  • precision – Defines precision of the reconstruction, True for complex128 for the resulting complex wave, otherwise results will be complex64

  • xp – Pass in either the numpy or cupy module to select CPU or GPU processing

Returns:

wav – the reconstructed complex image

Return type:

complex array

libertem_holo.base.reconstr.reconstruct_frame(frame: ~numpy.ndarray, sb_pos: tuple[float, float], aperture: ~numpy.ndarray, slice_fft: tuple[slice, slice], *, precision: bool = True, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Reconstruct a single hologram.

Parameters:
  • frame – A numpy or cupy array containing the input hologram

  • sb_pos – The sideband position, for example as returned by estimate_sideband_position

  • aperture – A numpy or cupy array containing the aperture to apply in fourier space

  • slice_fft – A slice to crop to the selected output shape in fourier space, as returned by get_slice_fft.

  • precision – Defines precision of the reconstruction, True for complex128 for the resulting complex wave, otherwise results will be complex64

  • xp – Pass in either the numpy or cupy module to select CPU or GPU processing

Image filtering and aperture building

Useful image filtering helpers.

libertem_holo.base.filters.butterworth_disk(shape: tuple[int, int], radius: float, order: int = 12, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Generate a filered disk-shaped aperture.

The edges are filtered with a butterworth filter of the given order.

Parameters:
  • shape – output shape of the aperture

  • radius – radius in pixels

  • order – order of the butterworth filter

  • xp – Either numpy or cupy

libertem_holo.base.filters.butterworth_line(shape: tuple[int, int], width: float, sb_position: tuple[int, int], length_ratio: float = 0.9, order: int = 12, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Generate a line filter.

This is useful to remove fresnel fringes.

Parameters:
  • shape – output shape of the aperture

  • width – width of the line in pixels

  • order – order of the butterworth filter

  • xp – Either numpy or cupy

libertem_holo.base.filters.central_line_filter(sb_position, out_shape, orig_shape, length_ratio=0.9, width=20, crop_to_out_shape=False)[source]

Return a line filter for the central band, that can be applied by multiplying it with the aperture.

libertem_holo.base.filters.clipped(img: ndarray, sigma: float = 6)[source]

Mask out outliers.

Return img, but with pixels deviating more than sigma from the mean masked out.

Useful for plotting:

>>> plt.imshow(img, vmax=np.max(clipped(img)))  
libertem_holo.base.filters.disk_aperture(out_shape: tuple[int, int], radius: float, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Generate a disk-shaped aperture

Parameters:
  • out_shape (interger) – Shape of the output array

  • radius (float) – Radius of the disk

  • xp – Either numpy or cupy

Returns:

2d array containing aperture

Return type:

aperture

libertem_holo.base.filters.exclusion_mask(img: ndarray, sigma: float = 6) ndarray[source]

Generate outlier mask.

Return a mask with True entries for pixels deviating more than sigma from the mean.

libertem_holo.base.filters.highpass(img: ndarray, sigma: float = 2) ndarray[source]

Return highpass by subtracting a gaussian lowpass filter.

libertem_holo.base.filters.line_filter(sb_position, out_shape, orig_shape, length_ratio=0.9, width=20, crop_to_out_shape=True)[source]

Return a line filter for the sideband.

It can be applied by multiplying it with the aperture. The filter has the zero frequency at the center of the image.

libertem_holo.base.filters.phase_unwrap(image)[source]

A phase_unwrap function that is unwrap the complex / wrapped phase image.

Parameters:

image (2d nd array) – Complex or Wrapped phase image

Return type:

2d nd array of the unwrapped phase image

libertem_holo.base.filters.remove_dead_pixels(img, sigma_lowpass=2.0, sigma_exclusion=6.0)[source]

Remove dead pixels.

Parameters:
  • img (np.array) – Input array

  • sigma_lowpass (float) – How much of the low frequencies should be removed before finding bad pixels

  • sigma_exclusion (float) – Pixels deviating more than this value from the mean will be removed

libertem_holo.base.filters.window_filter(input_array, window_type, window_shape)[source]

Apply window-based filter.

Return a filtered array with the same size of the input array

Parameters:
  • input_array (array) – Input array

  • window_type (string, float or tuple) – The type of window to be created. Any window type supported by scipy.signal.get_window is allowed here. See notes below for a current list, or the SciPy documentation for the version of SciPy on your machine.

  • window_shape (tuple of int or int) – The shape of the window. If an integer is provided, a 2D window is generated.

Notes

This function is based on scipy.signal.get_window and thus can access all of the window types available to that function (e.g., "hann", "boxcar"). Note that certain window types require parameters that have to be supplied with the window name as a tuple (e.g., ("tukey", 0.8)). If only a float is supplied, it is interpreted as the beta parameter of the Kaiser window. https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.windows.get_window.html

it is recommended to check the that after fft shift, the input array has value of 0 at the border.

Stack alignment

class libertem_holo.base.align.BiprismDeletionCorrelator(mask: ~numpy.ndarray, upsample_factor: int = 1, normalization: ~typing.Literal['phase'] | None = 'phase', xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation on low magnification while removing biprism.

classmethod get_masked(img, coords)[source]

Uses coordinates from plot_get_coords to create a mask of biprism.

classmethod plot_get_coords(img, coords_out)[source]

At low magnification, plot image of area with biprism visible. Click on edges of biprism to create the coordinates to mask it out, for cross correlation. First, click one edge of biprism from left side, then same edge, right side. Then, click other edge of biprism from left side, then right side. —–1——-2—– —–3——-4—–

class libertem_holo.base.align.BrightFieldCorrelator(holoparams: ~libertem_holo.base.utils.HoloParams, upsample_factor: int = 1, normalization: ~typing.Literal['phase'] | None = 'phase', xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation on bright field of hologram.

class libertem_holo.base.align.GradAngleCorrelator(holoparams: ~libertem_holo.base.utils.HoloParams, upsample_factor: int = 1, normalization: ~typing.Literal['phase'] | None = 'phase', xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation on gradient angle of phase image.

class libertem_holo.base.align.GradXYCorrelator(holoparams: ~libertem_holo.base.utils.HoloParams, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation on gradient x and Y, correlation maps summed.

class libertem_holo.base.align.ImageCorrelator(upsample_factor: int = 1, normalization: ~typing.Literal['phase'] | None = 'phase', hanning: bool = True, binning: int = 1, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation based image registration with some pre-filtering. Assumes real-space images as input.

class libertem_holo.base.align.NoopCorrelator[source]

Do nothing, successfully.

class libertem_holo.base.align.PhaseImageCorrelator(holoparams: ~libertem_holo.base.utils.HoloParams, upsample_factor: int = 1, normalization: ~typing.Literal['phase'] | None = 'phase', xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Cross correlation on reconstructed phase image.

class libertem_holo.base.align.RegResult(maximum, shift, corrmap)[source]
corrmap: ndarray

Alias for field number 2

maximum: tuple[float, float]

Alias for field number 0

shift: tuple[float, float]

Alias for field number 1

libertem_holo.base.align.align_stack(stack: ~numpy.ndarray, wave_stack: ~numpy.ndarray, static: ~numpy.ndarray | None, correlator: ~libertem_holo.base.align.Correlator | None = None, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) tuple[ndarray, ndarray, ndarray, ndarray][source]

Align stacks of N holograms.

Parameters:
  • stack – Stack of images with shape (N, h, w) which is used for the alignment. It should be of dtype float32 or float64, so for holograms, you may want to align on the abs(wave). You can also slice the original data, if you want to align to a region of interest.

  • wave_stack – The (complex) result of the reconstruction. This should have shape (N, h’, w’), meaning it should have the same number of complex images as the stack argument, but the 2d items in the stack can have a different shape.

  • static – A reference image to align against. If this is not given, the first image of the stack will be taken.

  • correlator – A Correlator instance. By default, an ImageCorrelator is used, with reasonable default parameters. If you need control over these, you can construct your own Correlator and pass it in. A requirement is that it has to work on the value given as the stack parameter.

Returns:

  • aligned_stack – The aligned stack. Same shape and dtype as wave_stack

  • shifts – The shifts that were applied to the stack. Shape (N, 2)

  • reference – The pre-processed reference image, as used in the registration. Useful for debugging, for example if the correlator filters its input.

  • corrs – The correlation maps, as returned from the correlator. Useful for debugging.

libertem_holo.base.align.cross_correlate(src, target, plot: bool = False, plot_title: str = '', normalization: ~typing.Literal['phase'] | None = 'phase', upsample_factor=1, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) tuple[ndarray, ndarray][source]

Rigid image registration by cross-correlation.

Supports optional phase normalization. Based on the phase_cross_correlation function of scikit-image, but with added GPU support via cupy, and some debugging facilities built in.

Parameters:
  • src – The static image, either a numpy or cupy array (if cupy, you should set xp=cp, too)

  • target – The moving image, either a numpy or cupy array (if cupy, you should set xp=cp, too)

  • normalization – ‘phase’ or None, same as for phase_cross_correlation

  • upsample_factor – Subpixel scaling factor, same as for phase_cross_correlation. Note that we test with a factor of up to 10; with larger values the results may not improve.

  • xp – numpy or cupy

libertem_holo.base.align.get_grad_angle(image, scale=3)[source]

From an image, get the angle of the gradient.

libertem_holo.base.align.get_grad_xy(image, scale=3)[source]

From an image, get the angle of the gradient.

libertem_holo.base.align.is_left(a: ndarray, b: ndarray, c: ndarray) ndarray[source]

Points a and b are points on a line and result is an array of True or False if c is left or right of line resp.

libertem_holo.base.align.stack_alignment_quality(wave_stack: ndarray, shifts)[source]

Stack quality ,judged by standard deviation on the stacking axis.

This should be mostly noise, if not, there may be issues from the alignment.

Simulation

libertem_holo.base.generate.hologram_frame(amp, phi, counts=1000.0, sampling=5.0, visibility=1.0, f_angle=30.0, gaussian_noise=None, poisson_noise=None, xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>)[source]

Generates holograms using phase and amplitude as an input

See Usage for detailed application example

New in version 0.3.0.

Notes

Theoretical basis for hologram simulations see in: Lichte, H., and M. Lehmann. Rep. Prog. Phys. 71 (2008): 016102. doi:10.1088/0034-4885/71/1/016102 [LL08]

Parameters:
  • amp (np.ndarray, 2d) – normalized amplitude and phase images of the same shape

  • phi (np.ndarray, 2d) – normalized amplitude and phase images of the same shape

  • counts (float, default: 1000.) – Number of electron counts in vacuum

  • sampling (float, default: 5.) – Hologram fringe sampling (number of pixels per fringe)

  • visibility (float, default: 1.) – Hologram fringe visibility (aka fringe contrast)

  • f_angle (float, default: 30.) – Angle in degrees of hologram fringes with respect to X-axis

  • gaussian_noise (float or int or None, default: None.) – Amount of Gaussian smoothing determined by sigma parameter applied to the hologram simulating effect of focus spread or PSF of the detector.

  • poisson_noise (float or int or None, default: None.) – Amount of Poisson applied to the hologram.

Returns:

holo – hologram image

Return type:

np.ndarray, 2d

Utility functions

Utility functions for working with holography data.

class libertem_holo.base.utils.HoloParams(sb_size: tuple[float, float], sb_position: tuple[float, float], aperture: np.ndarray, orig_shape: tuple[int, int], out_shape: tuple[int, int], scale_factor: float, xp: XPType)[source]

HoloParams class contians all parameters necessary for reconstruction.

aperture: ndarray

Alias for field number 2

classmethod from_hologram(hologram: ~numpy.ndarray, *, central_band_mask_radius: float | None = None, out_shape: tuple = None, sb_size: float | None = None, sb_position: tuple[float, float] | None = None, circle_filter_order: int = 20, line_filter_length: float = 0.9, line_filter_width: float | None = 3, line_filter_order: int = 2, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) HoloParams[source]

Determine reconstruction parameters from a hologram.

Automatically estimates sideband position and size, and returns the main parameters needed for holography reconstruction.

Parameters:
  • hologram – A single hologram, can be either a numpy or a cupy array

  • central_band_mask_radius – When estimating the sideband position, use a mask of this size to remove the central band

  • out_shape – The reconstruction shape, should be larger than the sideband size

  • sb_size – Override the sideband size; determined automatically if not given

  • sb_position – Override the sideband position; determined automatically if not given

  • circle_filter_order – Order of the butterworth filter applied to the circular part

  • line_filter_length – Length ratio of the line filter; as a fraction of the distance between the central band and the sideband

  • line_filter_width – Width of the line filter, in pixels. Passing in None will disable the line filter completely

  • line_filter_order – Order of the butterworth filter applied to the line part

  • xp – Pass in either the numpy or cupy module to select CPU or GPU processing

orig_shape: tuple[int, int]

Alias for field number 3

out_shape: tuple[int, int]

Alias for field number 4

sb_position: tuple[float, float]

Alias for field number 1

property sb_position_int: tuple[int, int]

Sideband position from float to int.

sb_size: tuple[float, float]

Alias for field number 0

scale_factor: float

Alias for field number 5

xp: Any

Alias for field number 6

libertem_holo.base.utils.estimate_sideband_position(holo_data: ~numpy.ndarray, holo_sampling: tuple[float, float], central_band_mask_radius: float | None = None, sb: ~typing.Literal['lower', 'upper'] = 'lower', xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) tuple[float, float][source]

Find the position of the sideband and return its position.

Parameters:
  • holo_data (ndarray) – The data of the hologram.

  • holo_sampling (tuple) – The sampling rate in both image directions.

  • central_band_mask_radius (float, optional) – The aperture radius used to mask out the centerband.

  • sb (str, optional) – Chooses which sideband is taken. ‘lower’ or ‘upper’

  • xp – Pass in either the numpy or cupy module to select CPU or GPU processing

Return type:

Tuple of the sideband position (y, x), referred to the unshifted FFT.

libertem_holo.base.utils.estimate_sideband_size(sb_position: tuple[float, float], holo_shape: tuple[int, int], sb_size_ratio: float = 0.5, xp: ~typing.Any = <module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) float[source]

Estimate the size of sideband filter.

Parameters:
  • holo_shape (array_like) – Holographic data array

  • sb_position (tuple) – The sideband position (y, x), referred to the non-shifted FFT.

  • sb_size_ratio (float, optional) – Size of sideband as a fraction of the distance to central band

  • xp – Pass in either the numpy or cupy module to select CPU or GPU processing

Returns:

sb_size – Size of sideband filter

Return type:

float

libertem_holo.base.utils.freq_array(shape: tuple[int, int], sampling: tuple[float, float] = (1.0, 1.0), xp=<module 'numpy' from '/home/runner/work/LiberTEM-holo/LiberTEM-holo/.tox/docs-build-ci/lib/python3.11/site-packages/numpy/__init__.py'>) ndarray[source]

Generate a frequency array.

Parameters:
  • shape ((int, int)) – The shape of the array.

  • sampling ((float, float), optional, (Default: (1., 1.))) – The sampling rates of the array.

Return type:

Array of the frequencies.

libertem_holo.base.utils.get_slice_fft(out_shape: tuple[int, int], sig_shape: tuple[int, int]) tuple[slice, slice][source]

Get a slice in fourier space to achieve the given output shape.

libertem_holo.base.utils.other_sb(sb_position, shape)[source]

Given the sb_position (as from the estimate function in fft coordinates), and the shape of the hologram, calculate the position of the other sideband position (also in fft coordinates)

libertem_holo.base.utils.remove_phase_ramp(img: ndarray, *, roi=None, method: Literal['gradient'] | Literal['fit'] = 'gradient') tuple[ndarray, ndarray, tuple[float, float]][source]

Remove a phase ramp from img.

Returns both the compensated image and the ramp that was removed.

Parameters:
  • img – The (phase) input image, has to be already unwrapped

  • roi – Either a slice (as returned by np.s_ for example), or an array of the region of interest in img. If not specified, the whole img is used. In any case, the ramp is subtraced from the whole image.

  • method

    • ‘gradient’: the average gradient in the specified region of interest

    • ’fit’: a least-square fit of a linear gradient to the data in the region of interest

Returns:

  • corrected_img – The image with the ramp removed

  • ramp – The ramp that was found, as a 2D gradient

  • (ramp_y, ramp_x) – The ramp slopes