Phonon

Mass-spring models from Quantum ESPRESSO.

class elphmod.ph.Model(flfrc=None, quadrupole_fmt=None, apply_asr=False, apply_asr_simple=False, apply_zasr=False, apply_rsr=False, lr=True, lr2d=None, amass=None, at=None, tau=None, atom_order=None, alph=None, epsil=None, zeu=None, Q=None, L=None, perp=True, divide_mass=True, divide_ndegen=True, ifc=None)[source]

Mass-spring model for the phonons.

Parameters:
flfrcstr

File with interatomic force constants from q2r.x or common filename part (without appended/inserted q-point number) of dynamical matrices from ph.x. Both the plain-text and XML format are supported.

quadrupole_fmtstr

File with quadrupole tensors in format suitable for epw.x.

apply_asrbool

Apply acoustic sum rule correction to force constants?

apply_asr_simplebool

Apply simple acoustic sum rule correction to force constants? This sets the self force constant to minus the sum of all other force constants.

apply_zasrbool

Apply acoustic sum rule correction to Born effective charges?

apply_rsrbool

Apply rotation sum rule correction to force constants?

lrbool

Compute long-range terms in case of polar material?

lr2dbool

Compute long-range terms for two-dimensional system if lr?

amassndarray

Atomic masses if flfrc is omitted.

atndarray

Bravais lattice vectors if flfrc is omitted.

taundarray

Positions of basis atoms if flfrc is omitted.

atom_orderlist of str

Ordered list of atoms if flfrc is omitted.

alphfloat

Ewald parameter if flfrc is omitted.

epsilndarray

Dielectric tensor if flfrc is omitted.

zeundarray

Born effective charges if flfrc is omitted.

Qndarray

Quadrupole tensors if quadrupole_fmt is omitted.

Lfloat

Range-separation parameter for two-dimensional electrostatics.

perpbool

Yield out-of-plane long-range terms if L is nonzero?

divide_massbool

Divide force constants and Born effective charges by atomic masses?

divide_ndegenbool

Divide force constants by degeneracy of Wigner-Seitz point? Only True yields correct phonons. False should only be used for debugging.

ifcstr

Name of file with interatomic force constants to be created if flfrc is prefix of files with dynamical matrices. Used to emulate q2r.x.

Attributes:
Rndarray

Lattice vectors \(\vec R\) of Wigner-Seitz supercell.

datandarray

Corresponding force constants divided by atomic masses in Ry2.

\[D_{\vec R i j} = \frac{\hbar^2}{\sqrt{M_i M_j}} \frac{\partial^2 E}{\partial u_{0 i} \partial u_{\vec R j}}\]

If not divide_mass, the prefactor \(\hbar^2 / \sqrt{M_i M_j}\) is absent and the units are Ry/bohr2 instead. If lr, this is only the short-range component of the force constants.

Mndarray

Atomic masses \(M_i\).

andarray

Bravais lattice vectors.

rndarray

Positions of basis atoms.

lndarray.

Bond lengths.

atom_orderlist of str

Ordered list of atoms.

alphafloat

Ewald parameter.

epsndarray

Dielectric tensor.

Zndarray

Born effective charges.

zndarray

Born effective charges divided by square root of atomic masses.

Qndarray

Quadrupole tensors.

qndarray

Quadrupole tensors divided by square root of atomic masses.

Lfloat

Range-separation parameter for two-dimensional electrostatics.

perpbool

Yield out-of-plane long-range terms if L is nonzero?

divide_massbool

Have force constants and Born charges been divided by atomic masses?

divide_ndegenbool

Have force constants been divided by degeneracy of Wigner-Seitz point?

sizeint

Number of displacement directions/bands.

natint

Number of atoms.

nqtuple of int

Shape of original q-point mesh.

q0ndarray

Original q-point mesh.

D0ndarray

Dynamical matrices on original q-point mesh.

cellslist of tuple of int, optional

Lattice vectors of unit cells if the model describes a supercell.

Nlist of tuple of int, optional

Primitive vectors of supercell if the model describes a supercell.

lrbool

Compute long-range terms?

lr2dbool

Compute long-range terms for two-dimensional system if lr?

bndarray

Reciprocal primitive vectors if lr.

prefactorfloat

Prefactor of long-range terms if lr.

r_efffloat

Effective thickness if lr2d.

scalefloat

Relevant scaling factor if lr.

D0_lrndarray

Constant part of long-range correction to dynamical matrix if lr.

C(R1=0, R2=0, R3=0)[source]

Get interatomic force constants for arbitrary lattice vector.

Parameters:
R1, R2, R3int, default 0

Lattice vector in units of primitive vectors.

Returns:
ndarray

Element of data or zero.

D(q1=0, q2=0, q3=0)[source]

Set up dynamical matrix for arbitrary q point.

Parameters:
q1, q2, q3float, default 0.0

q point in crystal coordinates with period \(2 \pi\).

Returns:
ndarray

Fourier transform of data, possibly plus a long-range term.

D_lr(q1=0, q2=0, q3=0)[source]

Calculate long-range part of dynamical matrix.

clear()[source]

Delete all lattice vectors and associated matrix elements.

decay()[source]

Plot force constants as a function of bond length.

Returns:
ndarray

Bond lengths.

ndarray

Frobenius norm of force-constant matrices.

generate_long_range(q1=0, q2=0, q3=0, perp=None, eps=1e-10)[source]

Generate long-range terms.

Parameters:
qndarray

q point in reciprocal lattice units \(q_i \in [0, 2 \pi)\).

perpbool

Yield out-of-plane terms? Defaults to attribute perp.

epsfloat

Tolerance for vanishing lattice vectors.

Returns:
ndarray

Scalar prefactor for relevant reciprocal lattice vectors.

ndarray

Direction-dependent term for relevant reciprocal lattice vectors.

order_atoms(*order)[source]

Reorder atoms.

Together with shift_atoms(), this function helps reconcile inconsistent definitions of the basis/motif of the Bravais lattice.

Parameters:
*orderint

New atom order.

prepare_long_range(G_max=28.0, G_2d=False)[source]

Prepare calculation of long-range terms for polar materials.

The following two routines are based on rgd_blk and rgd_blk_epw of Quantum ESPRESSO and the EPW code.

Copyright (C) 2010-2016 S. Ponce’, R. Margine, C. Verdi, F. Giustino Copyright (C) 2001-2012 Quantum ESPRESSO group

Please refer to:

  • Phonons: Gonze et al., Phys. Rev. B 50, 13035 (1994)

  • Coupling: Verdi and Giustino, Phys. Rev. Lett. 115, 176401 (2015)

  • 2D case: Sohier, Calandra, and Mauri, Phys. Rev. B 94, 085415 (2016)

  • Quadrupoles: Ponce’ et al., Phys. Rev. Research 3, 043022 (2021)

  • Improved 2D phonons: Royo and Stengel, Phys. Rev. X 11, 041027 (2021)

  • Improved 2D coupling: Ponce’ et al., Phys. Rev. B 107, 155424 (2023)

Parameters:
G_maxfloat

Cutoff for reciprocal lattice vectors.

G_2dbool, default False

Do not sample out-of-plane reciprocal lattice vectors regardless of lr2d?

sample_orig()[source]

Sample dynamical matrix on original q-point mesh.

shift_atoms(s, S)[source]

Move selected atoms across unit-cell boundary.

Together with order_atoms(), this function helps reconcile inconsistent definitions of the basis/motif of the Bravais lattice.

Parameters:
sslice

Slice of atom indices corresponding to selected basis atom(s).

Stuple of int

Shift of as multiple of primitive lattice vectors.

standardize(eps=0.0, symmetrize=False)[source]

Standardize mass-spring data.

  • Keep only nonzero force-constant matrices.

  • Sum over repeated lattice vectors.

  • Sort lattice vectors.

  • Optionally symmetrize force constants:

\[D_{\vec q} = D_{\vec q}^\dagger, D_{\vec R} = D_{-\vec R}^\dagger\]
sum_force_constants()[source]

Calculate sum of absolute values of short-range force constants.

For the optimal range-separation parameter this value is minimal.

supercell(N1=1, N2=1, N3=1, sparse=False)[source]

Map mass-spring model onto supercell.

Parameters:
N1, N2, N3tuple of int or int, default 1

Supercell lattice vectors in units of primitive lattice vectors.

sparsebool, default False

Only calculate q = 0 dynamical matrix as a sparse matrix to save memory? The result is stored in the attribute Ds. Consider using standardize() with nonzero eps and symmetrize before.

Returns:
object

Mass-spring model for supercell.

symmetrize()[source]

Symmetrize dynamical matrix.

to_flfrc(flfrc, nr1=None, nr2=None, nr3=None)[source]

Save mass-spring model to force-constants file.

Parameters:
flfrcstr

Filename.

nr1, nr2, nr3int, optional

Mesh dimensions. If not given, nq is assumed.

to_pwi(pwi, **kwargs)[source]

Save atomic positions etc. to PWscf input file.

Parameters:
pwistr

Filename.

**kwargs

Keyword arguments with further parameters to be written.

unit_cell()[source]

Map mass-spring model back to unit cell.

See also

supercell
update_short_range(flfrc=None, apply_asr_simple=False)[source]

Update short-range part of interatomic force constants.

Parameters:
flfrcstr

Filename where short-range force constants are written.

apply_asr_simplebool

Apply simple acoustic sum rule correction to short-range force constants?

elphmod.ph.asr(phid)[source]

Apply simple acoustic sum rule correction to force constants.

elphmod.ph.divide_by_mass(dyn, amass)[source]

Divide dynamical matrices by atomic masses.

Parameters:
dynndarray

Dynamical matrices.

amasslist of float

Atomic masses.

elphmod.ph.fildyn_freq(fildyn='matdyn')[source]

Create fildyn.freq as created by Quantum ESPRESSO’s ph.x.

Parameters:
fildynstr, default ‘matdyn’

Prefix of files with dynamical matrices.

elphmod.ph.group(n, size=3)[source]

Create slice of dynamical matrix belonging to n-th atom.

elphmod.ph.interpolate_dynamical_matrices(D, q, nq, fildyn_template, fildyn, flfrc, write_fildyn0=True, apply_asr=False, apply_asr_simple=False, apply_rsr=False, divide_mass=True, qe_prefix='', clean=False)[source]

Interpolate dynamical matrices given for irreducible wedge of q points.

This function still uses the Quantum ESPRESSO executables q2qstar.x and q2r.x. They are called in serial by each MPI process, which leads to problems if they have been compiled for parallel execution. If you want to run this function in parallel, you have two choices:

  1. Configure Quantum ESPRESSO for compilation of serial executables via ./configure --disable-parallel and run make ph. If you do not want to make them available through the environmental variable PATH, you can also set the parameter qe-prefix to '/path/to/serial/q-e/bin/'. The trailing slash is required.

  2. If your MPI implementation supports nested calls to mpirun, you may try to set qe_prefix to 'mpirun -np 1 '. The trailing space is required.

Parameters:
Dlist of square arrays

Dynamical matrices for all irreducible q points.

qlist of 2-tuples

Irreducible q points in crystal coordinates with period \(2 \pi\).

nqint

Number of q points per dimension, i.e., size of uniform mesh.

fildyn_templatestr

Complete name of fildyn file from which to take header information.

fildynstr

Prefix for written files with dynamical matrices.

flfrcstr

Name of written file with interatomic force constants.

write_fildyn0bool

Write fildyn0 needed by q2r.x? Otherwise the file must be present.

apply_asrbool

Enforce acoustic sum rule by overwriting self force constants?

apply_asr_simplebool

Apply simple acoustic sum rule correction to force constants? This sets the self force constant to minus the sum of all other force constants.

apply_rsrbool

Enforce rotation sum rule by overwriting self force constants?

divide_massbool

Have input dynamical matrices been divided by atomic mass? This is independent of ph.divide_mass, which is always respected.

qe_prefixstr

String to prepend to names of Quantum ESPRESSO executables.

cleanbool

Delete all temporary files afterwards?

Returns:
function

Fourier-interpolant (via force constants) for dynamical matrices.

elphmod.ph.polarization(e, path, angle=60)[source]

Characterize as in-plane longitudinal/transverse or out-of-plane.

elphmod.ph.q2r(ph, D_irr=None, q_irr=None, nq=None, D_full=None, angle=60, apply_asr=False, apply_asr_simple=False, apply_rsr=False, flfrc=None, divide_mass=True)[source]

Interpolate dynamical matrices given for irreducible wedge of q points.

This function replaces interpolate_dynamical_matrices, which depends on Quantum ESPRESSO. For 2D lattices, it is sufficient to provide dynamical matrices D_irr for the irreducible q points q_irr. Here, for the square lattice, the rotation symmetry (90 degrees) is currently disabled! In turn, for 1D and 2D lattices, dynamical matrices D_full on the complete uniform q-point mesh must be given.

Parameters:
phModel

Mass-spring model.

D_irrlist of square arrays

Dynamical matrices for all irreducible q points.

q_irrlist of 2-tuples

Irreducible q points in crystal coordinates with period \(2 \pi\).

nqint or tuple of int

Number of q points per dimension, i.e., size of uniform mesh. Different numbers of q points along different axes can be specified via a tuple. Alternatively, nq is inferred from the shape of D_full.

D_fullndarray

Dynamical matrices on complete uniform q-point mesh.

anglefloat

Angle between mesh axes in degrees.

apply_asrbool

Enforce acoustic sum rule by overwriting self force constants?

apply_asr_simplebool

Apply simple acoustic sum rule correction to force constants? This sets the self force constant to minus the sum of all other force constants.

apply_rsrbool

Enforce rotation sum rule by overwriting self force constants?

flfrcstr

Filename. If given, save file with force constants as generated by q2r.x. This is done before acoustic and rotation sum are applied.

divide_massbool

Have input dynamical matrices been divided by atomic mass? This is independent of ph.divide_mass, which is always respected.

elphmod.ph.read_flfrc(flfrc)[source]

Read force constants in real or reciprocal space.

flfrc can be either the force-constants file generated by q2r.x or one of the dynamical-matrix files generated by ph.x.

Returns:
ndarray or tuple of list of ndarray

Force constants for input generated by q2r.x or equivalent q points and corresponding dynamical matrices for input generated by ph.x.

ndarray

Atomic masses.

ndarray

Bravais lattice vectors.

ndarray

Atomic positions/basis vectors.

list of str

Atomic symbols.

float

Ewald parameter.

ndarray

Dielectric tensor.

ndarray

Born effective charges.

See also

write_flfrc

Inverse function.

elphmod.ph.read_flfrc_xml(flfrc)[source]

Read force constants in real or reciprocal space from XML files.

See also

read_flfrc

Equivalent function for non-XML files.

elphmod.ph.read_q(fildyn0)[source]

Read list of irreducible q points from fildyn0.

elphmod.ph.read_quadrupole_fmt(quadrupole_fmt)[source]

Read file quadrupole.fmt suitable for epw.x.

Parameters:
quadrupole_fmtstr

Filename.

Returns:
Qndarray

Quadrupole tensor.

See also

write_quadrupole_fmt

Inverse function.

elphmod.ph.sgnsqrt(w2)[source]

Calculate signed square root.

elphmod.ph.spectral_function(D, omega, eta)[source]

Calculate phonon spectral function.

See Eq. 76 by Monacelli et al., J. Phys. Condens. Matter 33, 363001 (2021). We use an additional prefactor of 2 to normalize the frequency integral (from zero to infinity) of the phonon spectral function to the number of modes (for each q point).

Parameters:
Dndarray

Dynamical matrix. The last three axes belong to the two displacement directions and the frequency argument.

omegandarray

Frequency arguments.

etafloat

Broadening parameter for all phonon modes.

Returns:
ndarray

Spectral function. The first axes belongs to the frequency argument.

elphmod.ph.sum_rule_correction(ph, asr=True, rsr=True, eps=1e-15, report=True)[source]

Apply sum rule correction to force constants.

Unlike asr() called by the argument apply_asr_simple, the corrected force constants may not fulfill the point-symmetries of the crystal anymore. In turn, the total deviation from the original force constants is minimized.

Parameters:
phModel

Mass-spring model for the phonons.

asrbool

Enforce acoustic sum rule?

rsrbool

Enforce Born-Huang rotation sum rule?

epsfloat

Smallest safe absolute value of divisor.

reportbool

Print sums before and after correction?

elphmod.ph.write_flfrc(flfrc, phid, amass, at, tau, atom_order, alph=None, epsil=None, zeu=None)[source]

Write force constants in real or reciprocal space.

Parameters:
flfrcstr

Filename.

phidndarray or tuple of list of ndarray

Force constants or equivalent q points and corresponding dynamical matrices.

amassndarray

Atomic masses.

atndarray

Bravais lattice vectors.

taundarray

Atomic positions/basis vectors.

atom_orderlist of str

Atomic symbols.

alphfloat

Ewald parameter.

epsilndarray

Dielectric tensor.

zeundarray

Born effective charges.

See also

read_flfrc

Inverse function.

elphmod.ph.write_q(fildyn0, q, nq)[source]

Write list of irreducible q points to fildyn0.

elphmod.ph.write_quadrupole_fmt(quadrupole_fmt, Q)[source]

Write file quadrupole.fmt suitable for epw.x.

Parameters:
quadrupole_fmtstr

Filename.

Qndarray

Quadrupole tensor.

See also

read_quadrupole_fmt

Inverse function.

elphmod.ph.zasr(Z)[source]

Apply acoustic sum rule correction to Born effective charges.