Raman spectroscopy

GPAW offers two ways of calculating Raman intensities. One can use the ASE Raman utility together with the GPAW LRTDDFT module as shown in the Resonant Raman tutorial (Resonant-)Raman spectra of gas-phase water.

GPAW also implements Raman spectroscopy for zone-center phonons of extended systems using the electron-phonon coupling (see Electron-phonon coupling) within 3rd order perturbation theory Taghizadeh et a. 1 , which is discussed here. This method is currently only implementated for the LCAO mode.

The Stokes Raman intensity can be written as

\[I(\omega) = I_0 \sum_\nu \frac{n_\nu+1}{\omega_\nu} \vert \sum_{\alpha, \beta} u_{in}^\alpha R_{\alpha \beta}^\nu u_{out}^\beta \vert^2 \delta(\omega-\omega_\nu)\]

where \(\nu\) denotes phonon modes and \(\alpha\), \(\beta\) denote polarisations of the incoming and outgoing laser light. The Raman tensor \(R_{\alpha \beta}^\nu\) has six terms and is given by Ref. 1 Eq. (10)

\[\begin{split}R_{\alpha \beta}^\nu \equiv \sum_{ijmn \mathbf{k}} \left[ \frac{p_{ij}^\alpha (g_{jm}^\nu \delta_{in} - g_{ni}^\nu \delta_{jm})p_{mn}^\beta}{(\hbar \omega_{in}-\varepsilon_{ji})(\hbar \omega_{out}-\varepsilon_{mn})} + \frac{p_{ij}^\alpha (p_{jm}^\beta \delta_{in} - p_{ni}^\beta \delta_{jm})g_{mn}^\nu}{(\hbar \omega_{in}-\varepsilon_{ji})(\hbar \omega_{\nu}-\varepsilon_{mn})} + \\ \frac{p_{ij}^\beta (g_{jm}^\nu \delta_{in} - g_{ni}^\nu \delta_{jm})p_{mn}^\alpha}{(-\hbar \omega_{out}-\varepsilon_{ji})(-\hbar \omega_{in}-\varepsilon_{mn})} + \frac{p_{ij}^\beta (p_{jm}^\alpha \delta_{in} - p_{ni}^\alpha \delta_{jm})g_{mn}^\nu}{(-\hbar \omega_{out}-\varepsilon_{ji})(\hbar \omega_{\nu}-\varepsilon_{mn})} + \\ \frac{g_{ij}^\nu (p_{jm}^\alpha \delta_{in} - p_{ni}^\alpha \delta_{jm})p_{mn}^\beta}{(-\hbar \omega_{\nu}-\varepsilon_{ji})(\hbar \omega_{out}-\varepsilon_{mn})} + \frac{g_{ij}^\nu (p_{jm}^\beta \delta_{in} - p_{ni}^\beta \delta_{jm})p_{mn}^\alpha}{(-\hbar \omega_{\nu}-\varepsilon_{ji})(-\hbar \omega_{in}-\varepsilon_{mn})} \right] f_i(1-f_j)f_n(1-f_m)\end{split}\]

The first term is considered to be the resonant term of the expression, the other terms represent different time orderings of the interaction in the Feynman diagrams.

To compute the Raman intensity we need these ingredients: The momentum matrix elements \(p_{ij}^\alpha=\langle i \mathbf{k} | \hat p^\alpha| j \mathbf{k} \rangle\), the electron-phonon matrix \(g_{ij}^\nu = \langle i \mathbf{k} \vert \partial_{\nu{q=0}} V^{KS} \vert j \mathbf{k} \rangle\) in the optical limit \(\mathbf{q}=0\) and of course knowledge of the electronic states and phonon modes throughout the Brillouin zone. For these calculations we can employ in the get_momentum_transitions() method, the GPAW electron-phonon module Electron-phonon coupling and the ASE phonon module, respectively.

The get_momentum_transitions() method is currently not aware of symmetry. It is therefore required to switch off point-group symmetry in GPAW, so that matrix elements for all k-points and not just the irreducible ones are calculated. The momentum matrix elements transform as \(p \rightarrow -p^*\) under time-reversal symmetry \(\mathbf{k} \rightarrow -\mathbf{k}\). Accordingly, it is possible to enable time-reversal symmetry. By default the routine saves a file called mom_skvnm.npy containing the momentum matrix. This can be deactivated using the savetofile switch. The matrix is always the return value of get_momentum_transitions().

Energy changes for phonons and potential changes for electron-phonon couplings are both computed using a finite displacement technique. Both quantities can be obtained simultaenously. In principle the phonon modes can be obtained in any fashion, which yields an ASE phonon object though. For small systems the finite displacement method has the disadvantage of leading to an interaction of a displaced atom with its periodic images. This can lead to large errors especially in the electron-phonon matrix. This can be avoided by using a sufficiently large supercell for the finite displacement simulations.

If phonon and effective potential are calculated simultaenously, results are saved in the same cache files with default name \(elph\).

The Raman module offers a wrapper, EPC(), around the ElectronPhononCoupling() class to facilitate easier computation of the electron-phonon matrix.

Some more details are elaborated in the related tutorial Raman spectroscopy for extended systems.

References

1(1,2)

A. Taghizadeh, U. Leffers, T.G. Pedersen, K.S. Thygesen, “A library of ab initio Raman spectra for automated identification of 2D materials”, Nature Communications 11, 3011 (2020).

Code

gpaw.raman.dipoletransition.get_momentum_transitions(wfs, savetofile: bool = True) numpy.ndarray[source]

Finds the momentum matrix elements: <nk|p|mk> = k delta_nm - i <nk|nabla|mk>

Parameters
  • wfs – LCAO WaveFunctions object

  • savetofile (bool) – Determines whether matrix is written to the file mom_skvnm.npy (default=True)

class gpaw.raman.elph.EPC(atoms, calc=None, supercell=(1, 1, 1), name='elph', delta=0.01, calculate_forces=False)[source]

Modify ElectronPhononCoupling to fit Raman stuff better.

Primarily, always save the supercell matrix in separate files, so that calculation of the elph matrix can be better parallelised.

Initialize with base class args and kwargs.

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

  • 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.

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

gpaw.raman.elph.EPC.calculate_supercell_matrix(self, calc, name=None, filter=None, include_pseudo=True)

Calculate elph supercell matrix.

This is a necessary intermediary step before calculating the electron- phonon matrix.

The Raman version always uses dump=2 when calling the ElectronPhononCoupling routine.

Parameters
  • calc (GPAW) – Converged ground-state calculation. Same grid as before.

  • name (str) – User specified name of the generated pickle file(s). If not provided, the string in the name attribute is used.

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

  • include_pseudo (bool) – Include the contribution from the psedupotential in the atomic gradients. If False, only the gradient of the effective potential is included (default: True).

gpaw.raman.elph.EPC.get_elph_matrix(self, calc, phonon, savetofile=True) numpy.ndarray

Calculate the electronphonon matrix in Bloch states.

Always uses q=0.

Parameters
  • calc (GPAW) – Converged ground-state calculation. NOT supercell.

  • phonon (ase.phonons.Phonons) – Phonon object

  • savetofile (bool) – Switch for saving to gsqklnn.npy file

gpaw.raman.raman.calculate_raman(calc, w_ph, w_in, d_i, d_o, resonant_only=False, gamma_l=0.1, momname='mom_skvnm.npy', elphname='gsqklnn.npy', gridspacing=1.0, suffix=None)[source]

Calculates the first order Raman tensor contribution for a given polarization.

Parameters
  • calc (GPAW) – Converged ground state calculation

  • w_ph (ndarray, str) – Array of phonon frequencies in eV, or name of file with them

  • w_in (float) – Laser energy in eV

  • d_i (int) – Incoming polarization

  • d_o (int) – Outgoing polarization

  • gamma_l (float) – Line broadening in eV

  • momname (str) – Name of momentum file

  • elphname (str) – Name of electron-phonon file

  • gridspacing (float) – grid spacing in cm^-1

  • suffix (str) – suffix for Rlab file name. Defaults to XXXnm if not given.

gpaw.raman.raman.calculate_raman_intensity(w_ph, d_i, d_o, suffix, T=300)[source]

Calculates Raman intensities from Raman tensor.

Parameters
  • w_ph (ndarray, str) – Array of phonon frequencies in eV, or name of file with them

  • d_i (int) – Incoming polarization

  • d_o (int) – Outgoing polarization

  • suffix (str) – Suffix for Rlab and RI files. Usually XXXnm

gpaw.raman.raman.plot_raman(figname, RIsuffix, relative=True, w_min=None, w_max=None)[source]

Plots a given Raman spectrum.

Parameters
  • figname (str) – Filename for figure.

  • RIsuffix (str, list) – Suffix of Raman intensity files to use for plotting. For example “0_1_455nm” for RI_0_1_455nm.npy