spinguin.core.relaxation#

This module provides functions for calculating relaxation superoperators.

G0(tensor1: ndarray, tensor2: ndarray, l: int) float[source]#

Computes the time correlation function at t = 0, G(0), for two Cartesian tensors.

This is the multiplicative factor in front of the exponential decay for the isotropic rotational diffusion model.

Source: Eq. 70 from Hilla & Vaara: Rela2x: Analytic and automatic NMR relaxation theory https://doi.org/10.1016/j.jmr.2024.107828

Parameters:
  • tensor1 (ndarray) – Cartesian tensor 1.

  • tensor2 (ndarray) – Cartesian tensor 2.

  • l (int) – Common rank of the tensors.

Returns:

G_0 – Time correlation function evaluated at t = 0.

Return type:

float

Q_constant(S: float, Q_moment: float) float[source]#

Calculates the nuclear quadrupolar coupling constant in (rad/s) / (V/m^2).

Parameters:
  • S (float) – Spin quantum number.

  • Q_moment (float) – Nuclear quadrupole moment (in units of m^2).

Returns:

Q_const – Quadrupolar coupling constant.

Return type:

float

Q_intr_tensors(efg: ndarray, spins: ndarray, quad: ndarray) ndarray[source]#

Calculates the quadrupolar interaction tensors for a spin system.

Parameters:
  • efg (ndarray) – A 3-dimensional array specifying the electric field gradient tensors. Must be given in atomic units.

  • spins (ndarray) – A 1-dimensional array specifying the spin quantum numbers for each spin.

  • quad (ndarray) – A 1-dimensional array specifying the quadrupolar moments. Must be given in the units of m^2.

Returns:

Q_tensors – Quadrupolar interaction tensors.

Return type:

ndarray

dd_constant(y1: float, y2: float) float[source]#

Calculates the dipole-dipole coupling constant (excluding the distance).

Parameters:
  • y1 (float) – Gyromagnetic ratio of the first spin in units of rad/s/T.

  • y2 (float) – Gyromagnetic ratio of the second spin in units of rad/s/T.

Returns:

dd_const – Dipole-dipole coupling constant in units of rad/s * m^3.

Return type:

float

dd_coupling_tensors(xyz: ndarray, gammas: ndarray) ndarray[source]#

Calculates the dipole-dipole coupling tensor between all nuclei in the spin system.

Parameters:
  • xyz (ndarray) – A 2-dimensional array specifying the cartesian coordinates in the XYZ format for each nucleus in the spin system. Must be given in the units of Å.

  • gammas (ndarray) – A 1-dimensional array specifying the gyromagnetic ratios for each nucleus in the spin system. Must be given in the units of rad/s/T.

Returns:

dd_tensors – Array of dimensions (N, N, 3, 3) containing the 3x3 tensors between all nuclei.

Return type:

ndarray

get_sop_T(basis: ndarray, spins: ndarray, l: int, q: int, interaction_type: Literal['CSA', 'Q', 'DD'], spin_1: int, spin_2: int = None, sparse: bool = True) ndarray | csc_array[source]#

Helper function for the relaxation module. Calculates the coupled product superoperators for different interaction types.

Parameters:
  • basis (ndarray) – A 2-dimensional array containing the basis set that contains sequences of integers describing the Kronecker products of irreducible spherical tensors.

  • spins (ndarray) – A 1-dimensional array specifying the spin quantum numbers for each spin in the system.

  • l (int) – Operator rank.

  • q (int) – Operator projection.

  • interaction_type ({'CSA', 'Q', 'DD'}) – Describes the interaction type. Possible options are “CSA”, “Q”, and “DD”, which stand for chemical shift anisotropy, quadrupolar coupling, and dipole-dipole coupling, respectively.

  • spin_1 (int) – Index of the first spin.

  • spin_2 (int, optional) – Index of the second spin. Leave empty for single-spin interactions (e.g., CSA).

  • sparse (bool, default=True) – Specifies whether to return the superoperator as a sparse or dense array.

Returns:

sop – Coupled spherical tensor superoperator of rank l and projection q.

Return type:

ndarray or csc_array

ldb_thermalization(R: ndarray | csc_array, H_left: ndarray | csc_array, T: float, zero_value: float = 1e-18) ndarray | csc_array[source]#

Applies the Levitt-Di Bari thermalization to the relaxation superoperator.

Parameters:
  • R (ndarray or csc_array) – Relaxation superoperator to be thermalized.

  • H_left (ndarray or csc_array) – Left-side coherent Hamiltonian superoperator.

  • T (float) – Temperature of the spin bath in Kelvin.

  • zero_value (float, default=1e-18) – This threshold is used to estimate the convergence in the matrix exponential and to eliminate small values from the array.

Returns:

R – Thermalized relaxation superoperator.

Return type:

ndarray or csc_array

process_interactions(intrs: dict, zero_value: float) dict[source]#

Processes all incoherent interactions and organizes them by rank. Disregards interactions below a specified threshold.

Parameters:
  • intrs (dict) – A dictionary where the keys represent the interaction type, and the values contain the interaction tensors and the ranks.

  • zero_value (float) – If the eigenvalues of the interaction tensor, estimated using the 1-norm, are smaller than this threshold, the interaction is ignored.

Returns:

interactions – A dictionary where the interactions are organized by rank. The values contain all interactions with meaningful strength. The interactions are tuples in the format (“interaction”, spin_1, spin_2, tensor).

Return type:

dict

shielding_intr_tensors(shielding: ndarray, gammas: ndarray, B: float) ndarray[source]#

Calculates the shielding interaction tensors for a spin system.

Parameters:
  • shielding (ndarray) – A 3-dimensional array specifying the nuclear shielding tensors for each nucleus. The tensors must be given in the units of ppm.

  • gammas (ndarray) – A 1-dimensional array specifying the gyromagnetic ratios for each nucleus in the spin system. Must be given in the units of rad/s/T.

  • B (float) – External magnetic field in units of T.

Returns:

shielding_tensors – Array of shielding tensors.

Return type:

ndarray

sop_R_phenomenological(basis: ndarray, R1: ndarray, R2: ndarray, sparse: bool = True) ndarray | csc_array[source]#

Constructs the relaxation superoperator from given R1 and R2 values for each spin.

Parameters:
  • basis (ndarray) – A 2-dimensional array containing the basis set that contains sequences of integers describing the Kronecker products of irreducible spherical tensors.

  • R1 (ndarray) – A one dimensional array containing the longitudinal relaxation rates in 1/s for each spin. For example: np.array([1.0, 2.0, 2.5])

  • R2 (ndarray) – A one dimensional array containing the transverse relaxation rates in 1/s for each spin. For example: np.array([2.0, 4.0, 5.0])

  • sparse (bool, default=True) – Specifies whether to construct the relaxation superoperator as sparse or dense array.

Returns:

sop_R – Relaxation superoperator.

Return type:

ndarray or csc_array

sop_R_random_field()[source]#

TODO PERTTU?

sop_R_redfield(basis: ndarray, sop_H: csc_array, tau_c: float, spins: ndarray, B: float = None, gammas: ndarray = None, quad: ndarray = None, xyz: ndarray = None, shielding: ndarray = None, efg: ndarray = None, include_antisymmetric: bool = False, include_dynamic_frequency_shift: bool = False, relative_error: float = 1e-06, interaction_zero: float = 1e-09, aux_zero: float = 1e-18, relaxation_zero: float = 1e-12, parallel_dim: int = 1000, sparse: bool = True) ndarray | csc_array[source]#

Calculates the relaxation superoperator using Redfield relaxation theory.

Sources:

Eq. 54 from Hilla & Vaara: Rela2x: Analytic and automatic NMR relaxation theory https://doi.org/10.1016/j.jmr.2024.107828

Eq. 24 and 25 from Goodwin & Kuprov: Auxiliary matrix formalism for interaction representation transformations, optimal control, and spin relaxation theories https://doi.org/10.1063/1.4928978

Parameters:
  • basis (ndarray) – A 2-dimensional array containing the basis set that contains sequences of integers describing the Kronecker products of irreducible spherical tensors.

  • sop_H (ndarray or csc_array) – Coherent part of the Hamiltonian superoperator.

  • tau_c (float) – Rotational correlation time in seconds.

  • spins (ndarray) – A 1-dimensional array specifying the spin quantum numbers of the system.

  • B (float) – External magnetic field in Tesla.

  • gammas (ndarray, default=None) – A 1-dimensional array specifying the gyromagnetic ratios for each spin. Must be defined in the units of rad/s/T.

  • quad (ndarray, default=None) – A 1-dimensional array specifying the quadrupolar moments for each spin. Must be defined in the units of m^2.

  • xyz (ndarray, default=None) – A 2-dimensional array where the rows contain the Cartesian coordinates for each spin in the units of Å.

  • shielding (ndarray, default=None) – A 3-dimensional array where the shielding tensors are specified for each spin in the units of Å.

  • efg (ndarray, default=None) – A 3-dimensional array where the electric field gradient tensors are specified for each spin in the units of Å.

  • include_antisymmetric (bool, default=False) – Specifies whether the antisymmetric component of the CSA is included. This is usually very small and can be neglected.

  • include_dynamic_frequency_shift (bool, default=False) – Specifies whether the dynamic frequency shifts are included. This corresponds to the imaginary part of the relaxation superoperator that causes small shifts to the resonance frequencies. This effect is usually very small and can be neglected.

  • relative_error (float=1e-6) – Relative error for the Redfield integral that is calculated using auxiliary matrix method. Smaller values correspond to more accurate integrals.

  • interaction_zero (float=1e-9) – If the eigenvalues of the interaction tensor, estimated using the 1-norm, are smaller than this threshold, the interaction is ignored.

  • aux_zero (float=1e-18) – This threshold is used to estimate the convergence of the Taylor series when exponentiating the auxiliary matrix, and also to eliminate small values from the arrays in the matrix exponential squaring step.

  • relaxation_zero (float=1e-12) – Smaller values than this threshold are eliminated from the relaxation superoperator before returning the array.

  • parallel_dim (int=1000) – If the basis set dimension is larger than this value, the Redfield integrals are calculated in parallel. Otherwise, the integrals are calculated in serial.

  • sparse (bool=True) – Specifies whether to calculate the relaxation superoperator as sparse or dense array.

Returns:

sop_R – Relaxation superoperator.

Return type:

ndarray or csc_array

sop_R_redfield_term(l: int, q: int, type_r: str, spin_r1: int, spin_r2: int, tensor_r: ndarray, top_l_shared: dict, top_r_shared: dict, bottom_r_shared: dict, t_max: float, aux_zero: float, relaxation_zero: float, sop_Ts: dict, interactions: dict) tuple[int, int, str, int, int, csc_array][source]#

Helper function for the Redfield relaxation theory. This function calculates one term of the relaxation superoperator and enables the use of parallel computation.

NOTE: This function returns some of the input parameters to display the progress in the computation of the total Redfield relaxation superoperator.

Parameters:
  • l (int) – Operator rank.

  • q (int) – Operator projection.

  • type_r (str) – Interaction type. Possible options are “CSA”, “Q”, and “DD”.

  • spin_r1 (int) – Index of the first spin in the interaction.

  • spin_r2 (int) – Index of the second spin in the interaction. Leave empty for single-spin interactions (e.g., CSA).

  • tensor_r (np.ndarray) – Interaction tensor for the right-hand interaction.

  • top_l_shared (dict) – Dictionary containing the shared top left block of the auxiliary matrix.

  • top_r_shared (dict) – Dictionary containing the shared top right block of the auxiliary matrix.

  • bottom_r_shared (dict) – Dictionary containing the shared bottom right block of the auxiliary matrix.

  • t_max (float) – Integration limit for the auxiliary matrix method.

  • aux_zero (float) – Threshold for the convergence of the Taylor series when exponentiating the auxiliary matrix.

  • relaxation_zero (float) – Values below this threshold are disregarded in the construction of the relaxation superoperator term.

  • sop_Ts (dict) – Dictionary containing the shared coupled T superoperators for different interactions.

  • interactions (dict) – Dictionary containing the interactions organized by rank.

Returns:

  • l (int) – Operator rank.

  • q (int) – Operator projection.

  • type_r (str) – Interaction type.

  • spin_r1 (int) – Index of the first spin.

  • spin_r2 (int) – Index of the second spin.

  • sop_R_term (csc_array) – Relaxation superoperator term for the given interaction.

sop_R_sr2k(basis: ndarray, spins: ndarray, gammas: ndarray, chemical_shifts: ndarray, J_couplings: ndarray, sop_R: csc_array, B: float, sparse: bool = True) ndarray | csc_array[source]#

Calculates the scalar relaxation of the second kind (SR2K) based on Abragam’s formula.

Parameters:
  • basis (ndarray) – A 2-dimensional array containing the basis set that contains sequences of integers describing the Kronecker products of irreducible spherical tensors.

  • spins (ndarray) – A 1-dimensional array specifying the spin quantum numbers for each spin in the system.

  • gammas (ndarray) – A 1-dimensional array specifying the gyromagnetic ratios for each nucleus in the spin system. Must be given in the units of rad/s/T.

  • chemical_shifts (ndarray) – A 1-dimensional array containing the chemical shifts of each spin in the units of ppm.

  • J_couplings (ndarray) – A 2-dimensional array containing the scalar J-couplings between each spin in the units of Hz. Only the bottom triangle is considered.

  • sop_R (ndarray or csc_array) – Relaxation superoperator without scalar relaxation of the second kind.

  • B (float) – Magnetic field in units of T.

  • sparse (bool, default=True) – Specifies whether to return the relaxation superoperator as dense or sparse array.

Returns:

sop_R – Relaxation superoperator containing the contribution from scalar relaxation of the second kind.

Return type:

ndarray or csc_array

tau_c_l(tau_c: float, l: int) float[source]#

Calculates the rotational correlation time for a given rank l. Applies only for anisotropic rotationally modulated interactions (l > 0).

Source: Eq. 70 from Hilla & Vaara: Rela2x: Analytic and automatic NMR relaxation theory https://doi.org/10.1016/j.jmr.2024.107828

Parameters:
  • tau_c (float) – Rotational correlation time.

  • l (int) – Interaction rank.

Returns:

t_cl – Rotational correlation time for the given rank.

Return type:

float