pyepri.multisrc
This module contains low-level operators related to multisources 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 EPR projections of from a sequence of 2D source images. |
|
Compute EPR projections of from a sequence of 2D source images (output in Fourier domain). |
|
Compute EPR projections of from a sequence of 2D source images (output in Fourier domain, half of the full spectrum). |
|
Perform EPR backprojection from 2D multisources EPR projections (adjoint of the |
|
Perform EPR backprojection from 2D multisources EPR projections provided in Fourier domain. |
|
Perform EPR backprojection from 2D multisources EPR projections provided in Fourier domain. |
|
Compute 2D Toeplitz kernels allowing fast computation of a |
|
Perform a |
|
Compute EPR projections of from a sequence of 3D source images. |
|
Compute EPR projections of from a sequence of 3D source images (output in Fourier domain). |
|
Compute EPR projections of from a sequence of 3D source images (output in Fourier domain, half of the full spectrum). |
|
Perform EPR backprojection from 3D multisources EPR projections (adjoint of the |
|
Perform EPR backprojection from 3D multisources EPR projections provided in Fourier domain. |
|
Perform EPR backprojection from 3D multisources EPR projections provided in Fourier domain. |
|
Compute 3D Toeplitz kernels allowing fast computation of a |
|
Perform a |
|
Factorized consistency checks for functions in the |
Module Contents
- pyepri.multisrc.proj2d(u, delta, B, h, fgrad, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 2D source images.
This function can be used to simulate EPR projections from a mixture of 2D sources, in different experimental conditions (e.g., different microwave power).
In the following, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed using
K = len(u)
andL = max(len(h), len(fgrad))
.- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 2D array corresponding to the j-th source image of the mixture 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.
h (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The EPR projections are evaluated in the frequency domain (through their dicrete 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
proj – A sequence with length L containing the EPR projections synthesized for each experimental setting. More precisely,
proj[i]
is an array with shape(fgrad[i].shape[1], len(B))
andproj[i][k,:]
corresponds the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.proj2d_fft(u, delta, B, fft_h, fgrad, backend=None, eps=1e-06, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 2D source images (output in Fourier domain).
- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 2D array corresponding to the j-th source image of the mixture 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.
fft_h (sequence of sequence of complex array_like (with type backend.cls)) –
Contains the discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,fft_h[i][j]
is a monodimensional array corresponding to the FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(fft_h) = 1
, we assume thatfft_h[i][j] = fft_h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
fft_proj – A sequence with length L containing the discrete Fourier coefficients of the EPR projections synthesized for each experimental setting. More precisely,
fft_proj[i]
is an array with shape(fgrad[i].shape[1], len(B))
andfft_proj[i][k,:]
corresponds the FFT of the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of complex array_like (with type backend.cls)
See also
- pyepri.multisrc.proj2d_rfft(u, delta, B, rfft_h, fgrad, backend=None, eps=1e-06, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 2D source images (output in Fourier domain, half of the full spectrum).
- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 2D array corresponding to the j-th source image of the mixture 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.
rfft_h (sequence of sequence of complex array_like (with type backend.cls)) –
Contains half of the discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,rfft_h[i][j]
is a monodimensional array corresponding to the real input FFT (rfft) ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(rfft_h) = 1
, we assume thatrfft_h[i][j] = rfft_h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
rfft_proj – A sequence with length L containing half of the discrete Fourier coefficients of the EPR projections synthesized for each experimental setting. More precisely,
rfft_proj[i]
is an array with shape(fgrad[i].shape[1], 1+len(B)//2)
andrfft_proj[i][k,:]
corresponds the real input FFT (rfft) of the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of complex array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections (adjoint of the
proj2d
operation).In the following, and as we did in the documentation of the
pyepri.multisrc.proj2d()
function, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed usingK = len(u)
andL = len(proj)
.- Parameters:
proj (sequence of array_like (with type backend.cls)) – A sequence with length L containing the EPR projections associated to each experimental setup. More precisely,
proj[i][k,:]
corresponds to the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds the the shape of the j-th source image.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The backprojection process involves the computation of the 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes associated to the input projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections provided in Fourier domain.
- Parameters:
fft_proj (sequence of complex array_like (with type backend.cls)) – A sequence with length L containing the discrete Fourier coefficients of the EPR projections associated to each experimental setup. More precisely,
fft_proj[i][k,:]
corresponds to the FFT of the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of complex array_like (with type backend.cls)) –
Contains the conjugated discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,fft_h_conj[i][j]
is a monodimensional array corresponding to the conjugate of the FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(fft_h_conj) = 1
, we assume thatfft_h_conj[i][j] = fft_h_conj[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds to the shape of the j-th source image. optional input is in fact mandatory when no preallocated array is given (i.e., whenout=None
).out (sequence of array_like (with type backend.cls), optional) – Preallocated sequence of output arrays (with shape
out[j].shape = out_shape[j]
forj in range(K)
).eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes associated to the input projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections provided in Fourier domain.
- Parameters:
rfft_proj (sequence of complex array_like (with type backend.cls)) – A sequence with length L containing half of the discrete Fourier coefficients of the EPR projections associated to each experimental setup. More precisely,
rfft_proj[i][k,:]
corresponds to the real input FFT (rfft) of the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of complex array_like (with type backend.cls)) –
Contains half of the conjugated discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,rfft_h_conj[i][j]
is a monodimensional array corresponding to the conjugate of the real input FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(rfft_h_conj) = 1
, we assume thatrfft_h_conj[i][j] = rfft_h_conj[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds to the shape of the j-th source image. optional input is in fact mandatory when no preallocated array is given (i.e., whenout=None
).out (sequence of array_like (with type backend.cls), optional) – Preallocated sequence of output arrays (with shape
out[j].shape = out_shape[j]
forj in range(K)
).eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes associated to the input projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.compute_2d_toeplitz_kernels(B, h, delta, fgrad, src_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]
Compute 2D Toeplitz kernels allowing fast computation of a
proj2d
followed by abackproj2d
operation.In the following, and as we did in the documentation of the
pyepri.multisrc.proj2d()
function, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed usingK = len(src_shape)
andL = max(len(h), len(fgrad))
.- 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.
h (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).
fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.src_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
src_shape[j]
corresponds the the shape of the j-th source image.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The evaluation of the sequence of output Toeplitz kernels involves the computation Fourier coefficients of real-valued signals. 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
phi – Output sequence of 2D cross source kernels (
phi[k][j]
is the cross source kernel associated to source k (backward) and source j (forward)).- Return type:
sequence of sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.apply_2d_toeplitz_kernels(u, rfft2_phi, backend=None, notest=False)[source]
Perform a
proj2d
followed by abackproj2d
operation using precomputed Toeplitz kernels provided in Fourier domain.- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 2D array corresponding to the j-th source image of the mixture to be projected.rfft2_phi (sequence of sequence of complex array_like (with type backend.cls)) – Sequence of real input FFT of the 2D cross sources Toeplitz kernels computed using
pyepri.multisrc.compute_2d_toeplitz_kernels()
.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
u[0]
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – Sequence of output sources (stored in the same order as in
u
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.proj3d(u, delta, B, h, fgrad, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 3D source images.
This function can be used to simulate EPR projections from a mixture of 3D sources, in different experimental conditions (e.g., different microwave power).
In the following, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed using
K = len(u)
andL = max(len(h), len(fgrad))
.- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 3D array corresponding to the j-th source image of the mixture 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.
h (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The EPR projections are evaluated in the frequency domain (through their dicrete 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
proj – A sequence with length L containing the EPR projections synthesized for each experimental setting. More precisely,
proj[i]
is an array with shape(fgrad[i].shape[1], len(B))
andproj[i][k,:]
corresponds the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.proj3d_fft(u, delta, B, fft_h, fgrad, backend=None, eps=1e-06, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 3D source images (output in Fourier domain).
- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 3D array corresponding to the j-th source image of the mixture 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.
fft_h (sequence of sequence of complex array_like (with type backend.cls)) –
Contains the discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,fft_h[i][j]
is a monodimensional array corresponding to the FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(fft_h) = 1
, we assume thatfft_h[i][j] = fft_h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
fft_proj – A sequence with length L containing the discrete Fourier coefficients of the EPR projections synthesized for each experimental setting. More precisely,
fft_proj[i]
is an array with shape(fgrad[i].shape[1], len(B))
andfft_proj[i][k,:]
corresponds the FFT of the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of complex array_like (with type backend.cls)
See also
- pyepri.multisrc.proj3d_rfft(u, delta, B, rfft_h, fgrad, backend=None, eps=1e-06, nodes=None, notest=False)[source]
Compute EPR projections of from a sequence of 3D source images (output in Fourier domain, half of the full spectrum).
- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 3D array corresponding to the j-th source image of the mixture 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.
rfft_h (sequence of sequence of complex array_like (with type backend.cls)) –
Contains half of the discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,rfft_h[i][j]
is a monodimensional array corresponding to the real input FFT (rfft) ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(rfft_h) = 1
, we assume thatrfft_h[i][j] = rfft_h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
rfft_proj – A sequence with length L containing half of the discrete Fourier coefficients of the EPR projections synthesized for each experimental setting. More precisely,
rfft_proj[i]
is an array with shape(fgrad[i].shape[1], 1+len(B)//2)
andrfft_proj[i][k,:]
corresponds the real input FFT (rfft) of the EPR projection of the mixture u with field gradientfgrad[i][:,k]
sampled over the grid B.- Return type:
sequence of complex array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections (adjoint of the
proj3d
operation).In the following, and as we did in the documentation of the
pyepri.multisrc.proj2d()
function, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed usingK = len(u)
andL = len(proj)
.- Parameters:
proj (sequence of array_like (with type backend.cls)) – A sequence with length L containing the EPR projections associated to each experimental setup. More precisely,
proj[i][k,:]
corresponds to the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds the the shape of the j-th source image.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The backprojection process involves the computation of the 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the frequency nodes associated to the input projections. If not given, the nodes sequence is automatically inferred from B, h, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections provided in Fourier domain.
- Parameters:
fft_proj (sequence of complex array_like (with type backend.cls)) – A sequence with length L containing the discrete Fourier coefficients of the EPR projections associated to each experimental setup. More precisely,
fft_proj[i][k,:]
corresponds to the FFT of the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of complex array_like (with type backend.cls)) –
Contains the conjugated discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,fft_h_conj[i][j]
is a monodimensional array corresponding to the conjugate of the FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(fft_h_conj) = 1
, we assume thatfft_h_conj[i][j] = fft_h_conj[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds to the shape of the j-th source image. optional input is in fact mandatory when no preallocated array is given (i.e., whenout=None
).out (sequence of array_like (with type backend.cls), optional) – Preallocated sequence of output arrays (with shape
out[j].shape = out_shape[j]
forj in range(K)
).eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes associated to the input projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.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 multisources EPR projections provided in Fourier domain.
- Parameters:
rfft_proj (sequence of complex array_like (with type backend.cls)) – A sequence with length L containing half of the discrete Fourier coefficients of the EPR projections associated to each experimental setup. More precisely,
rfft_proj[i][k,:]
corresponds to the real input FFT (rfft) of the EPR projection of the multisources mixture acquired with field gradientfgrad[:,k]
and i-th experimental setup.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 (sequence of sequence of complex array_like (with type backend.cls)) –
Contains half of the conjugated discrete Fourier coefficients of the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,rfft_h_conj[i][j]
is a monodimensional array corresponding to the conjugate of the real input FFT ofh[i][j]
(the reference spectrum of the j-th source in the i-th experimental setting).When
L > 1
andlen(rfft_h_conj) = 1
, we assume thatrfft_h_conj[i][j] = rfft_h_conj[0][j]`
for alli in range(L)
and allj in range(K)
.fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.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
.out_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
out_shape[j]
corresponds to the shape of the j-th source image. optional input is in fact mandatory when no preallocated array is given (i.e., whenout=None
).out (sequence of array_like (with type backend.cls), optional) – Preallocated sequence of output arrays (with shape
out[j].shape = out_shape[j]
forj in range(K)
).eps (float, optional) – Precision requested (>1e-16).
nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes associated to the input projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – A sequence with lenght K containing the backprojected source images (with shape
out[j].shape = out_shape[j]
forj in range(K)
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.compute_3d_toeplitz_kernels(B, h, delta, fgrad, src_shape, backend=None, eps=1e-06, rfft_mode=True, nodes=None, notest=False)[source]
Compute 3D Toeplitz kernels allowing fast computation of a
proj3d
followed by abackproj3d
operation.In the following, and as we did in the documentation of the
pyepri.multisrc.proj3d()
function, the index j shall refer to the j-th source image, while the index i shall refer to the i-th experimental setup. We shall denote by K the number of sources and by L the number of different experimental setup (those numbers are computed usingK = len(src_shape)
andL = max(len(h), len(fgrad))
.- 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.
h (sequence of sequence of array_like (with type backend.cls)) –
Contains the reference spectra (sampled over the grid
B
) of each individual source (j) for each experimental setup (i). More precisely,h[i][j]
is a monodimensional array corresponding to the sampling overB
of the reference spectrum of the j-th source in the i-th experimental setting.When
L > 1
andlen(h) = 1
, we assume thath[i][j] = h[0][j]`
for alli in range(L)
and allj in range(K)
.delta (float) – Pixel size given in a length unit denoted below as [length-unit] (can be centimeter (cm), millimeter (mm), …).
fgrad (sequence of array_like (with type backend.cls)) –
A sequence with length L containing the coordinates of the field gradient vector used for each experimental setting. More precisely,
fgrad[i][:,k]
corresponds to the (X,Y,Z) coordinates of the field gradient vector associated to the k-th EPR projection in the i-th experimental setting.When
L > 1
andlen(fgrad) = 1
, we assume thatfgrad[i][j] = fgrad[0][i]`
for alli in range(L)
and allj in range(K)
.src_shape (sequence of sequence of int) – Sequence made of each source shape. More precisely,
src_shape[j]
corresponds the the shape of the j-th source image.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
.eps (float, optional) – Precision requested (>1e-16).
rfft_mode (bool, optional) – The evaluation of the sequence of output Toeplitz kernels involves the computation Fourier coefficients of real-valued signals. 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, userfft_mode=False
, to enable standard (complex) FFT mode and compute all the Fourier coefficients.nodes (sequence of dict, optional) –
A sequence of length L containing the precomputed frequency nodes used to evaluate the output projections for each experimental setting. If not given, the nodes sequence is automatically inferred from B, delta and fgrad.
When
L > 1
andlen(nodes) = 1
, we assume thatnodes[i] = nodes[0]`
for alli in range(L)
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
phi – Output sequence of 3D cross source kernels (
phi[k][j]
is the cross source kernel associated to source k (backward) and source j (forward)).- Return type:
sequence of sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc.apply_3d_toeplitz_kernels(u, rfft3_phi, backend=None, notest=False)[source]
Perform a
proj3d
followed by abackproj3d
operation using precomputed Toeplitz kernels provided in Fourier domain.- Parameters:
u (sequence of array_like (with type backend.cls)) – A sequence with length K containing the source images. More precisely,
u[j]
must be a 3D array corresponding to the j-th source image of the mixture to be projected.rfft3_phi (sequence of sequence of complex array_like (with type backend.cls)) – Sequence of real input FFT of the 3D cross sources Toeplitz kernels computed using
pyepri.multisrc.compute_3d_toeplitz_kernels()
.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
u[0]
.notest (bool, optional) – Set
notest=True
to disable consistency checks.
- Returns:
out – Sequence of output sources (stored in the same order as in
u
).- Return type:
sequence of array_like (with type backend.cls)
See also
- pyepri.multisrc._check_nd_inputs_(ndims, backend, B=None, delta=None, fgrad=None, u=None, h=None, proj=None, eps=None, nodes=None, rfft_mode=None, out_shape=None, fft_h=None, fft_h_conj=None, rfft_h=None, rfft_h_conj=None, fft_proj=None, rfft_proj=None, src_shape=None, rfft2_phi=None, rfft3_phi=None, out=None)[source]
Factorized consistency checks for functions in the
pyepri.multisrc
submodule.