Electron-Phonon Coupling Theory

Phonons can interact with the electrons in a variety of ways. For example, when an electron moves through the crystal, it can scatter off of a phonon, thereby transferring some of its energy to the lattice. Conversely, when a phonon vibrates, it can create an oscillating electric field that can interact with the electrons, inducing a change in their energies and momenta. The coupling between electrons and lattice vibrations is responsible for a range of interesting and important phenomena, from electrical and thermal conductivity to superconductivity.

The first order electron-phonon coupling matrix \(g_{mn}^\nu(\mathbf{k}, \mathbf{q})\) couples the electronic states \(m(\mathbf{k}+ \mathbf{q}),n(\mathbf{k})\) via phonons \(\nu\) at wave vectors \(\mathbf{q}\) and frequencies \(\omega_\nu\):[1]

\[g_{mn}^\nu(\mathbf{k}, \mathbf{q}) = \sqrt{ \frac{\hbar}{2 m_0 \omega_\nu}} M_{mn}^\nu(\mathbf{k}, \mathbf{q}) .\]

with

\[M_{mn}^\nu(\mathbf{k}, \mathbf{q}) = \langle \psi_{m \mathbf{k}+ \mathbf{q}} \vert \nabla_u V^{KS} \cdot \mathbf{e}_\nu \vert \psi_{n\mathbf{k}} \rangle.\]

Here \(m_0\) is the sum of the masses of all the atoms in the unit cell and \(\nabla_u\) denotes the gradient with respect to atomic displacements. For the three translational modes at \(\vert \mathbf{q} \vert = 0\) the matrix elements \(g_{mn}^\nu = 0\), as a consequence of the acoustic sum rule.

Implementation

Within the PAW framework to Kohn-Sham potential can be split into a local part \(V(\mathbf{r})\) represented on a regular grid and a nonlocal part \(\Delta H^a_{i_1 i_2}\):

\[V^{KS} = V + \Delta H^a_{i_1 i_2}.\]

In GPAW \(\nabla_u V^{KS}(\mathbf{r})\) is determined using the finite difference method in a supercell. The potential at the displaced coordinates is computed by the DisplacementRunner() class, which is based on ASEs ase.phonon.Displacement class. The central difference derivative, as evaluted in the Supercell() class, consists of four contributions:

\[\langle \psi_{i} \vert \nabla_u V^{KS} \vert \psi_{j} \rangle = \langle \tilde \psi_{i} \vert \nabla_u V(\mathbf{r}) \vert \tilde \psi_{j} \rangle + \sum_{a, ij} \langle \tilde \psi_{i} \vert \tilde p^a_i \rangle (\nabla_u \Delta H^a_{i_1 i_2}) \langle \tilde p^a_j \vert \tilde \psi_{j} \rangle + \sum_{a, ij} \langle \tilde \psi_{i} \vert \nabla_u \tilde p^a_i \rangle \Delta H^a_{i_1 i_2} \langle \tilde p^a_j \vert \tilde \psi_{j} \rangle + \sum_{a, ij} \langle \tilde \psi_{i} \vert \tilde p^a_i \rangle \Delta H^a_{i_1 i_2} \langle \sum_{a, ij} \langle \tilde \psi_{i} \vert \nabla_u \tilde p^a_i \rangle \Delta H^a_{i_1 i_2} \langle \tilde p^a_j \vert \tilde \psi_{j} \rangle \tilde p^a_j \vert \tilde \psi_{j} \rangle\]

Here we do not project the derivatives onto electronics states actually, not rather onto LCAO orbitals \(\Psi_{NM}\), where \(N\) denotes the cell index and \(M\) the orbital index. We use a The Fourier transform from the \(\mathbf{k}\)-space Bloch to the real space representation so that we can can later to compute \(M_{mn}^\nu\) for arbitrary \(\mathbf{q}\):

\[\begin{split}\mathbf{g}_{\substack{N M\\ N^\prime M^\prime}}^{sc} = FFT\left[ \langle \Psi_{NM}(\mathbf{k}) \vert \nabla_u V^{KS} \vert \Psi_{N^\prime M^\prime}(\mathbf{k}) \rangle \right].\end{split}\]

Finally, the electron-phonon coupling matrix is obtained by projecting the supercell matrix into the primitive unit cell bands \(m, n\) and phonon modes \(\nu\) in ElectronPhononMatrix():

\[\begin{split}M_{mn}^\nu(\mathbf{k}, \mathbf{q}) = \sum_{\substack{N M\\ N^\prime M^\prime}} C_{mM}^{\star} C_{nM^\prime} \mathbf{g}_{\mu}^{sc} \cdot \mathbf{u}_{q \nu} e^{2\pi i [(\mathbf{k}+\mathbf{q})\cdot \mathbf{R_N} - \mathbf{k}\cdot \mathbf{R_N^\prime}]},\end{split}\]

where \(C_{nM}\) are the LCAO coefficients and \(\mathbf{u}_{q \nu}\) are the mass-scaled phonon displacement vectors.

Checkout Electron-phonon coupling for an example and exerice and Raman spectroscopy and Raman spectroscopy for extended systems for an application.

References

Code

class gpaw.elph.DisplacementRunner(atoms, calc, supercell=(1, 1, 1), name='elph', delta=0.01, calculate_forces=True)[source]

Class for calculating the changes in effective potential.

The derivative of the effective potential wrt atomic displacements is obtained from a finite difference approximation to the derivative by doing a self-consistent calculation for atomic displacements in the +/- directions. These calculations are carried out in the run member function.

Initialize with base class args and kwargs.

Parameters:
  • atoms (Atoms) – The atoms to work on. Primitive cell.

  • calc (GPAW) – Calculator for the supercell finite displacement calculation.

  • supercell (tuple, list) – Size of supercell given by the number of repetitions (l, m, n) of the small unit cell in each direction.

  • name (str) – Name to use for files (default: ‘elph’).

  • delta (float) – Magnitude of displacements. (default: 0.01 A)

  • calculate_forces (bool) – If true, also calculate and store the dynamical matrix.

run()[source]

Run the calculations for the required displacements.

class gpaw.elph.Supercell(atoms, supercell_name='supercell', supercell=(1, 1, 1), indices=None)[source]

Class for supercell-related stuff.

Initialize supercell class.

Parameters:
  • atoms (Atoms) – The atoms to work on. Primitive cell.

  • supercell_name (str) – User specified name of the generated JSON cache. Default is ‘supercell’.

  • supercell (tuple) – Size of supercell given by the number of repetitions (l, m, n) of the small unit cell in each direction.

classmethod calculate_gradient(fd_name, indices=None)[source]

Calculate gradient of effective potential and projector coefs.

This function loads the generated json files and calculates finite-difference derivatives.

Parameters:

fd_name (str) – name of finite difference JSON cache

calculate_supercell_matrix(calc, fd_name='elph', filter=None)[source]

Calculate matrix elements of the el-ph coupling in the LCAO basis.

This function calculates the matrix elements between LCAOs and local atomic gradients of the effective potential. The matrix elements are calculated for the supercell used to obtain finite-difference approximations to the derivatives of the effective potential wrt to atomic displacements.

The resulting g_xsNNMM is saved into a JSON cache.

Parameters:
  • calc (GPAW) – LCAO calculator for the calculation of the supercell matrix.

  • fd_name (str) – User specified name of the finite difference JSON cache. Default is ‘elph’.

  • filter (str) – Fourier filter atomic gradients of the effective potential. The specified components (normal or umklapp) are removed (default: None).

classmethod load_supercell_matrix(name='supercell')[source]

Load supercell matrix from cache.

Parameters:

name (str) – User specified name of the cache.

set_basis_info(*args)[source]

Store LCAO basis info for atoms in reference cell in attribute.

Parameters:

args (tuple) – If the LCAO calculator is not available (e.g. if the supercell is loaded from file), the load_supercell_matrix member function provides the required info as arguments.

class gpaw.elph.ElectronPhononMatrix(atoms, supercell_cache, phonon, load_sc_as_needed=True, indices=None)[source]

Class for containing the electron-phonon matrix

Initialize with base class args and kwargs.

Parameters:
  • atoms (Atoms) – Primitive cell object

  • supercell_cache (str) – Name of JSON cache containing supercell matrix

  • phonon (str, dict, Phonons) – Can be either name of phonon cache generated with electron-phonon DisplacementRunner or dictonary of arguments used in Phonons run or Phonons object.

  • load_sc_as_needed (bool) – Load supercell matrix elements only as needed. Greatly reduces memory requirement for large systems, but introduces huge filesystem overhead

  • indices (list) – List of atoms (indices) to use. Default: Use all.

bloch_matrix(calc, k_qc=None, savetofile=True, prefactor=True, accoustic=True)[source]

Calculate el-ph coupling in the Bloch basis for the electrons.

This function calculates the electron-phonon coupling between the specified Bloch states, i.e.:

           ______
 mnl      / hbar               ^
g    =   /-------  < m k + q | e  . grad V  | n k >
 kq    \/ 2 M w                 ql        q
               ql

In case the prefactor=False is given, the bare matrix element (in units of eV / Ang) without the sqrt prefactor is returned.

Parameters:
  • calc (GPAW) – Converged calculator object containing the LCAO wavefuntions (don’t use point group symmetry)

  • k_qc (np.ndarray) – q-vectors of the phonons. Must only contain values comenserate with k-point sampling of calculator. Default: all kpoints used.

  • savetofile (bool) – If true (default), saves matrix to gsqklnn.npy

  • prefactor (bool) – if false, don’t multiply with sqrt prefactor (Default: True)

  • accoustic (bool) – if True, for 3 accoustic modes set g=0 for q=0 (Default: True)