Skip to content

nabu.preproc.phase

source module nabu.preproc.phase

Classes

  • PaganinPhaseRetrieval Paganin Phase Retrieval for an infinitely distant point source. Formula (10) in [1].

Functions

  • lmicron_to_db Utility to convert the "Lmicron" parameter of PyHST to a value of delta/beta.

  • compute_paganin_margin Compute the convolution margin to use when calling PaganinPhaseRetrieval class.

source lmicron_to_db(Lmicron, energy, distance)

Utility to convert the "Lmicron" parameter of PyHST to a value of delta/beta.

Parameters

  • Lmicron : float Length in microns, values of the parameter "PAGANIN_Lmicron" in PyHST2 parameter file.

  • energy : float Energy in keV.

  • distance : float Sample-detector distance in microns

Notes

The conversion is done using the formula

\[ L^2 = \pi \lambda D \frac{\delta}{\beta} \]

source class PaganinPhaseRetrieval(shape, distance=0.5, energy=20, delta_beta=250.0, pixel_size=1e-06, padding='edge', use_rfft=True, use_R2C=None, fftw_num_threads=None, fft_num_threads=None)

Paganin Phase Retrieval for an infinitely distant point source. Formula (10) in [1].

Parameters

  • shape : int or tuple Shape of each radio, in the format (num_rows, num_columns), i.e (size_vertical, size_horizontal). If an integer is provided, the shape is assumed to be square.

  • distance : float, optional Propagation distance in meters.

  • energy : float, optional Energy in keV.

  • delta_beta : float, optional delta/beta ratio, where n = (1 - delta) + i*beta is the complex refractive index of the sample.

  • pixel_size : float or tuple, optional Detector pixel size in meters. Default is 1e-6 (one micron) If a tuple is passed, the pixel size is set as (horizontal_size, vertical_size).

  • padding : str, optional Padding method. Available are "zeros", "mean", "edge", "sym", "reflect". Default is "edge". Please refer to the "Padding" section below for more details.

  • use_rfft : bool, optional Whether to use Real-to-Complex (R2C) transform instead of standard Complex-to-Complex transform, providing better performances

  • use_R2C : bool, optional DEPRECATED, use use_rfft instead

  • fftw_num_threads : bool or None or int, optional DEPRECATED - please use fft_num_threads

  • fft_num_threads : bool or None or int, optional Number of threads for FFT. Default is to use all available threads. You can pass a negative number to use N - fft_num_threads cores.

Important

Mind the units! Distance and pixel size are in meters, and energy is in keV.

References

[1] D. Paganin Et Al, "Simultaneous phase and amplitude extraction from a single defocused image of a homogeneous object", Journal of Microscopy, Vol 206, Part 1, 2002

Notes

Padding methods

The phase retrieval is a convolution done in Fourier domain using FFT, so the Fourier transform size has to be at least twice the size of the original data. Mathematically, the data should be padded with zeros before being Fourier transformed. However, in practice, this can lead to artefacts at the edges (Gibbs effect) if the data does not go to zero at the edges. Apart from applying an apodization (Hamming, Blackman, etc), a common strategy to avoid these artefacts is to pad the data. In tomography reconstruction, this is usually done by replicating the last(s) value(s) of the edges ; but one can think of other methods:

  • "zeros": the data is simply padded with zeros.
  • "mean": the upper side of extended data is padded with the mean of the first row, the lower side with the mean of the last row, etc.
  • "edge": the data is padded by replicating the edges. This is the default mode.
  • "sym": the data is padded by mirroring the data with respect to its edges. See numpy.pad().
  • "reflect": the data is padded by reflecting the data with respect to its edges, including the edges. See numpy.pad().

Formulas

The radio is divided, in the Fourier domain, by the original "Paganin filter" [1].

\[ F = 1 + \frac{\delta}{\beta} \lambda D \\i |k|^2 \]

where k is the wave vector.

Methods

source method PaganinPhaseRetrieval.compute_filter()

source method PaganinPhaseRetrieval.pad_with_values(data, top_val=0, bottom_val=0, left_val=0, right_val=0)

Pad the data into self.padded_data with values.

Parameters

  • data : numpy.ndarray data (radio)

  • top_val : float or numpy.ndarray, optional Value(s) to fill the top of the padded data with.

  • bottom_val : float or numpy.ndarray, optional Value(s) to fill the bottom of the padded data with.

  • left_val : float or numpy.ndarray, optional Value(s) to fill the left of the padded data with.

  • right_val : float or numpy.ndarray, optional Value(s) to fill the right of the padded data with.

source method PaganinPhaseRetrieval.pad_data(data, padding_method=None)

Raises

  • ValueError

source method PaganinPhaseRetrieval.apply_filter(radio, padding_method=None, output=None)

source method PaganinPhaseRetrieval.lmicron_to_db(Lmicron)

Utility to convert the "Lmicron" parameter of PyHST to a value of delta/beta. Please see the doc of nabu.preproc.phase.lmicron_to_db()

source compute_paganin_margin(shape, cutoff=1000.0, **pag_kwargs)

Compute the convolution margin to use when calling PaganinPhaseRetrieval class.

Parameters

  • shape : tuple Detector shape in the form (n_z, n_x)