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

>>> shape = tuple(dataset.shape.sig)
>>> sb_position = [2, 3]
>>> sb_size = 4.4
>>> aperture = disk_aperture(out_shape=shape, radius=sb_size)
>>> holo_udf = HoloReconstructUDF(
...     out_shape=shape,
...     sb_position=sb_position,
...     aperture=aperture,
... )
>>> wave = ctx.run_udf(dataset=dataset, udf=holo_udf)['wave'].data
classmethod with_default_aperture(*, out_shape: tuple[int, int], sb_size: float, sb_position: tuple[float, float], precision: bool = True) HoloReconstructUDF[source]

Instantiate with a default disk-shaped aperture.

Examples

>>> udf = HoloReconstructUDF.with_default_aperture(
...     out_shape=(128, 128),
...     sb_size=7.6,
...     sb_position=(32, 32),
... )

Base holography functions

Image filtering

Useful image filtering helpers.

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.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.phase_ramp_finding(img, order=1)[source]

Find a phase ramp in img.

A phase ramp finding function that is used to find the phase ramp across the field of view.

Parameters
  • img (2d nd array) – Complex image or phase image.

  • order (int) – Phase ramp, 1 (default) is linear.

  • ramp (2d tuple, ()) – Phase ramp in y, x, if not None.

Return type

ramp, order, tuple, float

libertem_holo.base.filters.phase_ramp_removal(size, order=1, ramp=None)[source]

Remove phase ramp.

A phase ramp removal function that finds and removes a phase ramp across the field of view.

Parameters
  • size (2d tuple, ()) – Size of the Complex image or phase image

  • order (int) – Phase ramp, 1 (default) is linear.

  • ramp (2d tuple, ()) – Phase ramp in y, x, if not None.

Return type

2d nd array of the corrected 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.ramp_compensation(image)[source]

A ramp or wedge compensation for a 2D image with a linear optimization methods.

Parameters

image (2D-Array) – Input array

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.

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)[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

Masks and apertures

Functions for building apertures.

libertem_holo.base.mask.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, fft-shifted.

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.mask.line_filter(orig_shape: tuple[int, int], sb_pos: tuple[int, int], width: int, length: int, slice_fft: slice) ndarray[source]

Remove Fresnel fringes from biprism with a line filter in fourier space.

The starting point is the sideband position. The end points depend on the length and in the direction to top right image.

Parameters
  • orig_shape – the shape of the image.

  • sb_pos – Position of the sideband that is used for reconstruction of holograms.

  • width – Width of the line (rectangle) in pixels.

  • length – Length of the line (rectangle) in pixels.

  • slice_fft – A slice in fft shifted coordinates

Return type

2d array containing line filter

Reconstruction

Core reconstruction functions of LiberTEM-holo.

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

libertem_holo.base.reconstr.display_fft_image(image, sb_position, slice_fft, mask=1, detail=True)[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.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.reconstr.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.reconstr.freq_array(shape: tuple[int, int], sampling: tuple[float, float] = (1.0, 1.0)) 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.reconstr.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.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