corrct.processing.post

Post-processing routines.

Created on Tue Mar 24 15:25:14 2020

@author: Nicola VIGANÒ, Computational Imaging group, CWI, The Netherlands, and ESRF - The European Synchrotron, Grenoble, France

Module Contents

Functions

com

Compute center-of-mass for given volume.

power_spectrum

Compute the power spectrum of a n-dimensional signal.

compute_frc

Compute the FRC/FSC (Fourier ring/shell correlation) between two images / volumes.

estimate_resolution

Estimate the resolution or bandwidth, given an FRC and a threshold curve.

plot_frcs

Compute and plot the FSCs / FRCs of some volumes.

fit_scale_bias

Fit the scale and bias of an image, against its projection in a different space.

Data

eps

API

corrct.processing.post.eps

None

corrct.processing.post.com(vol: numpy.typing.NDArray, axes: Optional[numpy.typing.ArrayLike] = None) numpy.typing.NDArray[source]

Compute center-of-mass for given volume.

Parameters

vol : NDArray The input volume. axes : ArrayLike, optional Axes on which to compute center-of-mass. The default is None.

Returns

NDArray The center-of-mass position.

corrct.processing.post.power_spectrum(img: numpy.typing.NDArray, axes: Optional[collections.abc.Sequence[int]] = None, smooth: Optional[int] = 5, taper_ratio: Optional[float] = 0.05, power: int = 2) numpy.typing.NDArray[source]

Compute the power spectrum of a n-dimensional signal.

Parameters

img : NDArray The n-dimensional signal. axes : Optional[Sequence[int]], optional The axes over which we want to compute the power spectrum, by default None smooth : Optional[int], optional The smoothing kernel size, by default 5 taper_ratio : Optional[float], optional Whether to taper the signal at the edges (for truncated signals), by default 0.05 power : int, optional The exponent to use, by default 2

Returns

NDArray The power spectrum

corrct.processing.post.compute_frc(img1: numpy.typing.NDArray, img2: Optional[numpy.typing.NDArray], snrt: float = 0.2071, axes: Optional[collections.abc.Sequence[int]] = None, smooth: Optional[int] = 5, taper_ratio: Optional[float] = 0.05, supersampling: int = 1, theo_threshold: bool = True) tuple[numpy.typing.NDArray, numpy.typing.NDArray][source]

Compute the FRC/FSC (Fourier ring/shell correlation) between two images / volumes.

Please refer to the following article for more information: M. van Heel and M. Schatz, “Fourier shell correlation threshold criteria,” J. Struct. Biol., vol. 151, no. 3, pp. 250–262, Sep. 2005.

Parameters

img1 : NDArray First image / volume. img2 : NDArray Second image / volume. snrt : float, optional SNR to be used for generating the threshold curve for resolution definition. The SNR value of 0.4142 corresponds to the hald-bit curve for a full dataset. When splitting datasets in two sub-datasets, that value needs to be halved. The default is 0.2071, which corresponds to the half-bit threashold for half dataset. axes : Sequence[int], optional The axes over which we want to compute the FRC/FSC. If None, all axes will be used The default is None. smooth : Optional[int], optional Size of the Hann smoothing window. The default is 5. taper_ratio : Optional[float], optional Ratio of the edge pixels to be tapered off. This is necessary when working with truncated volumes / local tomography, to avoid truncation artifacts. The default is 0.05. supersampling : int, optional Supersampling factor of the images. Larger values increase the high-frequency range of the FRC/FSC function. The default is 1, which corresponds to the Nyquist frequency.

Raises

ValueError Error returned when not passing images of the same shape.

Returns

NDArray The computed FRC/FSC. NDArray The threshold curve corresponding to the given threshod SNR.

corrct.processing.post.estimate_resolution(frc: numpy.typing.NDArray, t_hb: numpy.typing.NDArray) Optional[tuple[float, float]][source]

Estimate the resolution or bandwidth, given an FRC and a threshold curve.

Parameters

frc : NDArray The FRC curve t_hb : NDArray The threshold curve

Returns

tuple[float, float] | None The resolution or bandwidth, if a crossing point was found. Otherwise None.

corrct.processing.post.plot_frcs(volume_pairs: collections.abc.Sequence[tuple[numpy.typing.NDArray, numpy.typing.NDArray]], labels: collections.abc.Sequence[str], title: Optional[str] = None, smooth: Optional[int] = 5, snrt: float = 0.2071, axes: Optional[collections.abc.Sequence[int]] = None, supersampling: int = 1, verbose: bool = False) tuple[matplotlib.figure.Figure, matplotlib.axes._axes.Axes][source]

Compute and plot the FSCs / FRCs of some volumes.

Parameters

volume_pairs : Sequence[Tuple[NDArray, NDArray]] A list of pairs of volumes to compute the FRCs on. labels : Sequence[str] The labels associated with each pair. title : Optional[str], optional The axes title, by default None. smooth : Optional[int], optional The size of the smoothing window for the computed curves, by default 5. snrt : float, optional The SNR of the T curve, by default 0.2071 - as per half-dataset SNR. axes : Sequence[int] | None, optional The axes along which we want to compute the FRC. The unused axes will be averaged. The default is None. verbose : bool, optional Whether to display verbose output, by default False.

corrct.processing.post.fit_scale_bias(img_data: numpy.typing.NDArray, prj_data: numpy.typing.NDArray, prj: Optional[corrct.operators.BaseTransform] = None) tuple[float, float][source]

Fit the scale and bias of an image, against its projection in a different space.

Parameters

img_data : NDArray The image data prj_data : NDArray The projected data prj : BaseTransform | None, optional The projection operator. The default is None, which uses the identity (TransformIdentity)

Returns

tuple[float, float] The scale and bias