corrct.physics.phase
Phase contrast support functions.
@author: Nicola VIGANÒ, CEA-IRIG, Grenoble, France
Module Contents
Functions
Compute the delta-over-beta parameter for a specific compound. |
|
Compute the delta-over-beta curves for the listed compounds, in the requested energy range. |
|
Plot frequency response of the wave propagation. |
|
Compute the phase contrast propagation filter for the given parameters. |
|
Apply a requested propagation filter to an image or a stack of images. |
API
- corrct.physics.phase.get_delta_beta(cmp_name: str, energy_keV: float, density: Union[float, None] = None) float [source]
Compute the delta-over-beta parameter for a specific compound.
Parameters
cmp_name : str Molar composition of a compound energy_keV : float Energy at which the d/b value should be computed density : Optional[float], optional Density of the compound, by default None
Returns
float The computed delta-over-beta.
- corrct.physics.phase.get_delta_beta_curves(compounds: collections.abc.Sequence[str], energy_keV_range: tuple[float, float, int] = (1.0, 800.0, 500), plot: bool = True) collections.abc.Sequence[numpy.typing.NDArray] [source]
Compute the delta-over-beta curves for the listed compounds, in the requested energy range.
Parameters
compounds : Sequence[str] Sequence of compounds. energy_keV_range : tuple[float, float, int], optional The energy range in keV, by default (1.0, 800.0, 500) plot : bool, optional Whether to plot the result, by default True
Returns
Sequence[NDArray] The computed delta-over-beta valueas.
- corrct.physics.phase._tie_freq_response(k2: numpy.typing.NDArray, dist_um: float, wlength_um: float, delta_beta: float) numpy.typing.NDArray [source]
- corrct.physics.phase._ctf_freq_response(k2: numpy.typing.NDArray, dist_um: float, wlength_um: float, delta_beta: float) numpy.typing.NDArray [source]
- corrct.physics.phase.plot_filter_responses(filter_length: int, pix_size_um: float, dist_um: float, wlength_um: float, delta_beta: float, domain: str = 'fourier') tuple [source]
Plot frequency response of the wave propagation.
Parameters
filter_length : int Length of the filter in pixels pix_size_um : float Pixel-wise of the detector in microns dist_um : float Distance of the detector from sample in microns wlength_um : float Wavelength of the wave in microns delta_beta : float Radio between the refraction index decrement and the absorption coefficient. domain : str Whether to plot Fourier or direct-space responses. By default, “Fourier”.
Returns
tuple Figure and axes of the plot
- corrct.physics.phase.get_propagation_filter(img_shape: Union[collections.abc.Sequence[int], numpy.typing.NDArray], pix_size_um: float, dist_um: float, wlength_um: float, delta_beta: float, filter_type: str = 'ctf', use_rfft: bool = False, plot_result: bool = False, dtype: numpy.typing.DTypeLike = np.float32) tuple[numpy.typing.NDArray, numpy.typing.NDArray] [source]
Compute the phase contrast propagation filter for the given parameters.
Parameters
img_shape : Union[Sequence[int], NDArray] Shape of the target image pix_size_um : float Pixel size of the detector (in microns) dist_um : float Propagation distance (in microns) wlength_um : float Wavelength of the radiation (in microns) delta_beta : float Delta-over-beta value for the given material filter_type : str, optional Type of the filter, by default “ctf”. Options: “ctf” (Contrast Transfer Function) | “tie” (Transport of Intensity Equation) use_rfft : bool, optional Whether to use the rfft, by default False plot_result : bool, optional Whether to plot the result, by default False dtype : DTypeLike, optional Data type to use, by default np.float32
Returns
tuple[NDArray, NDArray] The Fourier-space and real-space filters, respectively
Raises
ValueError When choosing an incorrect filter type
- corrct.physics.phase.apply_propagation_filter(data_wvu: numpy.typing.NDArray, pix_size_um: float, dist_um: float, wlength_um: float, delta_beta: float, filter_type: str = 'tie') numpy.typing.NDArray [source]
Apply a requested propagation filter to an image or a stack of images.
Parameters
data_wvu : NDArray The Image or stack of images pix_size_um : float Pixel size of the detector (in microns) dist_um : float Propagation distance (in microns) wlength_um : float Wavelength of the radiation (in microns) delta_beta : float Delta-over-beta value for the given material filter_type : str, optional Type of the filter, by default “tie”. Options: “ctf” (Contrast Transfer Function) | “tie” (Transport of Intensity Equation)
Returns
NDArray The filtered image