Isotropic and anisotropic hyperfine coupling paramters

Python API and CLI

Use the hyperfine_parameters() function or the CLI tool:

$ python3 -m gpaw.hyperfine --help
gpaw.hyperfine.hyperfine_parameters(calc, exclude_core=False)[source]

Calculate isotropic and anisotropic hyperfine coupling parameters.

One tensor (\(A_{ij}\)) per atom is returned in eV units. In Hartree atomic units, we have the isotropic part \(a = \text{Tr}(\mathbf{A}) / 3\):

\[a = \frac{2 \alpha^2 g_e m_e}{3 m_p} \int \delta_T(\mathbf{r}) \rho_s(\mathbf{r}) d\mathbf{r},\]

and the anisotropic part \(\mathbf{A} - a\):

\[\frac{\alpha^2 g_e m_e}{4 \pi m_p} \int \frac{3 r_i r_j - \delta_{ij} r^2}{r^5} \rho_s(\mathbf{r}) d\mathbf{r}.\]

Remember to multiply each tensor by the g-factors of the nuclei.

Use exclude_core=True to exclude contribution from “frozen” core.

For details, see Peter E. Blöchl and Oleg V. Yazyev et al..

The results should be divided by the net mangetic moments averaged over an entire supercell. If one wants to calculate the localized HF effects on an atom or group of atoms in an anti-ferromagnetic material, one needs to divide the HF constants with the average magnetic moment of that atom or group of atoms. As, an anti-ferromagnetic system, overall should have a net magnetic moment of zero.

G-factors

Here is a list of g-factors (from Wikipedia):

Nucleus

g-factor

\(^{1}\)H

5.586

\(^{3}\)He

-4.255

\(^{7}\)Li

2.171

\(^{13}\)C

1.405

\(^{14}\)N

0.404

\(^{17}\)O

-0.757

\(^{19}\)F

5.254

\(^{23}\)Na

1.477

\(^{27}\)Al

1.457

\(^{29}\)Si

-1.111

\(^{31}\)P

2.261

\(^{57}\)Fe

0.181

\(^{63}\)Cu

1.485

\(^{67}\)Zn

0.350

\(^{129}\)Xe

-1.545

Hydrogen 21 cm line

Here is how to calculate the famous hydrogen spectral line of 21 cm:

import numpy as np
from ase import Atoms
import ase.units as units
from gpaw import GPAW, PW
from gpaw.hyperfine import hyperfine_parameters

h = Atoms('H', magmoms=[1])
h.center(vacuum=3)
h.calc = GPAW(mode=PW(400), txt=None)
e = h.get_potential_energy()
A = hyperfine_parameters(h.calc)[0] * 5.586
a = np.trace(A) / 3
frequency = a * units._e / units._hplanck  # Hz
wavelength = units._c / frequency  # meters
print(f'{wavelength * 100:.1f} cm')

The output will be 23.2 cm. It’s slightly off because the LDA spin-density at the position of the hydrogen nucleus is a bit too low (should be \(1/\pi\) in atomic units).