soprano.nmr.tensor#
Contains the NMRTensor class, simplifying the process of diagonalisation of an NMR tensor as well as its representation in multiple conventions
Classes
|
Class containing an electric field gradient tensor, a subclass of NMRTensor. |
|
Class containing a magnetic shielding tensor, a subclass of NMRTensor |
|
Class containing an NMR tensor, useful to access all its most important properties and representations. |
- class soprano.nmr.tensor.ElectricFieldGradient(data, species, order='n', quadrupole_moment=None, gamma=None)[source]#
Bases:
NMRTensor
Class containing an electric field gradient tensor, a subclass of NMRTensor.
The EFG tensor is a symmetric 3x3 matrix, with a trace of zero and has therefore only 5 independent degrees of freedom. These are often captured using the following:
largest absolute eigenvalue (\(V_{zz}\))
the asymmetry parameter ( \(\eta\) )
α, β, γ Euler angles describing the orientation of the tensor.
The principal components of the tensor are the eigenvalues of the tensor, and they are usually sorted according to the NQR convention ( \(|V_{zz}| \geq |V_{yy}| \geq |V_{xx}|\) ). Note however, that Simpson uses the convention \(|V_{zz}| \geq |V_{xx}| \geq |V_{yy}|\) . This is can be compensated for by using the order parameter when creating the tensor.
Note that some conventions use the reduced anisotropy ( \(\zeta\) ) instead of the asymmetry parameter.
Initialise the EFGTensor
Create an EFGTensor object from a 3x3 matrix.
- Parameters:
data (np.ndarray or tuple) – 3x3 matrix containing the tensor, or pair [evals, evecs] for the symmetric part alone.
species (str) – Isotope symbol of the atom the tensor refers to. e.g ‘2H’, ‘13C’
order (str) – Order to use for eigenvalues/eigenvectors. Can be ‘i’ (ORDER_INCREASING), ‘d’ (ORDER_DECREASING), ‘h’ (ORDER_HAEBERLEN) or ‘n’ (ORDER_NQR). Default is ‘h’ for EFG tensors.
quadrupole_moment (float) – Quadrupole moment of the nucleus in millibarn.
gamma (float) – Nuclear gyromagnetic ratio in rad/(s*T). If not provided, will be looked up using the species.
- property Cq: float#
Calculates the quadrupolar constant in Hz for this EFG tensor. The constant will be zero for non-quadrupole active nuclei. The quadrupole moment used is that for the nucleus of the isotope specified in the species attribute.
This property is defined as
\[C_Q = \frac{e^2qQ}{h}\]in Hz. It is important to keep in mind that therefore this represents a frequency; the corresponding ‘omega’ (pulsation) would be the same value multiplied by 2*pi. This is, for example, exactly the value required as input in Simpson’s SPINSYS section.
- Returns:
Quadrupolar constant value in Hz.
- Return type:
float
- property PAS#
Returns the principal axis system (PAS) of the tensor.
- property Pq: float#
Calculates the quadrupolar product in Hz for this EFG tensor.
\[P_Q = C_Q (1+\frac{\eta_Q^2}{3})^{1/2}\]- Returns:
Quadrupolar product value in Hz.
- Return type:
float
- property Vzz#
Returns the largest absolute eigenvalue of the tensor ( the principal component of the EFG tensor). This should be in atomic units (a.u.) if read in from e.g. a magres file.
- __repr__()#
Return a string representation of the tensor
- Return type:
str
- __str__()#
Return a string representation of the tensor. Nicely print out the 3x3 matrix (data), the conventions, and the derived properties.
- Return type:
str
- property degeneracy#
Returns the degeneracy of the tensor. For example, a tensor with eigenvalues [1, 1, 1] has a degeneracy of 3. A tensor with eigenvalues [1, 1, 2] has a degeneracy of 2. A tensor with eigenvalues [1, 2, 3] has a degeneracy of 1.
- equivalent_euler_angles(convention='zyz', passive=False)#
Returns the equivalent Euler angles of the tensor.
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
- Returns:
- np.array – Euler angles in radians.
Size of the array is (4, 3) as there are 4 equivalent sets of Euler angles.
- equivalent_euler_to(other, convention='zyz', passive=False)#
Returns the equivalent Euler angles that rotate the tensor to the other tensor.
- property eta#
Returns the asymmetry parameter of the tensor.
- euler_angles(convention='zyz', passive=False, degrees=False)#
Euler angles of the Principal Axis System
Return Euler angles of the PAS for this tensor in the required convention (currently supported: zyz, zxz).
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
convention (str)
passive (bool)
- Returns:
np.ndarray – Euler angles in radians. Size of the array is (3,)
- Return type:
ndarray
- euler_to(other, convention='zyz', passive=False, eps=1e-06)#
Returns the Euler angles that rotate the tensor to the other tensor.
- get_MAS_full_max_linewidth(Bext)[source]#
Returns the MAS full maximum line width (Hz) of the central transition.
Defined as:
\[\Delta\nu = a\frac{(6+\eta)^2}{504}\]where \(a\) is the quadrupolar perturbation and \(\eta\) is the asymmetry parameter.
- Parameters:
Bext (float) – External magnetic field in T.
- Returns:
MAS full maximum line width in Hz.
- Return type:
float
- get_larmor_frequency(Bext)[source]#
Returns the Larmor frequency of the nucleus in an external magnetic field in Hz.
\[\nu_L = \gamma B_{ext} / (2\pi)\]where \(\gamma\) is the gyromagnetic ratio of the nucleus in rad/(s*T) and \(B_{ext}\) is the external magnetic field in T.
- Parameters:
Bext (float) – External magnetic field in T.
- Returns:
Larmor frequency in Hz.
- Return type:
float
- get_quadrupolar_perturbation(Bext)[source]#
Returns the perturbation of the quadrupolar Hamiltonian due to an external magnetic field. The perturbation is given by:
\[a = \frac{\nu_Q^2}{\nu_L} (I(I+1) - 3/2)\]where \(\nu_Q\) is the quadrupolar frequency and \(\nu_L\) is the Larmor frequency.
- Parameters:
Bext (float) – External magnetic field in T.
- Returns:
Perturbation in Hz.
- Return type:
float
- static make_dipolar(a, i, j, cell=[0, 0, 0], isotopes={}, isotope_i=None, isotope_j=None, rotation_axis=None)#
Create a dipolar NMRTensor
Create a dipolar NMR tensor from an atoms object and the indices of two atoms. Values are in Hz.
Args:a (ase.Atoms): Atoms object of the structure to compute thetensor fori (int): index of first atomj (int): index of second atomcell (np.array): vector of the cell of the second atom, forcouplings between atoms in different cells.By default is [0,0,0].isotopes (dict): dictionary of specific isotopes to use, by elementsymbol. If the isotope doesn’t exist an error willbe raised.isotope_i (int): isotope of atom i. To be used ifdifferent atoms of the same element are supposedto be of different isotopes. If None it will fallback on the previous definitions. Otherwise itoverrides everything else.isotope_j (int): isotope of atom j. See above.rotation_axis (np.array): an axis around which the selected pairis rotating. If present, the tensorwill be averaged for infinitely fastrotation around it.Returns:diptens (NMRTensor): an NMRTensor object with the dipolarcoupling matrix as data.
- property nuq: float#
Calculates the quadrupolar frequency in Hz for this EFG tensor. This is also known as the quadrupolar splitting parameter.
\[\nu_Q = \frac{3C_{Q}}{2I(2I-1)}\]where I is the nuclear spin of the nucleus.
See this for conventions: http://anorganik.uni-tuebingen.de/klaus/nmr/index.php?p=conventions/efg/quadtools
- Returns:
Quadrupolar frequency value in Hz.
- Return type:
float
- rotation_to(other)#
Returns the rotation matrix that rotates the tensor to the other tensor.
- property spherical_repr#
Spherical representation of the tensor
Returns a 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor. The isotropic part is the average of the trace of the tensor, the antisymmetric part is the difference between the tensor and its transpose divided by 2, and the symmetric part is the sum of the tensor and its transpose divided by 2 minus the isotropic part. This construction is such that the sum of the components is equal to the original tensor.
\[\sigma = \sigma_{iso} + \sigma_{A} + \sigma_{S}\]where
\[ \begin{align}\begin{aligned}\sigma_{iso} = \frac{1}{3} \text{Tr}(\sigma) \mathbf{I}\\\sigma_{A} = \frac{1}{2} (\sigma - \sigma^T)\\\sigma_{S} = \frac{1}{2} (\sigma + \sigma^T) - \sigma_{iso}\end{aligned}\end{align} \]- Returns:
np.ndarray – 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor
- property zeta#
Returns the reduced anisotropy of the tensor.
- class soprano.nmr.tensor.MagneticShielding(data, species, order='h', reference=None, gradient=-1.0, tag=None)[source]#
Bases:
NMRTensor
Class containing a magnetic shielding tensor, a subclass of NMRTensor
It provides easy access to common representations of the tensor in common conventions such as: IPAC, Herzfeld-Berger, Haeberlen, Maryland, and Mehring. For more information on these conventions, see the documentation for the corresponding NamedTuple and also [2]
Initialise the MSTensor
Create an MSTensor object from a 3x3 matrix.
- Parameters:
species (str) – Element or isotope symbol of the atom the tensor refers to. e.g ‘H’, ‘C’, ‘13C’
data (np.ndarray or tuple) – 3x3 matrix containing the tensor, or pair [evals, evecs] for the symmetric part alone.
order (str) – Order to use for eigenvalues/eigenvectors. Can be ‘i’ (ORDER_INCREASING), ‘d’ (ORDER_DECREASING), ‘h’ (ORDER_HAEBERLEN) or ‘n’ (ORDER_NQR). Default is ‘h’ for MS tensors.
references (dict) – Dictionary of references for the magnetic shielding to chemical shift conversion. Keys are element symbols, values are the reference chemical shift in ppm.
gradients (dict) – Dictionary of gradients for the magnetic shielding tensor. Keys are element symbols, values are the gradient. Any unspecified gradients will be set to -1.
tag (str) – Optional tag to identify the tensor. In a magres file, this would be the ‘ms_sometag’ though for most magres file ms is not decomposed into contributions and the tag is simply ‘ms’. By default, this is None.
reference (float | None)
gradient (float)
- class HaeberlenNotation(sigma_iso, sigma, delta, eta)[source]#
Bases:
NamedTuple
Haeberlen convention [1] follows this ordering:
\[|\sigma_{zz} - \sigma_{iso} | \geq |\sigma_{xx} - \sigma_{iso} | \geq |\sigma_{yy} - \sigma_{iso} |\]and uses the following parameters to describe the shielding tensor:
\(\sigma_{iso}\) : isotropic magnetic shielding
\(\sigma = \sigma_{zz} - \sigma_{iso}\) : the reduced anisotropy
\(\Delta = \sigma_{zz} - (\sigma_{xx} + \sigma_{yy}) / 2 = 3\sigma / 2\) : the anisotropy
\(\eta = (\sigma_{yy} - \sigma_{xx}) / \sigma\) : the asymmetry ( \(0 \leq \eta \leq +1\) )
Create new instance of HaeberlenNotation(sigma_iso, sigma, delta, eta)
- Parameters:
sigma_iso (float)
sigma (float)
delta (float)
eta (float)
- _asdict()#
Return a new dict which maps field names to their values.
- classmethod _make(iterable)#
Make a new HaeberlenNotation object from a sequence or iterable
- _replace(**kwds)#
Return a new HaeberlenNotation object replacing specified fields with new values
- count(value, /)#
Return number of occurrences of value.
- delta: float#
Alias for field number 2
- eta: float#
Alias for field number 3
- index(value, start=0, stop=9223372036854775807, /)#
Return first index of value.
Raises ValueError if the value is not present.
- sigma: float#
Alias for field number 1
- sigma_iso: float#
Alias for field number 0
- class HerzfeldBergerNotation(sigma_iso, omega, kappa)[source]#
Bases:
NamedTuple
Herzfeld-Berger convention [3] uses the following parameters:
\(\sigma_{iso}\) : isotropic magnetic shielding
\(\Omega = \sigma_{33} - \sigma_{11}\) (where these are the max and min principal components respectively): the span
\(\kappa = 3(\sigma_{iso} - \sigma_{22}) / \Omega\) (where \(\sigma_{22}\) is the median principal component): the skew
Create new instance of HerzfeldBergerNotation(sigma_iso, omega, kappa)
- Parameters:
sigma_iso (float)
omega (float)
kappa (float)
- _asdict()#
Return a new dict which maps field names to their values.
- classmethod _make(iterable)#
Make a new HerzfeldBergerNotation object from a sequence or iterable
- _replace(**kwds)#
Return a new HerzfeldBergerNotation object replacing specified fields with new values
- count(value, /)#
Return number of occurrences of value.
- index(value, start=0, stop=9223372036854775807, /)#
Return first index of value.
Raises ValueError if the value is not present.
- kappa: float#
Alias for field number 2
- omega: float#
Alias for field number 1
- sigma_iso: float#
Alias for field number 0
- class IUPACNotaion(sigma_iso, sigma_11, sigma_22, sigma_33)[source]#
Bases:
NamedTuple
IUPAC convention [4] (equivalent to the Mehring convention [5]). Follows the high frequency-positive order. Therefore: \(\sigma_{11}\) corresponds to the direction of least shielding, with the highest frequency, \(\sigma_{33}\) corresponds to the direction of highest shielding, with the lowest frequency.
\[\sigma_{11} \leq \sigma_{22} \leq \sigma_{33}\]The isotropic value, \(\sigma_{iso}\), is the average values of the principal components, and corresponds to the center of gravity of the line shape.
Note that the IUPAC convention is equivalent to the Mehring convention.
Create new instance of IUPACNotaion(sigma_iso, sigma_11, sigma_22, sigma_33)
- Parameters:
sigma_iso (float)
sigma_11 (float)
sigma_22 (float)
sigma_33 (float)
- _asdict()#
Return a new dict which maps field names to their values.
- classmethod _make(iterable)#
Make a new IUPACNotaion object from a sequence or iterable
- _replace(**kwds)#
Return a new IUPACNotaion object replacing specified fields with new values
- count(value, /)#
Return number of occurrences of value.
- index(value, start=0, stop=9223372036854775807, /)#
Return first index of value.
Raises ValueError if the value is not present.
- sigma_11: float#
Alias for field number 1
- sigma_22: float#
Alias for field number 2
- sigma_33: float#
Alias for field number 3
- sigma_iso: float#
Alias for field number 0
- class MarylandNotation(sigma_iso, omega, kappa)[source]#
Bases:
HerzfeldBergerNotation
The same as the Herzfeld-Berger notation.
Create new instance of HerzfeldBergerNotation(sigma_iso, omega, kappa)
- Parameters:
sigma_iso (float)
omega (float)
kappa (float)
- _asdict()#
Return a new dict which maps field names to their values.
- classmethod _make(iterable)#
Make a new HerzfeldBergerNotation object from a sequence or iterable
- _replace(**kwds)#
Return a new HerzfeldBergerNotation object replacing specified fields with new values
- count(value, /)#
Return number of occurrences of value.
- index(value, start=0, stop=9223372036854775807, /)#
Return first index of value.
Raises ValueError if the value is not present.
- kappa: float#
Alias for field number 2
- omega: float#
Alias for field number 1
- sigma_iso: float#
Alias for field number 0
- class MehringNotation(sigma_iso, sigma_11, sigma_22, sigma_33)[source]#
Bases:
IUPACNotaion
Mehring convention [5] is equivalent to the IUPAC convention.
Create new instance of IUPACNotaion(sigma_iso, sigma_11, sigma_22, sigma_33)
- Parameters:
sigma_iso (float)
sigma_11 (float)
sigma_22 (float)
sigma_33 (float)
- _asdict()#
Return a new dict which maps field names to their values.
- classmethod _make(iterable)#
Make a new IUPACNotaion object from a sequence or iterable
- _replace(**kwds)#
Return a new IUPACNotaion object replacing specified fields with new values
- count(value, /)#
Return number of occurrences of value.
- index(value, start=0, stop=9223372036854775807, /)#
Return first index of value.
Raises ValueError if the value is not present.
- sigma_11: float#
Alias for field number 1
- sigma_22: float#
Alias for field number 2
- sigma_33: float#
Alias for field number 3
- sigma_iso: float#
Alias for field number 0
- property PAS#
Returns the principal axis system (PAS) of the tensor.
- __repr__()#
Return a string representation of the tensor
- Return type:
str
- property degeneracy#
Returns the degeneracy of the tensor. For example, a tensor with eigenvalues [1, 1, 1] has a degeneracy of 3. A tensor with eigenvalues [1, 1, 2] has a degeneracy of 2. A tensor with eigenvalues [1, 2, 3] has a degeneracy of 1.
- property element#
Returns the element of the tensor. Species could have the isotope info, but here we just want the element. e.g. 13C -> C 1H -> H
- equivalent_euler_angles(convention='zyz', passive=False)#
Returns the equivalent Euler angles of the tensor.
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
- Returns:
- np.array – Euler angles in radians.
Size of the array is (4, 3) as there are 4 equivalent sets of Euler angles.
- equivalent_euler_to(other, convention='zyz', passive=False)#
Returns the equivalent Euler angles that rotate the tensor to the other tensor.
- euler_angles(convention='zyz', passive=False, degrees=False)#
Euler angles of the Principal Axis System
Return Euler angles of the PAS for this tensor in the required convention (currently supported: zyz, zxz).
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
convention (str)
passive (bool)
- Returns:
np.ndarray – Euler angles in radians. Size of the array is (3,)
- Return type:
ndarray
- euler_to(other, convention='zyz', passive=False, eps=1e-06)#
Returns the Euler angles that rotate the tensor to the other tensor.
- property haeberlen_values#
The magnetic shielding tensor in Haeberlen Notation.
- property herzfeldberger_values#
The magnetic shielding tensor in Herzfeld-Berger Notation.
- property iupac_values#
The magnetic shielding tensor in IUPAC Notation.
- static make_dipolar(a, i, j, cell=[0, 0, 0], isotopes={}, isotope_i=None, isotope_j=None, rotation_axis=None)#
Create a dipolar NMRTensor
Create a dipolar NMR tensor from an atoms object and the indices of two atoms. Values are in Hz.
Args:a (ase.Atoms): Atoms object of the structure to compute thetensor fori (int): index of first atomj (int): index of second atomcell (np.array): vector of the cell of the second atom, forcouplings between atoms in different cells.By default is [0,0,0].isotopes (dict): dictionary of specific isotopes to use, by elementsymbol. If the isotope doesn’t exist an error willbe raised.isotope_i (int): isotope of atom i. To be used ifdifferent atoms of the same element are supposedto be of different isotopes. If None it will fallback on the previous definitions. Otherwise itoverrides everything else.isotope_j (int): isotope of atom j. See above.rotation_axis (np.array): an axis around which the selected pairis rotating. If present, the tensorwill be averaged for infinitely fastrotation around it.Returns:diptens (NMRTensor): an NMRTensor object with the dipolarcoupling matrix as data.
- property maryland_values#
The magnetic shielding tensor in Maryland Notation.
- property mehring_values#
The magnetic shielding tensor in Mehring Notation.
- rotation_to(other)#
Returns the rotation matrix that rotates the tensor to the other tensor.
- set_reference(reference)[source]#
Set the reference chemical shift for this tensor.
- Parameters:
reference (float)
- property shift#
Returns the isotropic chemical shift of the tensor (ppm). If the reference is not set, will raise a ValueError. You can set the reference using the reference attribute.
- property shift_tensor#
Returns the chemical shift tensor in ppm.
- property spherical_repr#
Spherical representation of the tensor
Returns a 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor. The isotropic part is the average of the trace of the tensor, the antisymmetric part is the difference between the tensor and its transpose divided by 2, and the symmetric part is the sum of the tensor and its transpose divided by 2 minus the isotropic part. This construction is such that the sum of the components is equal to the original tensor.
\[\sigma = \sigma_{iso} + \sigma_{A} + \sigma_{S}\]where
\[ \begin{align}\begin{aligned}\sigma_{iso} = \frac{1}{3} \text{Tr}(\sigma) \mathbf{I}\\\sigma_{A} = \frac{1}{2} (\sigma - \sigma^T)\\\sigma_{S} = \frac{1}{2} (\sigma + \sigma^T) - \sigma_{iso}\end{aligned}\end{align} \]- Returns:
np.ndarray – 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor
- class soprano.nmr.tensor.NMRTensor(data, order='i')[source]#
Bases:
object
Class containing an NMR tensor, useful to access all its most important properties and representations.
Initialise the NMRTensor
Create an NMRTensor object from a 3x3 matrix.
- Parameters:
data (np.ndarray or tuple) – 3x3 matrix containing the tensor, or pair [evals, evecs] for the symmetric part alone.
order (str) – Order to use for eigenvalues/eigenvectors. Can be ‘i’ (ORDER_INCREASING), ‘d’ (ORDER_DECREASING), ‘h’ (ORDER_HAEBERLEN) or ‘n’ (ORDER_NQR). Default is ‘i’.
- property PAS#
Returns the principal axis system (PAS) of the tensor.
- __str__()[source]#
Return a string representation of the tensor. Nicely print out the 3x3 matrix (data), the conventions, and the derived properties.
- Return type:
str
- property degeneracy#
Returns the degeneracy of the tensor. For example, a tensor with eigenvalues [1, 1, 1] has a degeneracy of 3. A tensor with eigenvalues [1, 1, 2] has a degeneracy of 2. A tensor with eigenvalues [1, 2, 3] has a degeneracy of 1.
- equivalent_euler_angles(convention='zyz', passive=False)[source]#
Returns the equivalent Euler angles of the tensor.
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
- Returns:
- np.array – Euler angles in radians.
Size of the array is (4, 3) as there are 4 equivalent sets of Euler angles.
- equivalent_euler_to(other, convention='zyz', passive=False)[source]#
Returns the equivalent Euler angles that rotate the tensor to the other tensor.
- euler_angles(convention='zyz', passive=False, degrees=False)[source]#
Euler angles of the Principal Axis System
Return Euler angles of the PAS for this tensor in the required convention (currently supported: zyz, zxz).
- Parameters:
use (convention {str} -- Euler angles convention to)
(default – {‘zyz’})
angles (passive {bool} -- Whether to return the passive Euler)
(default – {False})
convention (str)
passive (bool)
- Returns:
np.ndarray – Euler angles in radians. Size of the array is (3,)
- Return type:
ndarray
- euler_to(other, convention='zyz', passive=False, eps=1e-06)[source]#
Returns the Euler angles that rotate the tensor to the other tensor.
- static make_dipolar(a, i, j, cell=[0, 0, 0], isotopes={}, isotope_i=None, isotope_j=None, rotation_axis=None)[source]#
Create a dipolar NMRTensor
Create a dipolar NMR tensor from an atoms object and the indices of two atoms. Values are in Hz.
Args:a (ase.Atoms): Atoms object of the structure to compute thetensor fori (int): index of first atomj (int): index of second atomcell (np.array): vector of the cell of the second atom, forcouplings between atoms in different cells.By default is [0,0,0].isotopes (dict): dictionary of specific isotopes to use, by elementsymbol. If the isotope doesn’t exist an error willbe raised.isotope_i (int): isotope of atom i. To be used ifdifferent atoms of the same element are supposedto be of different isotopes. If None it will fallback on the previous definitions. Otherwise itoverrides everything else.isotope_j (int): isotope of atom j. See above.rotation_axis (np.array): an axis around which the selected pairis rotating. If present, the tensorwill be averaged for infinitely fastrotation around it.Returns:diptens (NMRTensor): an NMRTensor object with the dipolarcoupling matrix as data.
- rotation_to(other)[source]#
Returns the rotation matrix that rotates the tensor to the other tensor.
- property spherical_repr#
Spherical representation of the tensor
Returns a 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor. The isotropic part is the average of the trace of the tensor, the antisymmetric part is the difference between the tensor and its transpose divided by 2, and the symmetric part is the sum of the tensor and its transpose divided by 2 minus the isotropic part. This construction is such that the sum of the components is equal to the original tensor.
\[\sigma = \sigma_{iso} + \sigma_{A} + \sigma_{S}\]where
\[ \begin{align}\begin{aligned}\sigma_{iso} = \frac{1}{3} \text{Tr}(\sigma) \mathbf{I}\\\sigma_{A} = \frac{1}{2} (\sigma - \sigma^T)\\\sigma_{S} = \frac{1}{2} (\sigma + \sigma^T) - \sigma_{iso}\end{aligned}\end{align} \]- Returns:
np.ndarray – 3x3x3 array containing the isotropic, antisymmetric and symmetric parts of the tensor