Code reference

Filtering / Sideband isolation

qpretrieve.filter.get_filter_array(filter_name: str, filter_size: float, freq_pos: tuple[float, float], fft_shape: tuple[int, int]) ndarray[source]

Create a Fourier filter for holography

Parameters:
  • filter_name (str) –

    specifies the filter to use, one of

    • ”disk”: binary disk with radius filter_size

    • ”smooth disk”: disk with radius filter_size convolved with a radial gaussian (sigma=filter_size/5)

    • ”gauss”: radial gaussian (sigma=0.6*filter_size)

    • ”square”: binary square with side length 2*filter_size

    • ”smooth square”: square with side length 2*filter_size convolved with square gaussian (sigma=filter_size/5)

    • ”tukey”: a square tukey window of width 2*filter_size and alpha=0.1

  • filter_size (float) –

    Size of the filter in Fourier space.

    In this function, filter_size is interpreted in the same normalized frequency units as returned by numpy.fft.fftfreq() (and typically used together with numpy.fft.fftshift()), i.e. values are in the range [-0.5, 0.5) along each axis.

    Implications for each filter_name: - For the radial/square filters (e.g. “disk”, “gauss”, “square”, …),

    values up to about 0.5 correspond to radii within the Fourier domain. Values larger than 0.5 (e.g. 1) effectively yield an all-pass filter (no filtering).

    • The “tukey” filter is different and that scales with the array size.

    For “tukey”, filter_size=1 can exceed the array bounds.

  • freq_pos (tuple of floats) – The position of the filter in frequency coordinates as returned by numpy.fft.fftfreq().

  • fft_shape (tuple of int) – The shape of the Fourier transformed image (2d) for which the filter will be applied. The shape must be squared (two identical integers).

Returns:

filt_arr – The Fourier-shifted filtering array. For mask images, this is a boolean array. For more elaborate filters, this is a float array.

Return type:

2d ndarray, boolean of float

Fourier transform methods

Numpy

class qpretrieve.fourier.ff_numpy.FFTFilterNumpy(data: ndarray, subtract_mean: bool = True, padding: int = 2, copy: bool = True, dtype_conversion=None)[source]

Wraps the numpy Fourier transform

Parameters:
  • data

    The experimental input real-valued image. Allowed input shapes are:
    • 2d (y, x)

    • 3d (z, y, x)

    • rgb (y, x, 3) or rgba (y, x, 4)

  • subtract_mean (bool) – If True, subtract the mean of data before performing the Fourier transform. This setting is recommended as it can reduce artifacts from frequencies around the central band.

  • padding (int) –

    Boundary padding with zeros; This value determines how large the padding region should be. If set to zero, then no padding is performed. If set to a positive integer, the size is computed to the next power of two (square image):

    2 ** xp.ceil(xp.log(padding * max(data.shape) / xp.log(2)))
    

  • copy (bool) – If set to True, make sure that data is not edited. If you set this to False, then caching FFT results will not work anymore.

  • dtype_conversion – The dtype that should be used to convert the input data before preprocessing occurs. This defaults to complex if the input data is complex, otherwise to float (64-bit) for all other situations. For some use-cases, for example when using a GPU, you might want to be more specific e.g., cp.float32.

Notes

The initial Fourier transform of the input data is cached using weak references (only if copy is True).

backend_check()

Warn if the FFTFilter superclass doesn’t match expected backend. Raise an Error if the FFTFilterPyFFTW is used with numpy backend.

Added in version 0.6.1.

filter(filter_name: str, filter_size: float, freq_pos: float, float, scale_to_filter: bool | float = False) xp.ndarray
Parameters:
  • filter_name (str) –

    specifies the filter to use, one of

    • ”disk”: binary disk with radius filter_size

    • ”smooth disk”: disk with radius filter_size convolved with a radial gaussian (sigma=filter_size/5)

    • ”gauss”: radial gaussian (sigma=0.6*filter_size)

    • ”square”: binary square with side length 2*filter_size

    • ”smooth square”: square with side length 2*filter_size convolved with square gaussian (sigma=filter_size/5)

    • ”tukey”: a square tukey window of width 2*filter_size and alpha=0.1

  • filter_size (float) – Size of the filter in Fourier space. The filter size interpreted as a Fourier frequency index (“pixel size”) and must be between 0 and max(fft_shape)/2

  • freq_pos (tuple of floats) – The position of the filter in frequency coordinates as returned by numpy.fft.fftfreq().

  • scale_to_filter (bool or float) – Crop the image in Fourier space after applying the filter, effectively removing surplus (zero-padding) data and increasing the pixel size in the output image. If True is given, then the cropped area is defined by the filter size, if a float is given, the cropped area is defined by the filter size multiplied by scale_to_filter. You can safely set this to True for filters with a binary support. For filters such as “smooth square” or “gauss” (filter is not a boolean array but a floating-point array), the higher you set scale_to_filter, the more information will be included in the scaled image.

Notes

The FFT result is cached using weak references in FFTCache. If you call this function a lot of times with different arguments, then it might look like a memory leak. However, you just have to delete the FFTFilter instance and everything will get garbage-collected.

property shape: tuple

Shape of the Fourier transform data

origin

original data (with subtracted mean)

padding

whether padding is enabled

subtract_mean

whether the mean was subtracted

origin_padded

padded input data

fft_origin

frequency-shifted Fourier transform

fft_filtered

filtered Fourier transform

fft_used

used Fourier transform (can have a different shape)

PyFFTW

Cupy

Interference image analysis

Off-Axis Holography (DHM)

class qpretrieve.interfere.if_oah.OffAxisHologram(data: ndarray, fft_interface: str | Type[FFTFilter] = 'auto', subtract_mean=True, padding=2, copy=True, dtype_conversion=None, **pipeline_kws)[source]

Off-axis hologram data analysis

Parameters:
  • data

    The experimental input real-valued image. Allowed input shapes are:
    • 2d (y, x)

    • 3d (z, y, x)

    • rgb (y, x, 3) or rgba (y, x, 4)

  • fft_interface – A Fourier transform interface. See qpretrieve.fourier.get_available_interfaces() to get a list of implemented interfaces. Default is “auto”, which will use qpretrieve.fourier.get_best_interface(). This is in line with old behaviour. See Notes for more details.

  • subtract_mean (bool) – If True, remove the mean of the hologram before performing the Fourier transform. This setting is recommended as it can reduce artifacts from frequencies around the central band.

  • padding (bool) –

    Boundary padding with zeros; This value determines how large the padding region should be. If set to zero, then no padding is performed. If set to a positive integer, the size is computed to the next power of two (square image):

    2 ** xp.ceil(xp.log(padding * max(data.shape) / xp.log(2)))
    

  • copy (bool) – Whether to create a copy of the input data.

  • dtype_conversion – The dtype that should be used to convert the input data before preprocessing occurs. This defaults to complex if the input data is complex, otherwise to float (64-bit) for all other situations. For some use-cases, for example when using a GPU, you might want to be more specific e.g., cp.float32.

  • pipeline_kws – Any additional keyword arguments for run_pipeline() as defined in default_pipeline_kws.

Notes

For fft_interface, if you do not have the relevant package installed, then an error will be raised. For example, setting fft_interface=FFTFilterPyFFTW will fail if you do not have pyfftw installed.

default_pipeline_kws = {'filter_name': 'disk', 'filter_size': 0.3333333333333333, 'filter_size_interpretation': 'sideband distance', 'invert_phase': False, 'scale_to_filter': False, 'sideband_freq': None}

Default OAH pipeline keyword arguments

property phase: ndarray

Retrieved phase information

property amplitude: ndarray

Retrieved amplitude information

run_pipeline(**pipeline_kws) ndarray[source]

Run OAH analysis pipeline

Parameters:
  • filter_name (str) – specifies the filter to use, see qpretrieve.filter.get_filter_array().

  • filter_size (float) – Size of the filter in Fourier space. The interpretation of this value depends on filter_size_interpretation.

  • filter_size_interpretation (str) – If set to “sideband distance”, the filter size is interpreted as the relative distance between central band and sideband (this is the default). If set to “frequency index”, the filter size is interpreted as a Fourier frequency index (“pixel size”) and must be between 0 and max(hologram.shape)/2. If set to “physical radius”, the base radius is \(r \approx n dx NA / \lambda\) in Fourier pixels. In this mode filter_size is a scaling factor (1.0 means direct physical radius).

  • scale_to_filter (bool or float) – Crop the image in Fourier space after applying the filter, effectively removing surplus (zero-padding) data and increasing the pixel size in the output image. If True is given, then the cropped area is defined by the filter size, if a float is given, the cropped area is defined by the filter size multiplied by scale_to_filter. You can safely set this to True for filters with a binary support. For filters such as “smooth square” or “gauss” (filter is not a boolean array but a floating-point array), the higher you set scale_to_filter, the more information will be included in the scaled image.

  • sideband_freq (tuple of floats) – Frequency coordinates of the sideband to use. By default, a heuristic search for the sideband is done. If you pass a 3D array, the first hologram is used to determine the sideband frequencies.

  • pixel_size (float) – Sensor pixel size dx in meters for physical-radius mode.

  • numerical_aperture (float) – Collection NA for physical-radius mode.

  • wavelength (float) – Illumination wavelength in meters for physical-radius mode.

  • invert_phase (bool) – Invert the phase data.

compute_filter_size(filter_size: float, filter_size_interpretation: str, sideband_freq: tuple[float, float] = None, pixel_size: float | None = None, numerical_aperture: float | None = None, wavelength: float | None = None) float

Compute the actual filter size in Fourier space

property field: ndarray

Retrieved amplitude information

get_data_with_input_layout(data: ndarray | str) ndarray

Convert data to the original input array layout.

Parameters:

data – Either an array (xp.ndarray) or name (str) of the relevant data.

Returns:

data_out – array in the original input array layout

Return type:

xp.ndarray

Notes

If data is the RGBA array layout, then the alpha (A) channel will be an array of ones.

get_pipeline_kw(key)

Current pipeline keyword argument with fallback to defaults

process_like(other)

Process this dataset in the same way as other dataset

fft

qpretrieve Fourier transform interface class

fft_origin

originally computed Fourier transform

fft_filtered

filtered Fourier data from last run of run_pipeline

pipeline_kws

pipeline keyword arguments

qpretrieve.interfere.if_oah.find_peak_cosine(ft_data: ndarray, copy: bool = True) tuple[float, float][source]

Find the side band position of a 2d regular off-axis hologram

The Fourier transform of a cosine function (known as the striped fringe pattern in off-axis holography) results in two sidebands in Fourier space.

The hologram is Fourier-transformed and the side band is determined by finding the maximum amplitude in Fourier space.

Parameters:
  • ft_data (2d ndarray) – FFt-shifted Fourier transform of the hologram image

  • copy (bool) – copy ft_data before modification

Returns:

fsx, fsy – coordinates of the side band in Fourier space frequencies

Return type:

tuple of floats

Quadriwave lateral shearing interferometry (QLSI)

class qpretrieve.interfere.if_qlsi.QLSInterferogram(data, reference=None, *args, **kwargs)[source]

Interferometric analysis of quadri-wave lateral shearing holograms

Parameters:
  • data

    The experimental input real-valued image. Allowed input shapes are:
    • 2d (y, x)

    • 3d (z, y, x)

    • rgb (y, x, 3) or rgba (y, x, 4)

  • fft_interface – A Fourier transform interface. See qpretrieve.fourier.get_available_interfaces() to get a list of implemented interfaces. Default is “auto”, which will use qpretrieve.fourier.get_best_interface(). This is in line with old behaviour. See Notes for more details.

  • subtract_mean (bool) – If True, remove the mean of the hologram before performing the Fourier transform. This setting is recommended as it can reduce artifacts from frequencies around the central band.

  • padding (bool) –

    Boundary padding with zeros; This value determines how large the padding region should be. If set to zero, then no padding is performed. If set to a positive integer, the size is computed to the next power of two (square image):

    2 ** xp.ceil(xp.log(padding * max(data.shape) / xp.log(2)))
    

  • copy (bool) – Whether to create a copy of the input data.

  • dtype_conversion – The dtype that should be used to convert the input data before preprocessing occurs. This defaults to complex if the input data is complex, otherwise to float (64-bit) for all other situations. For some use-cases, for example when using a GPU, you might want to be more specific e.g., cp.float32.

  • pipeline_kws – Any additional keyword arguments for run_pipeline() as defined in default_pipeline_kws.

Notes

For fft_interface, if you do not have the relevant package installed, then an error will be raised. For example, setting fft_interface=FFTFilterPyFFTW will fail if you do not have pyfftw installed.

default_pipeline_kws = {'filter_name': 'square', 'filter_size': 0.5, 'filter_size_interpretation': 'sideband distance', 'invert_phase': False, 'qlsi_pitch_term': None, 'scale_to_filter': False, 'sideband_freq': None, 'wavelength': None}

Default QLSI pipeline keyword arguments

property amplitude: ndarray

Retrieved amplitude information

property field: ndarray

Retrieved amplitude information

property phase: ndarray

Retrieved phase information

run_pipeline(**pipeline_kws) ndarray[source]

Run QLSI analysis pipeline

Parameters:
  • filter_name (str) – specifies the filter to use, see qpretrieve.filter.get_filter_array().

  • filter_size (float) – Size of the filter in Fourier space. The interpretation of this value depends on filter_size_interpretation.

  • filter_size_interpretation (str) – If set to “sideband distance”, the filter size is interpreted as the relative distance between central band and sideband (this is the default). If set to “frequency index”, the filter size is interpreted as a Fourier frequency index (“pixel size”) and must be between 0 and max(hologram.shape)/2. If set to “physical radius”, the base radius is \(r \approx n dx NA / \lambda\) in Fourier pixels. In this mode filter_size is a scaling factor (1.0 means direct physical radius).

  • scale_to_filter (bool or float) – Crop the image in Fourier space after applying the filter, effectively removing surplus (zero-padding) data and increasing the pixel size in the output image. If True is given, then the cropped area is defined by the filter size, if a float is given, the cropped area is defined by the filter size multiplied by scale_to_filter. You can safely set this to True for filters with a binary support. For filters such as “smooth square” or “gauss” (filter is not a boolean array but a floating-point array), the higher you set scale_to_filter, the more information will be included in the scaled image.

  • sideband_freq (tuple of floats) – Frequency coordinates of the sideband to use. By default, a heuristic search for the sideband is done. If you pass a 3D array, the first hologram is used to determine the sideband frequencies.

  • pixel_size (float) – Sensor pixel size dx in meters for physical-radius mode.

  • numerical_aperture (float) – Collection NA for physical-radius mode.

  • invert_phase (bool) – Invert the phase data.

  • wavelength (float) – Wavelength to convert from the wavefront in meters to radians.

  • qlsi_pitch_term (float) –

    Scaling term converting the integrated gradient image to the unit meters. This term is computed from the lattice constant of the grating \(L\), the distance between the grating and the camera sensor \(d\) and the physical camera pixel width \(a\) according to

    \[\text{pitch_term} = \frac{La}{d}\]

    For the case where the lattice constant is four times the pixel width, this simplifies to \(4a^2/d\). Note that for a relay-lens system (grating not directly attached to the sensor) this factor is wavelength dependent due to chromatic aberrations introduced by the lenses. For gratings-on-a-camera configurations (e.g. Phasics SID4Bio), this is a device-specific quantity which has to be determined only once. E.g. for our SID4Bio camera, this value is 0.01887711 µm (1.87711e-08 m).

compute_filter_size(filter_size: float, filter_size_interpretation: str, sideband_freq: tuple[float, float] = None, pixel_size: float | None = None, numerical_aperture: float | None = None, wavelength: float | None = None) float

Compute the actual filter size in Fourier space

get_data_with_input_layout(data: ndarray | str) ndarray

Convert data to the original input array layout.

Parameters:

data – Either an array (xp.ndarray) or name (str) of the relevant data.

Returns:

data_out – array in the original input array layout

Return type:

xp.ndarray

Notes

If data is the RGBA array layout, then the alpha (A) channel will be an array of ones.

get_pipeline_kw(key)

Current pipeline keyword argument with fallback to defaults

process_like(other)

Process this dataset in the same way as other dataset

fft

qpretrieve Fourier transform interface class

fft_origin

originally computed Fourier transform

fft_filtered

filtered Fourier data from last run of run_pipeline

pipeline_kws

pipeline keyword arguments

qpretrieve.interfere.if_qlsi.find_peaks_qlsi(ft_data: ndarray, periodicity: int = 4, copy: bool = True) tuple[tuple[float, float], tuple[float, float]][source]

Find the two peaks in Fourier space for the x and y gradient

Parameters:
  • ft_data (2d complex ndarray) – FFT-shifted Fourier transform of the QLSI image

  • periodicity (float) – Grid size of the QLSI image. For the Phasics SID4Bio camera, this is 4 (i.e. the peak-to-peak distance of the individual foci in the QLSI image is four pixels)

  • copy (bool) – Set to False to perform operations in-place.

Returns:

  • (f1x, f1y) (tuple of floats) – Coordinates of the first gradient peak in frequency coordinates.

  • (f2x, f2y) (tuple of floats) – Coordinates of the second gradient peak in frequency coordinates.

Notes

At some point it might be necessary to add an angle keyword argument that gives the algorithm a hint about te rotation of the QLSI grid. Currently, peak detection is only done in the lower half of ft_data. If the peaks are exactly aligned with the pixel grid, then the current approach might not work. Also, setting angle=xp.pi would be equivalent to setting sideband to -1 in holo.py (would be a nice feature).

Data Array Layout

Module that provides convenience functions for converting data between array layouts.

Added in version 0.4.0.

qpretrieve.data_array_layout.convert_data_to_3d_array_layout(data: ndarray) tuple[ndarray, str][source]

Convert the data to the 3d array_layout

Returns:

  • data_out – 3d version of the data

  • array_layout – original array layout for future reference

Notes

If input is either a RGB or RGBA array layout as input, the first channel is taken as the image to process. In other words, it is assumed that all channels contain the same information, so the first channel is used. 3D RGB/RGBA array layouts, such as (50, 256, 256, 3), are not allowed (yet).

qpretrieve.data_array_layout.convert_3d_data_to_array_layout(data: ndarray, array_layout: str) ndarray[source]

Convert the 3d data to the desired array_layout.

Returns:

data_out – input data with the given array layout

Return type:

xp.ndarray

Notes

Currently, this function is limited to converting from 3d to other array layouts. Perhaps if there is demand in the future, this can be generalised for other conversions.

NDArray Backend

Module that controls and exposes the active ndarray backend (NumPy or CuPy).

Added in version 0.6.0.

class qpretrieve._ndarray_backend.NDArrayBackend[source]

Proxy object exposing the current ndarray backend.

get()[source]

Return the currently active backend module.

set(backend_name: str = 'numpy')[source]

Switch the backend between ‘numpy’ and ‘cupy’.

exception qpretrieve._ndarray_backend.NDArrayBackendWarning(message)[source]
add_note()

Exception.add_note(note) – add a note to the exception

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

qpretrieve._ndarray_backend.get_ndarray_backend()

Return the currently active backend module.

qpretrieve._ndarray_backend.set_ndarray_backend(backend_name: str = 'numpy')

Switch the backend between ‘numpy’ and ‘cupy’.