Input File Description

Program: ph.x / PHonon / Quantum ESPRESSO (version: 7.2)

TABLE OF CONTENTS

INTRODUCTION

Line-of-input: title_line

&INPUTPH

amass | outdir | prefix | niter_ph | tr2_ph | alpha_mix(niter) | nmix_ph | verbosity | reduce_io | max_seconds | dftd3_hess | fildyn | fildrho | fildvscf | epsil | lrpa | lnoloc | trans | lraman | eth_rps | eth_ns | dek | recover | low_directory_check | only_init | qplot | q2d | q_in_band_form | electron_phonon | el_ph_nsigma | el_ph_sigma | ahc_dir | ahc_nbnd | ahc_nbndskip | skip_upperfan | lshift_q | zeu | zue | elop | fpol | ldisp | nogg | asr | ldiag | lqdir | search_sym | nq1 | nq2 | nq3 | nk1 | nk2 | nk3 | k1 | k2 | k3 | diagonalization | read_dns_bare | ldvscf_interpolate | wpot_dir | do_long_range | do_charge_neutral | start_irr | last_irr | nat_todo | modenum | start_q | last_q | dvscf_star | drho_star

Line-of-input: xq(1) xq(2) xq(3)

qPointsSpecs

nqs | xq1 | xq2 | xq3 | nq

Line-of-input: atom(1) atom(2) ... atom(nat_todo)

ADDITIONAL INFORMATION

INTRODUCTION

Input data format: { } = optional, [ ] = it depends, # = comment

Structure of the input data:
===============================================================================

title_line

&INPUTPH
   ...
/

[ xq(1) xq(2) xq(3) ]                        # if ldisp != .true.  and  qplot != .true.

[ nqs                                        # if qplot == .true. 
  xq(1,i)    xq(2,i)    xq(3,1)    nq(1)
  ...
  xq(1,nqs)  xq(2,nqs)  xq(3,nqs)  nq(nqs) ]

[ atom(1)  atom(2)  ... atom(nat_todo) ]     # if nat_todo was specified
   

Line of input

Syntax:

title_line  

Description of items:

title_line CHARACTER
Title of the job, i.e., a line that is reprinted on output.
         

Namelist: &INPUTPH

amass(i), i=1,ntyp REAL
Default: 0.0
Atomic mass [amu] of each atomic type.
If not specified, masses are read from data file.
         
outdir CHARACTER
Default: value of the ESPRESSO_TMPDIR environment variable if set;
current directory ('./') otherwise
Directory containing input, output, and scratch files;
must be the same as specified in the calculation of
the unperturbed system.
         
prefix CHARACTER
Default: 'pwscf'
Prepended to input/output filenames; must be the same
used in the calculation of unperturbed system.
         
niter_ph INTEGER
Default: maxter=100
Maximum number of iterations in a scf step. If you want
more than 100, edit variable "maxter" in PH/phcom.f90
         
tr2_ph REAL
Default: 1e-12
 Threshold for self-consistency.
         
alpha_mix(niter) REAL
Default: alpha_mix(1)=0.7
Mixing factor (for each iteration) for updating
the scf potential:

vnew(in) = alpha_mix*vold(out) + (1-alpha_mix)*vold(in)
         
nmix_ph INTEGER
Default: 4
Number of iterations used in potential mixing. Using a larger value (8~20)
can significantly speed up convergence, at the cost of using more memory.
         
verbosity CHARACTER
Default: 'default'
 Options are:
            
'debug', 'high', 'medium' :
 verbose output
            
'low', 'default', 'minimal' :
 short output
            
reduce_io LOGICAL
Default: .false.
Reduce I/O to the strict minimum.

BEWARE: If the input flag reduce_io=.true. was
used, it is not allowed to restart from an interrupted
run.
         
max_seconds REAL
Default: 1.d7
 Maximum allowed run time before the job stops smoothly.
         
dftd3_hess CHARACTER
Default: 'prefix.hess'
 File where the D3 dispersion hessian matrix is read.
         
fildyn CHARACTER
Default: 'matdyn'
 File where the dynamical matrix is written.
         
fildrho CHARACTER
Default: ' '
File where the charge density responses are written. Note that the file
will actually be saved as ${outdir}/_ph0/${prefix}.${fildrho}1
where  ${outdir}, ${prefix} and ${fildrho} are the values of the
corresponding input variables
         
fildvscf CHARACTER
Default: ' '
File where the the potential variation is written
(for later use in electron-phonon calculation, see also fildrho).
         
epsil LOGICAL
Default: .false.
If .true. in a q=0 calculation for a non metal the
macroscopic dielectric constant of the system is
computed. Do not set epsil to .true. if you have a
metallic system or q/=0: the code will complain and stop.

Note: the input value of epsil will be ignored if ldisp=.true.
(the code will automatically set epsil to .false. for metals,
to .true. for insulators: see routine PHonon/PH/prepare_q.f90).
         
lrpa LOGICAL
Default: .false.
If .true. the dielectric constant is calculated at the
RPA level with DV_xc=0.
         
lnoloc LOGICAL
Default: .false.
If .true. the dielectric constant is calculated without
local fields, i.e. by setting DV_H=0 and DV_xc=0.
         
trans LOGICAL
Default: .true.
If .false. the phonons are not computed.
If trans .and. epsil are both .true.,
the effective charges are calculated.
If ldisp is .true., trans=.false. is overridden
(except for the case of electron-phonon calculations)
         
lraman LOGICAL
Default: .false.
If .true. calculate non-resonant Raman coefficients
using second-order response as in:
M. Lazzeri and F. Mauri, PRL 90, 036401 (2003).
         

Optional variables for Raman:

eth_rps REAL
Default: 1.0d-9
 Threshold for calculation of  Pc R |psi>.
            
eth_ns REAL
Default: 1.0e-12
 Threshold for non-scf wavefunction calculation.
            
dek REAL
Default: 1.0e-3
 Delta_xk used for wavefunction derivation wrt k.
            
recover LOGICAL
Default: .false.
 If .true. restart from an interrupted run.
         
low_directory_check LOGICAL
Default: .false.
If .true. search in the phsave directory only the
                 quantities requested in input.
         
only_init LOGICAL
Default: .false.
If .true. only the bands and other initialization quantities are calculated.
(used for GRID parallelization)
         
qplot LOGICAL
Default: .false.
 If .true. a list of q points is read from input.
         
q2d LOGICAL
Default: .false.
If .true. three q points and relative weights are
read from input. The three q points define the rectangle
q(:,1) + l (q(:,2)-q(:,1)) + m (q(:,3)-q(:,1)) where
0< l,m < 1. The weights are integer and those of points two
and three are the number of points in the two directions.
         
q_in_band_form LOGICAL
Default: .false.
This flag is used only when qplot is .true. and q2d is
.false.. When .true. each couple of q points q(:,i+1) and
q(:,i) define the line from q(:,i) to q(:,i+1) and nq
points are generated along that line. nq is the weigth of
q(:,i). When .false. only the list of q points given as
input is calculated. The weights are not used.
         
electron_phonon CHARACTER
Default: ' '
Options are:
            
'simple' :
Electron-phonon lambda coefficients are computed
for a given q and a grid of k-points specified by
the variables nk1, nk2, nk3, k1, k2, k3.
            
'interpolated' :
Electron-phonon is calculated by interpolation
over the Brillouin Zone as in M. Wierzbowska, et
al. arXiv:cond-mat/0504077
            
'lambda_tetra' :
The electron-phonon coefficient \lambda_{q \nu}
is calculated with the optimized tetrahedron method.
            
'gamma_tetra' :
The phonon linewidth \gamma_{q \nu} is calculated
from the electron-phonon interactions
using the optimized tetrahedron method.
            
'epa' :
Electron-phonon coupling matrix elements are written
to file prefix.epa.k for further processing by program
epa.x which implements electron-phonon averaged (EPA)
approximation as described in G. Samsonidze & B. Kozinsky,
Adv. Energy Mater. 2018, 1800246 doi:10.1002/aenm.201800246
arXiv:1511.08115
            
'ahc' :
Quantities required for the calculation of phonon-induced
electron self-energy are computed and written to the directory
ahc_dir. The output files can be read by postahc.x for
the calculation of electron self-energy.
Available for both metals and insulators.
trans=.false. is required.
            
For metals only, requires gaussian smearing (except for 'ahc').

If trans=.true., the lambdas are calculated in the same
run, using the same k-point grid for phonons and lambdas.
If trans=.false., the lambdas are calculated using
previously saved DeltaVscf in fildvscf, previously saved
dynamical matrix, and the present punch file. This allows
the use of a different (larger) k-point grid.
            
el_ph_nsigma INTEGER
Default: 10
The number of double-delta smearing values used in an
electron-phonon coupling calculation.
         
el_ph_sigma REAL
Default: 0.02
The spacing between double-delta smearing values used in
an electron-phonon coupling calculation.
         

Variables for electron_phonon = 'ahc':

ahc_dir CHARACTER
Default: outdir // 'ahc_dir/'
Directory where the output binary files are written.
            
ahc_nbnd INTEGER
Status: REQUIRED
Number of bands for which the electron self-energy is to be computed.
            
ahc_nbndskip INTEGER
Default: 0
Number of bands to exclude when computing the self-energy. Self-energy
is computed for bands with indices from ahc_nbndskip+1 to
ahc_nbndskip+ahc_nbnd. ahc_nbndskip+ahc_nbnd cannot
exceed nbnd of the preceding SCF or NSCF calculation.
            
skip_upperfan LOGICAL
Default: .false.
If .true., skip calculation of the upper Fan self-energy, which
involves solving the Sternheimer equation.
            
lshift_q LOGICAL
Default: .false.
Use a wave-vector grid displaced by half a grid step
in each direction - meaningful only when ldisp is .true.
When this option is set, the q2r.x code cannot be used.
         
zeu LOGICAL
Default: zeu=epsil
If .true. in a q=0 calculation for a non metal the
effective charges are computed from the dielectric
response. This is the default algorithm. If epsil=.true.
and zeu=.false. only the dielectric tensor is calculated.
         
zue LOGICAL
Default: .false.
If .true. in a q=0 calculation for a non metal the
effective charges are computed from the phonon
density responses. This is an alternative algorithm,
different from the default one (if trans .and. epsil )
The results should be the same within numerical noise.
         
elop LOGICAL
Default: .false.
If .true. calculate electro-optic tensor.
         
fpol LOGICAL
Default: .false.
If .true. calculate dynamic polarizabilities
Requires epsil=.true. ( experimental stage:
see example09 for calculation of methane ).
         
ldisp LOGICAL
Default: .false.
If .true. the run calculates phonons for a grid of
q-points specified by nq1, nq2, nq3 - for direct
calculation of the entire phonon dispersion.
         
nogg LOGICAL
Default: .false.
If .true. disable the "gamma_gamma" trick used to speed
up calculations at q=0 (phonon wavevector) if the sum over
the Brillouin Zone includes k=0 only. The gamma_gamma
trick exploits symmetry and acoustic sum rule to reduce
the number of linear response calculations to the strict
minimum, as it is done in code phcg.x.
         
asr LOGICAL
Default: .false.
Apply Acoustic Sum Rule to dynamical matrix, effective charges
Works only in conjunction with "gamma_gamma" tricks (see above)
         
ldiag LOGICAL
Default: .false.
If .true. forces the diagonalization of the dynamical
matrix also when only a part of the dynamical matrix
has been calculated. It is used together with start_irr
and last_irr. If all modes corresponding to a
given irreducible representation have been calculated,
the phonon frequencies of that representation are
correct. The others are zero or wrong. Use with care.
         
lqdir LOGICAL
Default: .false.
If .true. ph.x creates inside outdir a separate subdirectory
for each q vector. The flag is set to .true. when ldisp=.true.
and fildvscf /= ' ' or when an electron-phonon
calculation is performed. The induced potential is saved
separately for each q inside the subdirectories.
         
search_sym LOGICAL
Default: .true.
Set it to .false. if you want to disable the mode
symmetry analysis.
         
nq1, nq2, nq3 INTEGER
Default: 0,0,0
Parameters of the Monkhorst-Pack grid (no offset) used
when ldisp=.true. Same meaning as for nk1, nk2, nk3
in the input of pw.x.
         
nk1, nk2, nk3, k1, k2, k3 INTEGER
Default: 0,0,0,0,0,0
When these parameters are specified the phonon program
runs a pw non-self consistent calculation with a different
k-point grid thant that used for the charge density.
This occurs even in the Gamma case.
nk1, nk2, nk3 are the parameters of the Monkhorst-Pack grid
with offset determined by k1, k2, k3.
         
diagonalization CHARACTER
Default: 'david'
Diagonalization method for the non-SCF calculations.
            
'david' :
Davidson iterative diagonalization with overlap matrix
(default). Fast, may in some rare cases fail.
            
'cg' :
Conjugate-gradient-like band-by-band diagonalization.
Slower than 'david' but uses less memory and is
(a little bit) more robust.
            
read_dns_bare LOGICAL
Default: .false.
If .true. the PH code tries to read three files in the DFPT+U
calculation: dns_orth, dns_bare, d2ns_bare.
dns_orth and dns_bare are the first-order variations of
the occupation matrix, while d2ns_bare is the second-order
variation of the occupation matrix. These matrices are
computed only once during the DFPT+U calculation. However,
their calculation (especially of d2ns_bare) is computationally
expensive, this is why they are written to file and then can be
read (e.g. for restart) in order to save time.
         
ldvscf_interpolate LOGICAL
Default: .false.
If .true., use Fourier interpolation of phonon potential
to compute the induced part of phonon potential at each
q point. Results of a dvscf_q2r.x run is needed.
Requires trans = .false..
         

Optional variables for dvscf interpolation:

wpot_dir CHARACTER
Default: outdir // 'w_pot/'
Directory where the w_pot binary files are written.
Must be the same with wpot_dir used in dvscf_q2r.x.
The real space potential files are stored in wpot_dir
with names ${prefix}.wpot.irc${irc}//"1".
            
do_long_range LOGICAL
Default: .false.
If .true., add the long-range part of the potential
to the Fourier interpolated potential as in:
S. Ponce et al, J. Chem. Phys. 143, 102813 (2015).
Reads dielectric matrix and Born effective charges from
the ${wpot_dir}/tensors.dat file, written in dvscf_q2r.x.
Currently, only the dipole (Frohlich) part is implemented.
The quadrupole part is not implemented.
            
do_charge_neutral LOGICAL
Default: .false.
If .true., impose charge neutrality on the Born effective
charges. Used only if do_long_range = .true..
            

Specification of irreducible representation

start_irr INTEGER
Default: 1
See: last_irr
Perform calculations only from start_irr to last_irr
irreducible representations.

IMPORTANT:
   * start_irr must be <= 3*nat
   * do not specify nat_todo together with
     start_irr, last_irr
            
last_irr INTEGER
Default: 3*nat
See: start_irr
Perform calculations only from start_irr to last_irr
irreducible representations.

IMPORTANT:
   * start_irr must be <= 3*nat
   * do not specify nat_todo together with
     start_irr, last_irr
            
nat_todo INTEGER
Default: 0, i.e. displace all atoms
Choose the subset of atoms to be used in the linear response
calculation: nat_todo atoms, specified in input (see below)
are displaced. Can be used to estimate modes for a molecule
adsorbed over a surface without performing a full fledged
calculation. Use with care, at your own risk, and be aware
that this is an approximation and may not work.
IMPORTANT:
   * nat_todo <= nat
   * if linear-response is calculated for a given atom, it
     should also be done for all symmetry-equivalent atoms,
     or else you will get incorrect results
            
modenum INTEGER
Default: 0
For single-mode phonon calculation : modenum is the index of the
irreducible representation (irrep) into which the reducible
representation formed by the 3*nat atomic displacements are
decomposed in order to perform the phonon calculation.
Note that a single-mode calculation will not give you the
frequency of a single phonon mode: in general, the selected
"modenum" is not an eigenvector. What you get on output is
a column of the dynamical matrix.
            

q-point specification

start_q INTEGER
Default: 1
See: last_q
Used only when ldisp=.true..
Computes only the q points from start_q to last_q.

IMPORTANT:
   * start_q must be <= nqs (number of q points found)
   * do not specify nat_todo together with
     start_q, last_q
            
last_q INTEGER
Default: number of q points
See: start_q
Used only when ldisp=.true..
Computes only the q points from start_q to last_q.

IMPORTANT
   * last_q must be <= nqs (number of q points)
   * do not specify nat_todo together with
     start_q, last_q
            
dvscf_star STRUCTURE
Default: disabled
It contains the following components:

dvscf_star%open  (logical, default: .false.)
dvscf_star%dir   (character, default: outdir//"Rotated_DVSCF" or the
                  ESPRESSO_FILDVSCF_DIR environment variable)
dvscf_star%ext   (character, default: "dvscf") the extension to use
                  for the name of the output files, see below
dvscf_star%basis (character, default: "cartesian") the basis on which
                  the rotated dvscf will be saved
dvscf_star%pat   (logical, default: false) save an optional file with the
                  displacement patterns and q vector for each dvscf file

IF dvscf_star%open is .true. use symmetry to compute and store the variation
of the self-consistent potential on every q* in the star of the present q.

The rotated dvscf will then be stored in directory dvscf_star%dir with name
prefix.dvscf_star%ext.q_name//"1". Where q_name is derived from the coordinates
of the q-point, expressed as fractions in crystalline coordinates
(notice that ph.x reads q-points in cartesian coordinates).
E.g. q_cryst= (0, 0.5, -0.25) -> q_name = "0_1o2_-1o4"

The dvscf can be represented on a basis of cartesian 1-atom displacements
(dvscf_star%basis='cartesian') or on the basis of the modes at the rotated q-point
(dvscf_star%basis='modes'). Notice that the el-ph wannier code requires 'cartesian'.
Each dvscf file comes with a corresponding pattern file with an additional ".pat"
suffix; this file contains information about the basis and the q-point of the dvscf.

Note: rotating dvscf can require a large amount of RAM memory and can be i/o
      intensive; in its current implementation all the operations are done
      on a single processor.
Note2: this feature is currently untested with image parallelisation.
            
drho_star STRUCTURE
Default: disabled
See: dvscf_star
It contains the following components:

drho_star%open  (logical, default: .false.)
drho_star%dir   (character, default: outdir//"Rotated_DRHO" or the
                 ESPRESSO_FILDRHO_DIR environment variable)
drho_star%ext   (character, default: "drho") the extension to use
                 for the name of the output files, see below
drho_star%basis (character, default: "modes") the basis on which
                 the rotated drho will be saved
drho_star%pat   (logical, default: true) save an optional file with the
                 displacement patterns and q vector for each drho file

Like dvscf_star, but for the perturbation of the charge density.
Notice that the defaults are different.
            
IF ldisp != .true. and qplot != .true. :

Line of input

Syntax:

xq(1) xq(2) xq(3)   

Description of items:

xq(1) xq(2) xq(3) REAL
The phonon wavevector, in units of 2pi/a0
(a0 = lattice parameter).
Not used if ldisp=.true. or qplot=.true.
               
ELSEIF qplot == .true. :

Specification of q points when qplot == .true.

Card: qPointsSpecs

Syntax:

nqs  
 xq1(1)   xq2(1)   xq3(1)   nq(1) 
 xq1(2)   xq2(2)   xq3(2)   nq(2) 
 . . .
 xq1(nqs)   xq2(nqs)   xq3(nqs)   nq(nqs) 

Description of items:

nqs INTEGER
Number of q points in the list. Used only if qplot=.true.
                     
xq1, xq2, xq3 REAL
q-point coordinates; used only with ldisp=.true. and qplot=.true.
The phonon wavevector, in units of 2pi/a0 (a0 = lattice parameter).
The meaning of these q points and their weights nq depend on the
flags q2d and q_in_band_form. (NB: nq is integer)
                        
nq INTEGER
The weight of the q-point; the meaning of nq depends
on the flags q2d and q_in_band_form.
                        
IF nat_todo was specified :

Line of input

Syntax:

atom(1) atom(2) ... atom(nat_todo)   

Description of items:

atom(1) atom(2) ... atom(nat_todo) INTEGER
Contains the list of indices of atoms used in the
calculation if nat_todo is specified.
               

ADDITIONAL INFORMATION

NB: The program ph.x writes on the tmp_dir/_ph0/{prefix}.phsave directory
a file for each representation of each q point. This file is called
dynmat.#iq.#irr.xml where #iq is the number of the q point and #irr
is the number of the representation. These files contain the
contribution to the dynamical matrix of the irr representation for the
iq point.

If recover=.true. ph.x does not recalculate the
representations already saved in the tmp_dir/_ph0/{prefix}.phsave
directory.  Moreover ph.x writes on the files patterns.#iq.xml in the
tmp_dir/_ph0/{prefix}.phsave directory the displacement patterns that it
is using. If recover=.true. ph.x does not recalculate the
displacement patterns found in the tmp_dir/_ph0/{prefix}.phsave directory.

This mechanism allows:

  1) To recover part of the ph.x calculation even if the recover file
     or files are corrupted. You just remove the _ph0/{prefix}.recover
     files from the tmp_dir directory. You can also remove all the _ph0
     files and keep only the _ph0/{prefix}.phsave directory.

  2) To split a phonon calculation into several jobs for different
     machines (or set of nodes). Each machine calculates a subset of
     the representations and saves its dynmat.#iq.#irr.xml files on
     its tmp_dir/_ph0/{prefix}.phsave directory. Then you collect all the
     dynmat.#iq.#irr.xml files in one directory and run ph.x to
     collect all the dynamical matrices and diagonalize them.

NB: To split the q points in different machines, use the input
variables start_q and last_q. To split the irreducible
representations, use the input variables start_irr, last_irr. Please
note that different machines will use, in general, different
displacement patterns and it is not possible to recollect partial
dynamical matrices generated with different displacement patterns.  A
calculation split into different machines will run as follows: A
preparatory run of ph.x with start_irr=0, last_irr=0 produces the sets
of displacement patterns and save them on the patterns.#iq.xml files.
These files are copied in all the tmp_dir/_ph0/{prefix}.phsave directories
of the machines where you plan to run ph.x. ph.x is run in different
machines with complementary sets of start_q, last_q, start_irr and
last_irr variables.  All the files dynmat.#iq.#irr.xml are
collected on a single tmp_dir/_ph0/{prefix}.phsave directory (remember to
collect also dynmat.#iq.0.xml).  A final run of ph.x in this
machine collects all the data contained in the files and diagonalizes
the dynamical matrices.  This is done requesting a complete dispersion
calculation without using start_q, last_q, start_irr, or last_irr.
See an example in examples/GRID_example.

On parallel machines the q point and the irreps calculations can be split
automatically using the -nimage flag. See the phonon user guide for further
information.
      
This file has been created by helpdoc utility on Thu Mar 30 12:13:56 CEST 2023.