pyepri.spectralspatial

This module contains low-level operators related to spectral-spatial EPR imaging (projection, backprojection, projection-backprojection using Toeplitz kernels).

Functions

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

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

compute_4d_weights(nodes[, backend, nrm, isign, ...])

Precompute weights to accelerate further proj4d & backproj4d calls.

proj4d(u, delta, B, fgrad[, backend, weights, eps, ...])

Compute EPR projections of a 4D image (adjoint of the backproj4d operation).

proj4d_fft(u, delta, B, fgrad[, backend, weights, ...])

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

proj4d_rfft(u, delta, B, fgrad[, backend, weights, ...])

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

backproj4d(proj, delta, B, fgrad, out_shape[, ...])

Perform EPR backprojection from 4D EPR projections (adjoint of the proj4d operation).

backproj4d_fft(fft_proj, delta, B, fgrad[, backend, ...])

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

backproj4d_rfft(rfft_proj, delta, B, fgrad[, backend, ...])

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

compute_4d_toeplitz_kernel(B, delta, fgrad, out_shape)

Compute 4D Toeplitz kernel allowing fast computation of a proj4d followed by a backproj4d operation.

apply_4d_toeplitz_kernel(u, rfft4_phi[, backend, notest])

Perform a proj4d followed by a backproj4d 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.spectralspatial submodule.

Module Contents

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

Compute 4D irregular frequency nodes involved in 4D 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 (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, 'xi': xi, 't', t, 'lt': lt 'indexes': indexes, 'idt': idt, 'idt0': idt0, 'rfft_mode': rfft_mode} where

  • x, y, z are the 3D frequency nodes computed using pyepri.monosrc.compute_3d_frequency_nodes()

  • indexes is a one dimensional array, with same length as x, y and z, corresponding to the indexes where should be dispatched the computed Fourier coefficients in the (r)fft of the 4D 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).

  • t and idt are one dimensional arrays such that t[idt] has the same length as x, y and z and represents the frequency nodes along the B-axis of the 4D image (the t array has no duplicate values). Said differently, the 4D frequency nodes involved in the 4D projection and backprojection operations are the (x, y, z, t[idt]) nodes.

  • idt0 corresponds to the location of the t == 0 index

  • lt is a 2D array such that lt[l,k] = l * t[k] (those values are involved in the computation of weights for 4D projection and backprojection functions).

Return type:

dict

pyepri.spectralspatial.compute_4d_weights(nodes, backend=None, nrm=1 + 0j, isign=1, make_contiguous=True, notest=False)[source]

Precompute weights to accelerate further proj4d & backproj4d calls.

Precomputing the weights is useful to save computation time when multiple evaluations of proj4d() or backproj4d() with option memory_usage=0 are needed.

Parameters:
  • nodes (dict, optional) – Frequency nodes computed using compute_4d_frequency_nodes().

  • 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 stored into nodes.

  • nrm (complex, optional) – Normalization parameter involved in the weights definition (see below).

  • isign (float, optional) – Must be equal to +1 or -1, used as sign in the complex exponential defining the weights (see below). The user must use isign = 1 to compute the weigths involved in proj4d() (also proj4d_fft() and proj4d_rfft()) and isign = -1 to compute those involved in backproj4d() (also backproj4d_fft() and backproj4d_rfft()).

  • make_contiguous (bool, optional) – Set make_contiguous = True to make the output array contiguous and ordered in row-major order using weights = weigths.ravel().reshape(weights.shape).

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

Returns:

weights – A two-dimensional array equal to w = (nrm * (backend.exp(isign * 1j * nodes[`lt`]))[:, nodes['idt']].

Return type:

complex array_like (with type backend.cls)

pyepri.spectralspatial.proj4d(u, delta, B, fgrad, backend=None, weights=None, eps=1e-06, rfft_mode=True, nodes=None, memory_usage=1, notest=False)[source]

Compute EPR projections of a 4D image (adjoint of the backproj4d operation).

Parameters:
  • u (array_like (with type backend.cls)) – Four-dimensional array corresponding to the input spectral-spatial 4D image to be projected (axis 0 is the spectral axis, axes 1, 2 and 3 correspond to the Y, X and Z spatial axes).

  • 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.

  • 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, fgrad).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=1).

    Note that those weights are only used when memory_usage=0.

  • 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

  • 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.spectralspatial.proj4d_fft(u, delta, B, fgrad, backend=None, weights=None, eps=1e-06, out=None, nodes=None, memory_usage=1, notest=False)[source]

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

Parameters:
  • u (array_like (with type backend.cls)) – Four-dimensional array corresponding to the input spectral-spatial 4D image to be projected (axis 0 is the spectral axis, axes 1, 2 and 3 correspond to the Y, X and Z spatial axes).

  • 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.

  • 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, fgrad).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=1).

    Note that those weights are only used when memory_usage=0.

  • 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

  • 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.spectralspatial.proj4d_rfft(u, delta, B, fgrad, backend=None, weights=None, eps=1e-06, out=None, nodes=None, memory_usage=1, notest=False)[source]

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

Parameters:
  • u (array_like (with type backend.cls)) – Four-dimensional array corresponding to the input spectral-spatial 4D image to be projected (axis 0 is the spectral axis, axes 1, 2 and 3 correspond to the Y, X and Z spatial axes).

  • 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.

  • 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, fgrad).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=1).

    Note that those weights are only used when memory_usage=0.

  • 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

  • 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.spectralspatial.backproj4d(proj, delta, B, fgrad, out_shape, backend=None, weights=None, eps=1e-06, rfft_mode=True, nodes=None, memory_usage=1, notest=False)[source]

Perform EPR backprojection from 4D EPR projections (adjoint of the proj4d 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.

  • 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 4) –

    Shape of the output image out_shape = out.shape = (N0, N1, N2, N3).

    Note: N0 should be equal to len(B).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=-1).

    Note that those weights are only used when memory_usage=0.

  • 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

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

Returns:

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

Return type:

array_like (with type backend.cls)

pyepri.spectralspatial.backproj4d_fft(fft_proj, delta, B, fgrad, backend=None, out_shape=None, out=None, weights=None, eps=1e-06, nodes=None, memory_usage=1, notest=False)[source]

Perform EPR backprojection from 4D 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.

  • 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, …).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=-1).

    Note that those weights are only used when memory_usage=0.

  • 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, fgrad).

  • out_shape (integer or integer tuple of length 4, optional) –

    Shape of the output image out_shape = out.shape = (N0, N2, N3). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

    Note: N0 should be equal to len(B).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N0, 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

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

Returns:

out – Backprojected 4D 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.spectralspatial.backproj4d_rfft(rfft_proj, delta, B, fgrad, backend=None, out_shape=None, out=None, weights=None, eps=1e-06, nodes=None, preserve_input=False, memory_usage=1, notest=False)[source]

Perform EPR backprojection from 4D 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.

  • 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, …).

  • weights (complex array_like (with type backend.cls), optional) –

    A two dimensional array with shape (len(B), len(nodes['idt'])) of weights precomputed using compute_4d_weights() (with option isign=-1).

    Note that those weights are only used when memory_usage=0.

  • 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, fgrad).

  • out_shape (integer or integer tuple of length 4, optional Shape) –

    of the output image out_shape = out.shape = (N0, N1, N2, N3). This optional input is in fact mandatory when no preallocated array is given (i.e., when out=None).

    Note: N0 should be equal to len(B).

  • out (complex array_like (with type backend.cls), optional) – Preallocated output array with shape (N0, 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 compute_4d_frequency_nodes().

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: (recommended) pretty good tradeoff between speed and reduced memory usage

    • memory_usage = 2: slow computation but very light memory usage

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

Returns:

out – Backprojected 4D 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.spectralspatial.compute_4d_toeplitz_kernel(B, delta, fgrad, out_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, return_rfft4=False, notest=False, memory_usage=1)[source]

Compute 4D Toeplitz kernel allowing fast computation of a proj4d followed by a backproj4d 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.

  • 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 4) – Shape of the output kernel out_shape = phi.shape = (M0, M1, M2, M3). The kernel shape should be twice the EPR image shape (i.e., denoting by (N0, N1, N2, N3) the shape of the 4D EPR image, we should have (M0, M1, M2, M3) = (2*N0, 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, 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 compute_4d_frequency_nodes().

  • return_rfft4 (bool, optional) – Set return_rfft4 to return the real input FFT (rfft4) of the computed four-dimensional kernel (instead of the kernel itself).

  • memory_usage (int, optional) –

    Specify the computation strategy (depending on your available memory budget).

    • memory_usage = 0: fast computation but memory demanding

    • memory_usage = 1: same as memory_usage = 0 (for consistency with proj4d() and backproj4d())

    • memory_usage = 2: slow computation but very light memory usage

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

Returns:

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

Return type:

array_like (with type backend.cls)

pyepri.spectralspatial.apply_4d_toeplitz_kernel(u, rfft4_phi, backend=None, notest=False)[source]

Perform a proj4d followed by a backproj4d operation using a precomputed Toeplitz kernel provided in Fourier domain.

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

  • rfft4_phi (complex array_like (with type backend.cls)) – real input FFT of the 4D Toeplitz kernel computed using compute_4d_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, rfft4_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.spectralspatial._check_nd_inputs_(ndims, B, delta, fgrad, backend, u=None, proj=None, eps=None, nodes=None, rfft_mode=None, out_shape=None, fft_proj=None, rfft_proj=None, return_rfft4=None, nrm=None, isign=None, make_contiguous=None, memory_usage=None, weights=None, preserve_input=None)[source]

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