pymatgen.analysis.nmr module

class ChemicalShielding[source]

Bases: pymatgen.core.tensors.SquareTensor

This class extends the SquareTensor to perform extra analysis unique to NMR Chemical shielding tensors

Three notations to describe chemical shielding tensor (RK Harris; Magn. Reson. Chem. 2008, 46, 582–598; DOI: 10.1002/mrc.2225) are supported.

Authors: Shyam Dwaraknath, Xiaohui Qu

Create a Chemical Shielding tensor. Note that the constructor uses __new__ rather than __init__ according to the standard method of subclassing numpy ndarrays.

Parameters:
  • cs_matrix (1x3 or 3x3 array-like) – the 3x3 array-like representing the chemical shielding tensor or a 1x3 array of the primary sigma values corresponding to the principal axis system
  • vscale (6x1 array-like) – 6x1 array-like scaling the voigt-notation vector with the tensor entries
HaeberlenNotation

alias of HaeberlenNotion

class MarylandNotation(sigma_iso, omega, kappa)

Bases: tuple

Create new instance of MarylandNotation(sigma_iso, omega, kappa)

kappa

Alias for field number 2

omega

Alias for field number 1

sigma_iso

Alias for field number 0

class MehringNotation(sigma_iso, sigma_11, sigma_22, sigma_33)

Bases: tuple

Create new instance of MehringNotation(sigma_iso, sigma_11, sigma_22, sigma_33)

sigma_11

Alias for field number 1

sigma_22

Alias for field number 2

sigma_33

Alias for field number 3

sigma_iso

Alias for field number 0

classmethod from_maryland_notation(sigma_iso, omega, kappa)[source]
haeberlen_values

the Chemical shielding tensor in Haeberlen Notation

Type:Returns
maryland_values

the Chemical shielding tensor in Maryland Notation

Type:Returns
mehring_values

the Chemical shielding tensor in Mehring Notation

Type:Returns
principal_axis_system

Returns a chemical shielding tensor aligned to the principle axis system so that only the 3 diagnol components are non-zero

class ElectricFieldGradient[source]

Bases: pymatgen.core.tensors.SquareTensor

This class extends the SquareTensor to perform extra analysis unique to NMR Electric Field Gradient tensors in units of V/Angstrom^2

Authors: Shyam Dwaraknath, Xiaohui Qu

Create a Chemical Shielding tensor. Note that the constructor uses __new__ rather than __init__ according to the standard method of subclassing numpy ndarrays.

Parameters:
  • efg_matrix (1x3 or 3x3 array-like) – the 3x3 array-like representing the electric field tensor or a 1x3 array of the primary values corresponding to the principal axis system
  • vscale (6x1 array-like) – 6x1 array-like scaling the voigt-notation vector with the tensor entries
V_xx
V_yy
V_zz
asymmetry

(V_yy - V_xx)/V_zz

Type:Asymmetry of the electric field tensor defined as
coupling_constant(specie)[source]
Computes the couplling constant C_q as defined in:
Wasylishen R E, Ashbrook S E, Wimperis S. NMR of quadrupolar nuclei in solid materials[M]. John Wiley & Sons, 2012. (Chapter 3.2)
C_q for a specific atom type for this electric field tensor:
C_q=e*Q*V_zz/h

h: planck’s constant Q: nuclear electric quadrupole moment in mb (millibarn e: elementary proton charge

Parameters:specie – flexible input to specify the species at this site. Can take a isotope or element string, Specie object, or Site object
Returns:the coupling constant as a FloatWithUnit in MHz
principal_axis_system

Returns a electric field gradient tensor aligned to the principle axis system so that only the 3 diagnol components are non-zero