corrct.processing.misc
Miscellaneous processing routines.
@author: Nicola VIGANÒ, Computational Imaging group, CWI, The Netherlands, and ESRF - The European Synchrotron, Grenoble, France
Module Contents
Functions
Compute a circular mask for the reconstruction volume. |
|
Compute a ball with specified radius. |
|
Compute the azimuthal integration of a n-dimensional image or a stack of them. |
|
Compute the intersection point between two lines. |
|
Compute the normalized cross-correlation between two images. |
|
Display Fourier representation of the input image. |
Data
API
- corrct.processing.misc.eps
None
- corrct.processing.misc.NDArrayInt
None
- corrct.processing.misc.circular_mask(vol_shape_zxy: Union[collections.abc.Sequence[int], corrct.processing.misc.NDArrayInt], radius_offset: float = 0, coords_ball: Union[collections.abc.Sequence[int], corrct.processing.misc.NDArrayInt, None] = None, ball_norm: float = 2, vol_origin_zxy: Union[collections.abc.Sequence[float], numpy.typing.NDArray, None] = None, taper_func: Optional[str] = None, taper_target: Union[str, float] = 'edge', super_sampling: int = 1, squeeze: bool = True, dtype: numpy.typing.DTypeLike = np.float32) numpy.typing.NDArray [source]
Compute a circular mask for the reconstruction volume.
Parameters
vol_shape_zxy : Sequence[int] | NDArrayInt The size of the volume. radius_offset : float, optional The offset with respect to the volume edge. The default is 0. coords_ball : Sequence[int] | NDArrayInt | None, optional The coordinates to consider for the non-masked region. The default is None. ball_norm : float, optional The norm of the ball. The default is 2. vol_origin_zxy : Optional[Sequence[float]], optional The origin of the coordinates in voxels. The default is None. taper_func : str, optional The mask data type. Allowed types: “const” | “cos”. The default is “const”. taper_target : str | float, optional The size target. Allowed values: “edge” | “diagonal”. The default is “edge”. super_sampling : int, optional The pixel super sampling to be used for the mask. The default is 1. squeeze : bool, optional Whether to squeeze the mask. The default is True. dtype : DTypeLike, optional The type of mask. The default is np.float32.
Raises
ValueError In case of unknown taper_func value, or mismatching volume origin and shape.
Returns
NDArray The circular mask.
- corrct.processing.misc.ball(data_shape_vu: numpy.typing.ArrayLike, radius: Union[int, float], super_sampling: int = 5, dtype: numpy.typing.DTypeLike = np.float32, func: Optional[Callable] = None) numpy.typing.ArrayLike [source]
Compute a ball with specified radius.
Parameters
data_shape_vu : ArrayLike Shape of the output array. radius : int | float Radius of the ball. super_sampling : int, optional Super-sampling for having smoother ball edges. The default is 5. dtype : DTypeLike, optional Type of the output. The default is np.float32. func : Optional[Callable], optional Point-wise function for the local values. The default is None.
Returns
ArrayLike The ball.
- corrct.processing.misc.azimuthal_integration(img: numpy.typing.NDArray, axes: collections.abc.Sequence[int] = (-2, -1), domain: str = 'direct') numpy.typing.NDArray [source]
Compute the azimuthal integration of a n-dimensional image or a stack of them.
Parameters
img : NDArray The image or stack of images. axes : tuple(int, int), optional Axes of that need to be azimuthally integrated. The default is (-2, -1). domain : string, optional Domain of the integration. Options are: “direct” | “fourier”. Default is “direct”.
Raises
ValueError Error returned when not passing images or wrong axes.
Returns
NDArray The azimuthally integrated profile.
- corrct.processing.misc.lines_intersection(line_1: numpy.typing.NDArray, line_2: Union[float, numpy.typing.NDArray], position: str = 'first', x_lims: Optional[tuple[Optional[float], Optional[float]]] = None) Optional[tuple[float, float]] [source]
Compute the intersection point between two lines.
Parameters
line_1 : NDArray The first line. line_2 : float | NDArray The second line. It can be a scalar representing a horizontal line. position : str, optional The position of the point to select. Either “first” or “last”. The default is “first”.
Raises
ValueError If position is neither “first” nor “last”.
Returns
Tuple[float, float] | None It returns either the requested crossing point, or None in case the point was not found.
- corrct.processing.misc.norm_cross_corr(img1: numpy.typing.NDArray, img2: Optional[numpy.typing.NDArray] = None, axes: collections.abc.Sequence[int] = (-2, -1), t_match: bool = False, mode_full: bool = True, compute_profile: bool = True, plot: bool = True) Union[numpy.typing.NDArray, tuple[numpy.typing.NDArray, numpy.typing.NDArray]] [source]
Compute the normalized cross-correlation between two images.
Parameters
img1 : NDArray The first image. img2 : NDArray, optional The second images. If None, it computes the auto-correlation. The default is None. axes : Sequence[int], optional Axes along which to compute the cross-correlation. The default is (-2, -1). t_match : bool, optional Whether to perform the cross-correlation for template matching. The default is False. mode_full : bool, optional Whether to return the “full” or “same” convolution size. The default is True. compute_profile : bool, optional Whether to compute the azimuthal integration of the cross-correlation or not. The default is True. plot : bool, optional Whether to plot the profile of the cross-correlation curve. The default is True.
Returns
NDArray The one-dimensional cross-correlation curve.