pyepri.monosrc

This module contains low-level operators related to single source EPR imaging (projection, backprojection, projection-backprojection using Toeplitz kernels). Detailed mathematical definitions of the operators are provided in the Mathematical definitions section of the PyEPRI documentation.

Functions

compute_2d_frequency_nodes(B, delta, fgrad[, backend, ...])

Compute 2D irregular frequency nodes involved in 2D projection & backprojection operations.

proj2d(u, delta, B, h, fgrad[, backend, eps, ...])

Compute EPR projections of a 2D image (adjoint of the backproj2d operation).

proj2d_fft(u, delta, B, fft_h, fgrad[, backend, eps, ...])

Compute EPR projections of a 2D image (output in Fourier domain).

proj2d_rfft(u, delta, B, rfft_h, fgrad[, backend, ...])

Compute EPR projections of a 2D image (output in Fourier domain, half of the full spectrum).

backproj2d(proj, delta, B, h, fgrad, out_shape[, ...])

Perform EPR backprojection from 2D EPR projections (adjoint of the proj2d operation).

backproj2d_fft(fft_proj, delta, B, fft_h_conj, fgrad)

Perform EPR backprojection from 2D EPR projections provided in Fourier domain.

backproj2d_rfft(rfft_proj, delta, B, rfft_h_conj, fgrad)

Perform EPR backprojection from 2D EPR projections provided in Fourier domain (half of the full spectrum).

compute_2d_toeplitz_kernel(B, h1, h2, delta, fgrad, ...)

Compute 2D Toeplitz kernel allowing fast computation of a proj2d followed by a backproj2d operation.

apply_2d_toeplitz_kernel(u, rfft2_phi[, backend, notest])

Perform a proj2d followed by a backproj2d operation using a precomputed Toeplitz kernel provided in Fourier domain.

compute_3d_frequency_nodes(B, delta, fgrad[, backend, ...])

Compute 3D irregular frequency nodes involved in 3D projection & backprojection operations.

proj3d(u, delta, B, h, fgrad[, backend, eps, ...])

Compute EPR projections of a 3D image (adjoint of the backproj3d operation).

proj3d_fft(u, delta, B, fft_h, fgrad[, backend, eps, ...])

Compute EPR projections of a 3D image (output in Fourier domain).

proj3d_rfft(u, delta, B, rfft_h, fgrad[, backend, ...])

Compute EPR projections of a 3D image (output in Fourier domain, half of the full spectrum).

backproj3d(proj, delta, B, h, fgrad, out_shape[, ...])

Perform EPR backprojection from 3D EPR projections (adjoint of the proj3d operation).

backproj3d_fft(fft_proj, delta, B, fft_h_conj, fgrad)

Perform EPR backprojection from 3D EPR projections provided in Fourier domain.

backproj3d_rfft(rfft_proj, delta, B, rfft_h_conj, fgrad)

Perform EPR backprojection from 3D EPR projections provided in Fourier domain (half of the full spectrum).

compute_3d_toeplitz_kernel(B, h1, h2, delta, fgrad, ...)

Compute 3D Toeplitz kernel allowing fast computation of a proj3d followed by a backproj3d operation.

apply_3d_toeplitz_kernel(u, rfft3_phi[, backend, notest])

Perform a proj3d followed by a backproj3d operation using a precomputed Toeplitz kernel provided in Fourier domain.

_check_nd_inputs_(ndims, B, delta, fgrad, backend[, ...])

Factorized consistency checks for functions in the pyepri.monosrc submodule.

Module Contents

pyepri.monosrc.compute_2d_frequency_nodes(B, delta, fgrad, backend=None, rfft_mode=True, notest=False)[source]

Compute 2D irregular frequency nodes involved in 2D projection & backprojection operations.

Parameters:
  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (B, fgrad).

  • rfft_mode (bool, optional) – Set rfft_mode=True to compute only half of the frequency nodes (to be combined with the use of real FFT functions in further processing). Otherwise, set rfft_mode=False to compute all frequency nodes (to be combined with the use of complex FFT functions in further processing).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

nodes – A dictionary with content {'x': x, 'y': y, 'indexes': indexes, 'rfft_mode': rfft_mode} where

  • x is a one dimensional array containing the frequency nodes along the horizontal axis;

  • y is a one dimensional array, with same length as x, containing the frequency nodes along the vertical axis;

  • indexes is a one dimensional array, with same length as x, corresponding to the indexes where should be dispatched the computed Fourier coefficients in the rfft of the 2D projections.

  • rfft_mode a bool specifying whether the frequency nodes cover half of the frequency domain (rfft_mode=True) or the full frequency domain (rfft_mode=False).

Return type:

dict

pyepri.monosrc.proj2d(u, delta, B, h, fgrad, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]

Compute EPR projections of a 2D image (adjoint of the backproj2d operation).

Parameters:
  • u (array_like (with type backend.cls)) – Two-dimensional array corresponding to the input 2D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • h (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum sampled over the grid B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The EPR projections are evaluated in the frequency domain (through their discrete Fourier coefficients) before being transformed back to the B-domain. Set rfft_mode=True to enable real FFT mode (only half of the Fourier coefficients will be computed to speed-up computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode and compute all the Fourier coefficients.

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] corresponds the EPR projection of u with field gradient fgrad[:,k] sampled over the grid B.

Return type:

array_like (with type backend.cls)

pyepri.monosrc.proj2d_fft(u, delta, B, fft_h, fgrad, backend=None, eps=1e-06, out=None, nodes=None, notest=False)[source]

Compute EPR projections of a 2D image (output in Fourier domain).

Parameters:
  • u (array_like (with type backend.cls)) – Two-dimensional array corresponding to the input 2D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • fft_h (complex array_like (with type backend.cls)) – One dimensional array with length len(B) containing the discrete Fourier coefficients of the reference spectrum sampled over B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, rfft_h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • out (array_like (with type backend.cls), optional) – Preallocated output complex array with shape (fgrad.shape[1], len(B)).

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] contains the discrete Fourier coefficients of the EPR projection of u with field gradient fgrad[:,k].

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.proj2d_rfft(u, delta, B, rfft_h, fgrad, backend=None, eps=1e-06, out=None, nodes=None, notest=False)[source]

Compute EPR projections of a 2D image (output in Fourier domain, half of the full spectrum).

Parameters:
  • u (array_like (with type backend.cls)) – Two-dimensional array corresponding to the input 2D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • rfft_h (complex array_like (with type backend.cls)) – One dimensional array with length 1+len(B)//2 containing half of the discrete Fourier coefficients of the reference spectrum sampled over B (corresponds to the signal computed using the real input fft (rfft) of the reference spectrum).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, rfft_h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • out (array_like (with type backend.cls), optional) – Preallocated output array with shape (fgrad.shape[1], 1 + len(B)//2).

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, 1+len(B)//2) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] corresponds to half of the discrete Fourier coefficients of the EPR projection of u with field gradient fgrad[:,k].

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.backproj2d(proj, delta, B, h, fgrad, out_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]

Perform EPR backprojection from 2D EPR projections (adjoint of the proj2d operation).

Parameters:
  • proj (array_like (with type backend.cls)) – Two-dimensional array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1]) such that proj[k,:] corresponds to the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • h (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum sampled over the grid B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • out_shape (integer or integer tuple of length 2) – Shape of the output image out_shape = out.shape = (N1, N2) (or out_shape = N when N1 = N2 = N).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (proj, B, h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The backprojection process involves the computation of discrete Fourier coefficients of the input projections. Set rfft_mode=True to enable real FFT mode (only half of the Fourier coefficients will be computed to speed-up computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode and compute all the Fourier coefficients.

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – A two-dimensional array with specified shape corresponding to the backprojected image.

Return type:

array_like (with type backend.cls)

pyepri.monosrc.backproj2d_fft(fft_proj, delta, B, fft_h_conj, fgrad, backend=None, out_shape=None, out=None, eps=1e-06, nodes=None, notest=False)[source]

Perform EPR backprojection from 2D EPR projections provided in Fourier domain.

Parameters:
  • fft_proj (complex array_like (with type backend.cls)) –

    Two-dimensional array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1]) containing the EPR projections in Fourier domain.

    More precisely, fft_proj[k,:] corresponds to the FFT of the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • fft_h_conj (complex array_like (with type backend.cls)) – One dimensional array with length len(B) containing the conjugate of half of the discrete Fourier coefficients of the reference spectrum sampled over B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (fft_proj, B, fft_h_conj, fgrad).

  • out_shape (integer or integer tuple of length 2, optional) – Shape of the output image out_shape = out.shape = (N1, N2) (or out_shape = N when N1 = N2 = N). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N1, N2) and complex data type. If out_shape is specifed, the shape must match (i.e., we must have out.shape == out_shape), otherwise, out_shape is inferred from out.

  • eps (float, optional) – Precision requested (>1e-16).

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Backprojected image in complex format (imaginary part should be close to zero and can be thrown away).

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.backproj2d_rfft(rfft_proj, delta, B, rfft_h_conj, fgrad, backend=None, out_shape=None, out=None, eps=1e-06, nodes=None, notest=False)[source]

Perform EPR backprojection from 2D EPR projections provided in Fourier domain (half of the full spectrum).

Parameters:
  • rfft_proj (complex array_like (with type backend.cls)) –

    Two-dimensional array with shape (Nproj, 1+len(B)//2) (where Nproj = fgrad.shape[1]) containing the EPR projections in Fourier domain (half of the spectrum).

    More precisely, rfft_proj[k,:] corresponds to the real FFT (rfft) of the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • rfft_h_conj (complex array_like (with type backend.cls)) – One dimensional array with length 1+len(B)//2 containing the conjugate of half of the discrete Fourier coefficients of the reference spectrum sampled over B (corresponds to the conjugate of the signal computed using the real input FFT (rfft) of the reference spectrum).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (rfft_proj, B, rfft_h_conj, fgrad).

  • out_shape (integer or integer tuple of length 2, optional) – Shape of the output image out_shape = out.shape = (N1, N2) (or out_shape = N when N1 = N2 = N). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N1, N2) and complex data type. If out_shape is specifed, the shape must match (i.e., we must have out.shape == out_shape), otherwise, out_shape is inferred from out.

  • eps (float, optional) – Precision requested (>1e-16).

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Backprojected image in complex format (imaginary part should be close to zero and can be thrown away).

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.compute_2d_toeplitz_kernel(B, h1, h2, delta, fgrad, out_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, return_rfft2=False, notest=False)[source]

Compute 2D Toeplitz kernel allowing fast computation of a proj2d followed by a backproj2d operation.

Parameters:
  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • h1 (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum involved in the forward (proj2d) operation (and sampled over the grid B).

  • h2 (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum involved in the backward (backproj2d) operation (and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (2, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input array (B, h1, h2, fgrad).

  • out_shape (integer or integer tuple of length 2) – Shape of the output kernel out_shape = phi.shape = (M1, M2). The kernel shape should be twice the EPR image shape (i.e., denoting by (N1, N2) the shape of the EPR image, we should have (M1, M2) = (2*N1, 2*N2)).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The computation of the Toeplitz kernel involves the computation of discrete Fourier coefficients of real-valued signals. Set rfft_mode=True to enable real FFT mode (speed-up the computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode.

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output kernel. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_2d_frequency_nodes().

  • return_rfft2 (bool, optional) – Set return_rfft2 to return the real input FFT (rfft2) of the computed two-dimensional kernel (instead of the kernel itself).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

phi – Computed Toeplitz kernel (or its two-dimensional real input FFT when return_rfft2 is True).

Return type:

array_like (with type backend.cls)

pyepri.monosrc.apply_2d_toeplitz_kernel(u, rfft2_phi, backend=None, notest=False)[source]

Perform a proj2d followed by a backproj2d operation using a precomputed Toeplitz kernel provided in Fourier domain.

Parameters:
  • u (array_like (with type backend.cls)) – Two-dimensional array corresponding to the input 2D image to be projected-backprojected.

  • rfft2_phi (complex array_like (with type backend.cls)) – real input FFT of the 2D Toeplitz kernel computed using pyepri.monosrc.compute_2d_toeplitz_kernel().

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, rfft2_phi).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – output projected-backprojected image.

Return type:

array_like (with type backend.cls)

pyepri.monosrc.compute_3d_frequency_nodes(B, delta, fgrad, backend=None, rfft_mode=True, notest=False)[source]

Compute 3D irregular frequency nodes involved in 3D projection & backprojection operations.

Parameters:
  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • fgrad (array_like (with type backend.cls)) –

    Three dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (B, fgrad).

  • rfft_mode (bool, optional) – Set rfft_mode=True to compute only half of the frequency nodes (to be combined with the use of real FFT functions in further processing). Otherwise, set rfft_mode=False to compute all frequency nodes (to be combined with the use of complex FFT functions in further processing).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

nodes – A dictionary with content {'x': x, 'y': y, 'z': z, 'indexes': indexes, 'rfft_mode': rfft_mode} where

  • x is a one dimensional array containing the frequency nodes along the X-axis;

  • y is a one dimensional array, with same length as x, containing the frequency nodes along the Y-axis;

  • z is a one dimensional array, with same length as x, containing the frequency nodes along the Z-axis;

  • indexes is a one dimensional array, with same length as x, corresponding to the indexes where should be dispatched the computed Fourier coefficients in the rfft of the 3D projections.

  • rfft_mode a bool specifying whether the frequency nodes cover half of the frequency domain (rfft_mode=True) or the full frequency domain (rfft_mode=False).

Return type:

dict

See also

proj3d, backproj3d

pyepri.monosrc.proj3d(u, delta, B, h, fgrad, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]

Compute EPR projections of a 3D image (adjoint of the backproj3d operation).

Parameters:
  • u (array_like (with type backend.cls)) – Three-dimensional array corresponding to the input 3D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • h (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum sampled over the grid B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The EPR projections are evaluated in the frequency domain (through their discrete Fourier coefficients) before being transformed back to the B-domain. Set rfft_mode=True to enable real FFT mode (only half of the Fourier coefficients will be computed to speed-up computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode and compute all the Fourier coefficients.

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] corresponds the EPR projection of u with field gradient fgrad[:,k] sampled over the grid B.

Return type:

array_like (with type backend.cls)

pyepri.monosrc.proj3d_fft(u, delta, B, fft_h, fgrad, backend=None, eps=1e-06, out=None, nodes=None, notest=False)[source]

Compute EPR projections of a 3D image (output in Fourier domain).

Parameters:
  • u (array_like (with type backend.cls)) – Three-dimensional array corresponding to the input 2D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • fft_h (complex array_like (with type backend.cls)) – One dimensional array with length len(B) containing the discrete Fourier coefficients of the reference spectrum sampled over B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, fft_h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • out (array_like (with type backend.cls), optional) – Preallocated output complex array with shape (fgrad.shape[1], len(B)).

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] contains the discrete Fourier coefficients of the EPR projection of u with field gradient fgrad[:,k].

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.proj3d_rfft(u, delta, B, rfft_h, fgrad, backend=None, eps=1e-06, out=None, nodes=None, notest=False)[source]

Compute EPR projections of a 3D image (output in Fourier domain, half of the full spectrum).

Parameters:
  • u (array_like (with type backend.cls)) – Three-dimensional array corresponding to the input 2D image to be projected.

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) –

    One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

    WARNING: this function assumes that the range covered by B is large enough so that the computed EPR projections are fully supported by B. Using a too small range for B will result in unrealistic projections due to B-domain aliasing phenomena.

  • rfft_h (complex array_like (with type backend.cls)) – One dimensional array with length 1+len(B)//2 containing half of the discrete Fourier coefficients of the reference spectrum sampled over B (corresponds to the signal computed using the real input fft (rfft) of the reference spectrum).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, B, rfft_h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • out (array_like (with type backend.cls), optional) – Preallocated output array with shape (fgrad.shape[1], 1 + len(B)//2).

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Output array with shape (Nproj, 1+len(B)//2) (where Nproj = fgrad.shape[1] corresponds to the number of computed projections) such that out[k,:] corresponds to half of the discrete Fourier coefficients of the EPR projection of u with field gradient fgrad[:,k].

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.backproj3d(proj, delta, B, h, fgrad, out_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]

Perform EPR backprojection from 3D EPR projections (adjoint of the proj3d operation).

Parameters:
  • proj (array_like (with type backend.cls)) – Two-dimensional array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1]) such that proj[k,:] corresponds to the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • h (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum sampled over the grid B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • out_shape (integer or integer tuple of length 3) – Shape of the output image out_shape = out.shape = (N1, N2, N3) (or out_shape = N when N1 = N2 = N3 = N).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (proj, B, h, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The backprojection process involves the computation of discrete Fourier coefficients of the input projections. Set rfft_mode=True to enable real FFT mode (only half of the Fourier coefficients will be computed to speed-up computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode and compute all the Fourier coefficients.

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – A three-dimensional array with specified shape corresponding to the backprojected image.

Return type:

array_like (with type backend.cls)

pyepri.monosrc.backproj3d_fft(fft_proj, delta, B, fft_h_conj, fgrad, backend=None, out_shape=None, out=None, eps=1e-06, nodes=None, notest=False)[source]

Perform EPR backprojection from 3D EPR projections provided in Fourier domain.

Parameters:
  • fft_proj (complex array_like (with type backend.cls)) –

    Two-dimensional array with shape (Nproj, len(B)) (where Nproj = fgrad.shape[1]) containing the EPR projections in Fourier domain.

    More precisely, fft_proj[k,:] corresponds to the FFT of the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • fft_h_conj (complex array_like (with type backend.cls)) – One dimensional array with length len(B) containing the conjugate of half of the discrete Fourier coefficients of the reference spectrum sampled over B.

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (fft_proj, B, fft_h_conj, fgrad).

  • out_shape (integer or integer tuple of length 3, optional) – Shape of the output image out_shape = out.shape = (N1, N2, N3) (or out_shape = N when N1 = N2 = N3 = N). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N1, N2, N3) and complex data type. If out_shape is specifed, the shape must match (i.e., we must have out.shape == out_shape), otherwise, out_shape is inferred from out.

  • eps (float, optional) – Precision requested (>1e-16).

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Backprojected 3D image in complex format (imaginary part should be close to zero and can be thrown away).

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.backproj3d_rfft(rfft_proj, delta, B, rfft_h_conj, fgrad, backend=None, out_shape=None, out=None, eps=1e-06, nodes=None, notest=False)[source]

Perform EPR backprojection from 3D EPR projections provided in Fourier domain (half of the full spectrum).

Parameters:
  • rfft_proj (complex array_like (with type backend.cls)) –

    Two-dimensional array with shape (Nproj, 1+len(B)//2) (where Nproj = fgrad.shape[1]) containing the EPR projections in Fourier domain (half of the spectrum).

    More precisely, rfft_proj[k,:] corresponds to the real FFT (rfft) of the k-th EPR projection (acquired with field gradient fgrad[:,k] and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • rfft_h_conj (complex array_like (with type backend.cls)) – One dimensional array with length 1+len(B)//2 containing the conjugate of half of the discrete Fourier coefficients of the reference spectrum sampled over B (corresponds to the conjugate of the signal computed using the real input FFT (rfft) of the reference spectrum).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (rfft_proj, B, rfft_h_conj, fgrad).

  • out_shape (integer or integer tuple of length 2, optional) – Shape of the output image out_shape = out.shape = (N1, N2, N3) (or out_shape = N when N1 = N2 = N3 = N). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N1, N2, N3) and complex data type. If out_shape is specifed, the shape must match (i.e., we must have out.shape == out_shape), otherwise, out_shape is inferred from out.

  • eps (float, optional) – Precision requested (>1e-16).

  • nodes (dict, optional) – Precomputed frequency nodes associated to the input projections. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – Backprojected 3D image in complex format (imaginary part should be close to zero and can be thrown away).

Return type:

complex array_like (with type backend.cls)

pyepri.monosrc.compute_3d_toeplitz_kernel(B, h1, h2, delta, fgrad, out_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, return_rfft3=False, notest=False)[source]

Compute 3D Toeplitz kernel allowing fast computation of a proj3d followed by a backproj3d operation.

Parameters:
  • B (array_like (with type backend.cls)) – One dimensional array corresponding to the homogeneous magnetic field sampling grid, with unit denoted below as [B-unit] (can be Gauss (G), millitesla (mT), …), to use to compute the projections.

  • h1 (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum involved in the forward (proj3d) operation (and sampled over the grid B).

  • h2 (array_like (with type backend.cls)) – One dimensional array with same length as B corresponding to the reference spectrum involved in the backward (backproj3d) operation (and sampled over the grid B).

  • delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).

  • fgrad (array_like (with type backend.cls)) –

    Two dimensional array with shape (3, fgrad.shape[1]) such that fgrad[:,k] corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection to be computed.

    The physical unit of the field gradient should be consistent with that of B and delta, i.e., fgrad must be provided in [B-unit] / [length-unit] (e.g., G/cm, mT/cm, …).

  • out_shape (integer or integer tuple of length 3) – Shape of the output kernel out_shape = phi.shape = (M1, M2, M3). The kernel shape should be twice the EPR image shape (i.e., denoting by (N1, N2, N3) the shape of the EPR image, we should have (M1, M2, M3) = (2*N1, 2*N2, 2*N3)).

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (B, h1, h2, fgrad).

  • eps (float, optional) – Precision requested (>1e-16).

  • rfft_mode (bool, optional) – The computation of the Toeplitz kernel involves the computation of discrete Fourier coefficients of real-valued signals. Set rfft_mode=True to enable real FFT mode (speed-up the computation and reduce memory usage). Otherwise, use rfft_mode=False, to enable standard (complex) FFT mode.

  • nodes (dict, optional) – Precomputed frequency nodes used to evaluate the output kernel. If not given, nodes will be automatically inferred from B, delta and fgrad using pyepri.monosrc.compute_3d_frequency_nodes().

  • return_rfft3 (bool, optional) – Set return_rfft3 to return the real input FFT (rfft3) of the computed three-dimensional kernel (instead of the kernel itself).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

phi – Computed Toeplitz kernel (or its three-dimensional real input FFT when return_rfft3 is True).

Return type:

array_like (with type backend.cls)

pyepri.monosrc.apply_3d_toeplitz_kernel(u, rfft3_phi, backend=None, notest=False)[source]

Perform a proj3d followed by a backproj3d operation using a precomputed Toeplitz kernel provided in Fourier domain.

Parameters:
  • u (array_like (with type backend.cls)) – Three-dimensional array corresponding to the input 3D image to be projected-backprojected.

  • rfft3_phi (complex array_like (with type backend.cls)) – real input FFT of the 3D Toeplitz kernel computed using pyepri.monosrc.compute_3d_toeplitz_kernel().

  • backend (<class 'pyepri.backends.Backend'> or None, optional) –

    A numpy, cupy or torch backend (see pyepri.backends module).

    When backend is None, a default backend is inferred from the input arrays (u, rfft3_phi).

  • notest (bool, optional) – Set notest=True to disable consistency checks.

Returns:

out – output projected-backprojected image.

Return type:

array_like (with type backend.cls)

pyepri.monosrc._check_nd_inputs_(ndims, B, delta, fgrad, backend, u=None, h=None, h1=None, h2=None, proj=None, eps=None, nodes=None, rfft_mode=None, out_proj=None, out_im=None, out_shape=None, fft_h=None, fft_h_conj=None, rfft_h=None, rfft_h_conj=None, fft_proj=None, rfft_proj=None, return_rfft2=None, return_rfft3=None)[source]

Factorized consistency checks for functions in the pyepri.monosrc submodule.