deerlab.dipolarkernel#

dipolarkernel(t, r, *, pathways=None, mod=None, bg=None, method='fresnel', excbandwidth=inf, orisel=None, g=None, integralop=True, gridsize=1000, complex=False, clearcache=False, memorylimit=8, tinterp=None)[source]#

Compute the (multi-pathway) dipolar kernel operator.

The function computes the two-spin interaction dipolar kernel which accounts for the structure and shape of the echo modulation due to intra- and inter-molecular dipolar interactions over a sample. The dipolar kernel is constructed from the contributions arising from different dipolar pathways [1] during a pulse sequence. The two-spin dipolar kernel reads:

K(\mathbf{t},r) = \Lambda_0 + \sum_k \lambda_k 
\int_0^\pi \mathrm{d}\theta \sin(\theta) 
\exp\Bigg(-i \frac{D\mathbf{\delta}_k (\mathbf{t} - \mathbf{t}_{\mathrm{ref},k})}{r^3}(1-3\cos(\theta)^2) \Bigg)
\prod_k B(\lambda_k, \mathbf{\delta}_k (\mathbf{t} - \mathbf{t}_{\mathrm{ref},k}))

where \mathbf{t}=(t_1,\dots,t_D) are the experimental time incrementation vectors of a D-dimensional experiment and r represents the interspin distance. The dipolar pathways are characterized by their amplitudes \lambda_k, their refocusing times \mathbf{t}_{\mathrm{ref},k}, and their (sub)harmonic factors \mathbf{\delta}_k. The constant factor \Lambda_0 represents the total amplitude of all unmodulated contributions. The function B(\lambda,t) represents the background function due to the intermolecular contributions. The dipolar signal due to two-spin contributions given a distance distribution P(r) is given by the linear Fredholm integral

V(\mathbf{t}) = \int_0^\infty \mathrm{d}r P(r) K(\mathbf{t},r)

The function can also compute three-spin interaction dipolar kernels that account for the additional the structure and shape of the echo modulation in multi-spin systems [2]. The three-spin dipolar kernel is constructed from the contributions arising from the combinations of different pair dipolar pathways [2] during a pulse sequence. The three-spin dipolar kernel reads:

K^{(3)}(\mathbf{t},\mathbf{r}) = & \Lambda_0 + \sum_k \lambda_k \frac{1}{2\pi}
\int_0^{2\pi}\mathrm{d}\varphi 
\int_0^{\pi}\mathrm{d}\theta_1
\sin(\theta_1)
\exp\Bigg(-i \sum_q^3 \frac{D \delta_{k,q}(\mathbf{t} - \mathbf{t}_{\mathrm{ref},k,q})}{r_q^3}
(1 - 3\cos(\theta_{1q})^2 \Bigg) 
\prod_k\prod_q^3 B(\lambda_k, \mathbf{\delta}_{k,q} (\mathbf{t} - \mathbf{t}_{\mathrm{ref},k,q}))

with \cos(\theta_{1q}) = \sin\chi_{1q} \cos\varphi \sin\theta_1 + \cos\chi_{1q} \cos\theta_1, where \mathbf{r} = (r_1,r_2,r_3) represents the vector of the three interspin distances between the three interacting spins, and \chi_{1q} are the angles between the first and q-th interspin vector. The three-spin dipolar pathways are now characterized by an amplitude \lambda_k, and by a set of three refocusing times \mathbf{t}_{\mathrm{ref},k,q} and three (sub)harmonic factors \mathbf{\delta}_{k,q}.

Due to the higher dimensionality, the three-spin kernel is not constructed as a full matrix but folds the three distance dimensions into a single one, such that the contribution to the dipolar echo modulation is given directly by

V^{(3)}(\mathbf{t}) = \int_0^\infty\int_0^\infty \int_0^\infty  \mathrm{d}\mathbf{r} K^{(3)}(\mathbf{t},\mathbf{r})

without the specification of a distance distribution.

Parameters:
tndarray or list thereof

Experimental time incrementation vectors (time coordinates) t, in microseconds. For multi-dimensional experiments, specify a list [t1,...,tD] of the time coordinates \mathbf{t} along each of the D different experimental dimensions.

rndarray or list thereof

Distance vector r, in nanometers. For three-spin interactions, a list [r1,r2,r3] of three vectors of interspin distances must be specified.

pathwayslist of dictionaries or None, optional

List of dipolar pathways.

For two-spin dipolar interactions, each pathway is defined as a dictionary containing the pathway’s amplitude \lambda_k ('amp'), refocusing time t_{\mathrm{ref},k} ('reftime') in microseconds, and (sub)harmonic \delta_k ('harmonic'). For example, for a one-dimensional experiment {'amp:'lambda, 'reftime':tref, 'harmonic':delta} or {'amp:'lambda, 'reftime':tref}. If 'harmonic' is not specified, it is assumed to be 1 (modulated). For an unmodulated pathway, specify only the amplitude, i.e. {'amp':Lambda0}.

For multi-dimensional experiments, the pathway refocusing times mathbf{t}_{mathrm{ref},k}` and harmonics \mathbf{\delta}_k along the different dimensions must be specified as a list, e.g. for a pathway in a two-dimensional experiment {'amp:'lambda, 'reftime':[tref1,tref2], 'harmonic':[delta1,delta2]}

For three-spin interactions, each pathway is defined as a dictionary containing the pathway’s amplitude \lambda_k ('amp'), a tuple of refocusing times t_{\mathrm{ref},k,q} ('reftime'), and a tuple of (sub)harmonics \delta_{k,q} ('harmonic'). For example, for a one-dimensional experiment {'amp:'lambda, 'reftime':(tref1,tref2,tref3), 'harmonic':(delta1,delta2,delta3)}.

If neither pathways or mod are specified (or set to None), pathways=[{'amp':1,'reftime':0}] is used as the default.

modscalar or None, optional

Modulation depth for single-pathway dipolar models. This definition is equivalent to pathways=[{'amp':1-mod},{'amp':mod,'reftime':0}]. If neither pathways or mod are specified (or None), mod=1 is used as the default.

bgcallable or array_like or None, optional

Dipolar background basis function. Must be a callable function, accepting a time-vector array as first input and a pathway amplitude as a second, i.e. bg = lambda t,lam: bg_model(t,par,lam). If set to None, no background function is included. For a single-pathway model, the numerical background decay can be passed as an array.

methodstring, optional

Numerical method for kernel matrix calculation:

  • 'fresnel' - uses Fresnel integrals for the kernel (fast, accurate)

  • 'integral' - uses explicit integration function (slow, accurate)

  • 'grid' - powder average via explicit grid integration (slow, inaccurate)

The default is 'fresnel'. For all three-spin dipolar kernels, the 'grid' method is used automatically.

oriselcallable or None, optional

Probability distribution of possible orientations of the interspin vector to account for orientation selection. Must be a function taking a value of the angle θ∈[0,π/2] between the interspin vector and the external magnetic field and returning the corresponding probability density. If specified as None (by default), a uniform distribution is assumed. Requires the 'grid' or 'integral' methods.

excbandwidthscalar, optional

Excitation bandwidth of the pulses in MHz to account for limited excitation bandwidth. If not specified, an infinite excitation bandwidth is used. Requires the 'grid' or 'integral' methods.

gscalar, 2-element array, optional

Electron g-values of the spin centers [g1, g2]. If a single g is specified, [g, g] is used. If not specified, g = 2.002319… is used for both spins.

complexboolean, optional

If enabled, return the complex-valued dipolar kernel, otherwise return the real-valued dipolar kernel. Disabled by default. Requires the 'fresnel' or 'grid' methods.

integralopboolean, optional

Return the kernel as an integral operator (i.e K = K*dr) or not (K). Usage as an integral operator means that the matrix operation V=K@P with a normalized distance distribution (i.e. trapz(r,P)==1) leads to a signal V with amplitude V(t=0)=1. The default is True.

gridsizescalar, optional

Number of points on the grid of powder orientations to be used in the 'grid' kernel calculation method. By default set to 1000 points. If the expected memory costs of the combined sizes of the t vector(s), r vector(s) and grid exceed the memorylimit value, a MemoryError will be raised.

clearcacheboolean, optional

Clear the cached dipolar kernels at the beginning of the function to save memory. Disabled by default.

memorylimitscalar, optional

Memory limit to be allocated for the dipolar kernel. If the requested kernel exceeds this limit, the execution will be stopped. The default is 8GB.

tinterparray_like, optional

Effective dipolar evolution time vector used for interpolation. If specified, the elementary dipolar kernel is computed based on tinterp and all contributions from the different pathways are interpolated based on it. This results in a significant reduction of the computation time. The vector tinterp must cover all possible time shifts of the original t by the different refocusing times of the dipolar pathways.

Returns:
Kndarray

Dipolar kernel matrix.

Notes

For a multi-pathway DEER [1], [2], [3] signal (e.g, 4-pulse DEER with 2+1 contribution, 5-pulse DEER with 4-pulse DEER residual signal, and more complicated experiments), pathways contains a list of pathway amplitudes and refocusing times (in microseconds). The background function specified as bg is used as basis function, and the actual multipathway background included into the kernel is computed using deerlab.dipolarbackground. The background in included in the dipolar kernel definition [3]. Optionally, the harmonic (1 = fundamental, 2 = first harmonic, etc.) can be specified for modeling RIDME signals or subharmonics (1/2 = first subharmonic) can be specified for modeling SIFTER signals.

References

[1] (1,2)

L. Fábregas Ibáñez, M. H. Tessmer, G. Jeschke, and S. Stoll. Dipolar pathways in dipolar EPR spectroscopy, Phys. Chem. Chem. Phys., 24 2022, 2504-2520

[2] (1,2,3)

L. Fábregas Ibáñez, M. H. Tessmer, G. Jeschke, and S. Stoll. Dipolar pathways in multi-spin and multi-dimensional dipolar EPR spectroscopy, Phys. Chem. Chem. Phys., 24 2022, 22645-22660

[3] (1,2)

L. Fábregas Ibáñez, G. Jeschke, and S. Stoll. DeerLab: A comprehensive toolbox for analyzing dipolar EPR spectroscopy data, Magn. Reson., 1, 209–224, 2020

[4]

L. Fábregas Ibáñez, and G. Jeschke Optimal background treatment in dipolar spectroscopy, Physical Chemistry Chemical Physics, 22, 1855–1868, 2020.

Examples

To specify a model with an unmodulated offset and a single dipolar pathway that refocuses at time t=0, use:

lam = 0.4  # modulation depth
pathways = [
    {'amp': [1-lam]},
    {'amp': lam, 'reftime': 0}
    ]

K = dl.dipolarkernel(t,r,pathways=pathways)

A shorthand input syntax equivalent to this input:

lam = 0.4 # modulation
K = dl.dipolarkernel(t,r,mod=lam)

To specify a three-pathway kernel of the 4-pulse DEER experiment that, e.g that includes the 2+1 contribution, use:

tau1 = 0.5  # First interpulse delay (µs) 
tau2 = 4.0  # Second interpulse delay (µs)

lam1 = 0.40  # Pathway #1 amplitude
lam2 = 0.10  # Pathway #2 amplitude
lam3 = 0.10  # Pathway #3 amplitude

tref1 = tau1       # Pathway #1 refocusing time
tref2 = tau1+tau2  # Pathway #2 refocusing time
tref3 = 0          # Pathway #3 refocusing time

# Define the list of dipolar pathways
pathways = [
    {'amp': 1-lam1-lam2-lam3},
    {'amp': lam1, 'reftime': tref1},
    {'amp': lam2, 'reftime': tref2},
    {'amp': lam3, 'reftime': tref3},
]

K = dl.dipolarkernel(t,r,pathways=pathways)

To specify a dipolar kernel for an experiment on a three-spin system based on a single-pathway model use:

tau1 = 0.5  # First interpulse delay (µs)         
lam = 0.40  # Pathway #1 amplitude

# Three dipolar pathways needed, one for each dipolar interaction in the three-spin system 
pathways = [
    {'amp': 1-3*lam},
    {'amp': lam, 'reftime': (tau1,None,None), 'harmonic': (1,0,0)},
    {'amp': lam, 'reftime': (None,tau1,None), 'harmonic': (0,1,0)},
    {'amp': lam, 'reftime': (None,None,tau1), 'harmonic': (0,0,1)},
]

K = dl.dipolarkernel(t, [r1,r2,r3], pathways=pathways)

To construct a kernel for a three-spin system studied with a two-dimensional experiment, with a single dipolar pathway refocusing at, e.g., t=(\tau_1,\tau_2) use:

# Three dipolar pathways needed, one for each dipolar interaction in the three-spin system 
pathways = [
    {'amp': 1-3*lam},
    {'amp': lam, 'reftime': ([tau1,tau2],None,None), 'harmonic': ([1,1],0,0)},
    {'amp': lam, 'reftime': (None,[tau1,tau2],None), 'harmonic': (0,[1,1],0)},
    {'amp': lam, 'reftime': (None,None,[tau1,tau2]), 'harmonic': (0,0,[1,1])},
]

K = dl.dipolarkernel(t, [r1,r2,r3], pathways=pathways)