quantum-espresso/atomic/Doc/INPUT_LD1.def

1011 lines
27 KiB
Modula-2

input_description -distribution {Quantum Espresso} -program ld1.x {
toc {}
intro {
Input data cards for ld1.x program:
Always present:
1 namelist &input
1.1 optional cards for all-electron calculations
Needed for PP generation:
2 namelist &inputp
2.1 additional cards for PP generation
Needed for pseudo-potential (PP) test. optional for PP generation:
3 namelist &test
3.1 optional cards for PP test
}
namelist INPUT {
label { This namelist is always needed ! }
var title -type CHARACTER {
info { A string describing the job. }
status { OPTIONAL }
}
group {
label { Either: }
var zed -type REAL {
see { atom }
info {
The nuclear charge (1 < zed < 100).
@b IMPORTANT:
Specify either @ref zed OR @ref atom, not both!
}
}
label { Or: }
var atom -type CHARACTER {
see { zed }
info {
Atomic symbol: atom='H', 'He', 'Be', etc.
@b IMPORTANT:
Specify either @ref atom OR @ref zed, not both!
}
}
}
group {
label { Radial grid parameters: }
var xmin -type REAL {
see { dx }
default {
-7.0 if @ref iswitch>1 or @ref rel=0,
-8.0 otherwise
}
info { Radial grid parameter. }
}
var dx -type REAL {
info {
Radial grid parameter.
The radial grid is: r(i+1) = exp(xmin+i*dx)/zed a.u.
}
default {
0.0125 if @ref iswitch>1,
0.008 otherwise
}
}
var rmax -type REAL {
info { Outermost grid point. }
default { 100.0 a.u. }
}
}
var beta -type REAL {
info { parameter for potential mixing }
default { 0.2 }
}
var tr2 -type REAL {
info { convergence threshold for scf }
default { 1e-14 }
}
var iswitch -type INTEGER {
info {
1 all-electron calculation
2 PP test calculation
3 PP generation
4 LDA-1/2 correction, needs a previously generated PP file
}
default { 1 }
}
group {
label { Parameters for logarithmic derivatives: }
var nld -type INTEGER {
info { the number of logarithmic derivatives to be calculated }
}
var rlderiv -type REAL {
info { radius (a.u.) at which logarithmic derivatives are calculated }
}
vargroup -type REAL {
var eminld
var emaxld
info {
Energy range (min, max energy, in Ry) at which
logarithmic derivatives are calculated.
}
}
var deld -type REAL {
info {
Delta e (Ry) of energy for logarithmic derivatives.
}
}
var rpwe -type REAL {
info { radius (a.u.) at which partial wave expansions are calculated }
default { rlderiv }
}
message {
If the above parameters are not specified, logarithmic
derivatives and partial wave expansions are not calculated.
}
}
var rel -type INTEGER {
info {
0 ... non relativistic calculation
1 ... scalar relativistic calculation
2 ... full relativistic calculation with spin-orbit
}
default {
0 for Z <= 18;
1 for Z > 18
}
}
var lsmall -type LOGICAL {
default { .false. }
info {
if .true. writes on files the small component
}
}
var max_out_wfc -type INTEGER {
default { 7 }
info {
Maximum number of atomic wavefunctions written in the output
file.
}
}
var noscf -type LOGICAL {
default { .false. }
info {
if .true. the charge is not computed. The occupations are
not used and the eigenvalues and eigenfunctions are those
of a hydrogen-like atom.
}
}
var lsd -type INTEGER {
info {
0 ... non spin polarized calculation
1 ... spin-polarized calculation
@b BEWARE:
not allowed if iswitch=3 (PP generation) or with full
relativistic calculation
}
default { 0 }
}
var dft -type CHARACTER {
info {
Exchange-correlation functional.
Examples:
@b 'PZ' Perdew and Zunger formula for LDA
@b 'PW91' Perdew and Wang GGA
@b 'BP' Becke and Perdew GGA
@b 'PBE' Perdew, Becke and Ernzerhof GGA
@b 'BLYP' ...
For the complete list, see module "functionals" in ../Modules/
The default is @b 'PZ' for all-electron calculations,
it is read from the PP file in a PP calculation.
}
}
var latt -type INTEGER {
info {
0 ... no Latter correction
1 ... apply Latter correction
}
default { 0 }
}
var isic -type INTEGER {
info {
0 ... no Self-interaction correction
1 ... apply Self-interaction correction
}
default { 0 }
status {
only for all-electron calculation
}
}
var rytoev_fact -type REAL {
default { as specified in file Modules/constants.f90 }
info {
Factor used to convert Ry into eV.
}
}
var cau_fact -type REAL {
default { as specified in file Modules/constants.f90 }
info {
Speed of light in a.u..
(Be careful the default value is always used in the
relativistic exchange.)
}
}
var vdw -type LOGICAL {
default { .false. }
info {
If .true., the frequency dependent polarizability and van der
Waals coefficient C6 will be computed in Thomas-Fermi and
von Weizsaecker approximation(only for closed-shell ions).
}
status {
Gradient-corrected DFT not yet implemented.
}
}
var prefix -type CHARACTER {
default { 'ld1' }
info {
Prefix for file names - only for output file names
containing the orbitals, logarithmic derivatives, tests
See below for file names and the content of the file.
}
}
var verbosity -type CHARACTER {
default { 'low' }
info {
@b 'low' or @b 'high'
if @b 'high' with iswitch=2,3 prints separately core and
valence contributions to the energies. Print the
frozen-core energy.
}
}
var file_charge -type CHARACTER {
default { ' ' }
info {
Name of the file where the code writes the all-electron
total charge. No charge is written if file_charge=' '.
}
}
var config -type CHARACTER {
default { ' ' }
info {
A string with the electronic configuration.
Example:
'[Ar] 3d10 4s2 4p2.5'
* If @ref lsd=1, spin-up and spin-down state may appear twice
with the respective occupancy: 3p4 3p2 = 4 up,
2 down. Otherwise, the Hund's rule is assumed.
* If @ref rel=2, states with jj=l-1/2 are filled first.
If a state appears twice, the first one has jj=l-1/2,
the second one jj=l+1/2 (except S states)
(Use rel_dist if you want to average the electrons
over all available states.)
* If config='default' the code uses @ref zed to set the ground
state electronic configuration for the atom.
Negative occupancies are used to flag unbound states;
they are not actually used.
}
}
var relpert -type LOGICAL {
default { .false. }
info {
If .true. the relativistic corrections to the non-relativistic
Kohn-Sham energy levels (@ref rel=0 .and. @ref lsd=0) are computed using
first-order perturbation theory in all-electron calculations.
The corrections consist of the following terms:
E_vel: velocity (p^4) correction
E_Dar: Darwin term
E_S-O: spin-orbit coupling
The spin-orbit term vanishes for s-electron states and gives
rise to a splitting of (2*l+1)*E_S-O for the other states.
The separate contributions are printed only if verbosity='high'.
Formulas and notation are based on the Herman-Skillman book:
F. Herman and S. Skillman, "Atomic Structure Calculations",
Prentice-Hall, Inc., Englewood Cliffs, New Jersey, 1963
}
}
var rel_dist -type CHARACTER {
default { 'energy' }
info {
@b 'energy' or @b 'average'
* if @b 'energy' the relativistic l-1/2 states are filled first.
* if @b 'average' the electrons are uniformly distributed
among all the states with the given l.
}
}
var write_coulomb -type LOGICAL {
default { .false. }
info {
If .true., a fake pseudo-potential file with name X.UPF,
where X is the atomic symbol, is written. It contains
the radial grid and the wavefunctions as specified in input,
plus the info needed to build the Coulomb potential
for an all-electron calculation - for testing only.
}
}
}
card AllElectronCards -nameless 1 {
label {
If config is empty the electronic configuration is read from
the following cards:
}
choose {
when -test "rel < 2" {
syntax {
line {
var nwf -type INTEGER {
info {
number of wavefunctions
}
}
}
table AE_wfs {
rows -start 1 -end nwf {
col nl -type CHARACTER {
info { wavefunction label (e.g. 1s, 2s, etc.) }
}
col n -type INTEGER { info { principal quantum number } }
col l -type INTEGER { info { angular quantum number } }
col oc -type REAL { info { occupation number } }
col isw -type INTEGER {
info { the spin index (1-2) used only in the lsda case }
}
}
}
}
}
elsewhen -test "rel = 2" {
syntax {
line {
var nwf
}
table AE_wfs {
rows -start 1 -end nwf {
col nl
col n
col l
col oc
col jj -type REAL {
info {
The total angular momentum (0.0 is allowed for complete
shells: the codes fills 2l states with jj=l-1/2,
2l+2 with jj=l+1/2).
}
}
}
}
}
}
}
}
namelist INPUTP {
var zval -type REAL {
default { (calculated) }
info {
Valence charge.
zval is automatically calculated from available data.
If the value of zval is provided in input, it will be
checked versus the calculated value. The only case in
which you need to explicitly provide the value of zval
for noninteger zval (i.e. half core-hole pseudo-potentials).
}
}
var pseudotype -type INTEGER {
info {
1 ... norm-conserving, single-projector PP
IMPORTANT: if pseudotype=1 all calculations are done
using the SEMILOCAL form, not the separable nonlocal form
2 ... norm-conserving PP in separable form (obsolescent)
All calculations are done using SEPARABLE non-local form
IMPORTANT: multiple projectors allowed but not properly
implemented, use only if you know what you are doing
3 ... ultrasoft PP or PAW
}
}
var file_pseudopw -type CHARACTER {
status { REQUIRED }
info {
File where the generated PP is written.
* if the file name ends with "upf" or "UPF",
or in any case for spin-orbit PP (rel=2),
the file is written in UPF format;
* if the file name ends with 'psp' it is
written in native CPMD format (this is currently
an experimental feature); otherwise it is written
in the old "NC" format if pseudotype=1, or
in the old RRKJ format if pseudotype=2 or 3
(no default, must be specified).
}
}
var file_recon -type CHARACTER {
info {
File containing data needed for GIPAW reconstruction
of all-electron wavefunctions from PP results.
If you want to use additional states to perform the
reconstruction, add them at the end of the list
of all-electron states.
}
default { ' ' }
}
var lloc -type INTEGER {
default { -1 }
info {
Angular momentum of the local channel.
* lloc=-1 or lloc=-2 pseudizes the all-electron potential
if lloc=-2 the original recipe of Troullier-Martins
is used (zero first and second derivatives at r=0)
* lloc>-1 uses the corresponding channel as local PP
NB: if lloc>-1, the corresponding channel must be the last in the
list of wavefunctions appearing after the namelist &inputp
In the relativistic case, if lloc > 0 both the j=lloc-1/2 and
the j=lloc+1/2 wavefunctions must be at the end of the list.
}
}
var rcloc -type REAL {
status {
Must be specified only if @ref lloc=-1, otherwise the
corresponding value of @ref rcut is used.
}
info {
Matching radius (a.u.) for local pseudo-potential (no default).
}
}
var nlcc -type LOGICAL {
default { .false. }
info {
If .true. produce a PP with the nonlinear core
correction of Louie, Froyen, and Cohen
[PRB 26, 1738 (1982)].
}
}
var new_core_ps -type LOGICAL {
default { .false. }
status { requires nlcc=.true. }
info {
If .true. pseudizes the core charge with bessel functions.
}
}
var rcore -type REAL {
info {
Matching radius (a.u.) for the smoothing of the core charge.
If not specified, the matching radius is determined
by the condition: rho_core(rcore) = 2*rho_valence(rcore)
}
}
var tm -type LOGICAL {
default { .false. }
info {
* .true. for Troullier-Martins pseudization [PRB 43, 1993 (1991)]
* .false. for Rappe-Rabe-Kaxiras-Joannopoulos pseudization
[PRB 41, 1227 (1990), erratum PRB 44, 13175 (1991)]
}
}
var rho0 -type REAL {
info {
Charge at the origin: when the Rappe-Rabe-Kaxiras-Joannopoulos
method with 3 Bessel functions fails, specifying rho0 > 0
may allow to override the problem (using 4 Bessel functions).
Typical values are in the order of 0.01-0.02
}
default { 0.0 }
}
var lpaw -type LOGICAL {
info {
If .true. produce a PAW dataset, experimental feature
only for @ref pseudotype=3
}
default { .false. }
}
group {
var which_augfun -type CHARACTER {
default {
@tt 'AE' for Vanderbilt-Ultrasoft pseudo-potentials and 'BESSEL' for PAW datasets.
}
info {
If different from @b 'AE' the augmentation functions are pseudized
before @ref rmatch_augfun. The pseudization options are:
* @b 'PSQ' Use Bessel functions to pseudize Q
from the origin to rmatch_augfun.
These features are available only for PAW:
* @b 'BESSEL' Use Bessel functions to pseudize the Q.
* @b 'GAUSS' Use 2 Gaussian functions to pseudize the Q.
* @b 'BG' Use original Bloechl's recipe with a single gaussian.
Note: if lpaw is true and which_augfun is set to AE real all-
electron charge will be used, which will produce extremely
hard augmentation.
}
}
var rmatch_augfun -type REAL {
default { 0.5 a.u. }
status { Used only if which_augfun is different from 'AE'. }
info {
Pseudization radius for the augmentation functions. Presently
it has the same value for all L.
}
}
var rmatch_augfun_nc -type REAL {
default { .false. }
status { Used only if which_augfun is 'PSQ'. }
info {
If .true. the augmentation functions are pseudized
from the origin to min(rcut(ns),rcut(ns1)) where ns
and ns1 are the two channels for that Q. In this case
rmatch_augfun is not used.
}
}
}
var lsave_wfc -type LOGICAL {
default { .false. if .not. lpaw, otherwise .true. }
info {
Set it to .true. to save all-electron and pseudo wavefunctions
used in the pseudo-potential generation in the UPF file. Only
works for UPFv2 format.
}
}
var lgipaw_reconstruction -type LOGICAL {
default { .false. }
info {
Set it to .true. to generate pseudo-potentials containing the
additional info required for reconstruction of all-electron
orbitals, used by GIPAW. You will typically need to specify
additional projectors beyond those used in the generation of
pseudo-potentials. You should also specify @ref file_recon.
All projectors used in the reconstruction must be listed BOTH
in the test configuration after namelist &test AND in the
all-electron configuration (variable 'config', namelist &inputp,
Use negative occupancies for projectors on unbound states). The
core radii in the test configuration should be the same as in
the pseudo-potential generation section and will be used as the
radius of reconstruction. Projectors not used to generate the
pseudo-potential should have zero occupation number.
}
}
var use_paw_as_gipaw -type LOGICAL {
default { .false. }
info {
When generating a PAW dataset, setting this option to .true. will
save the core all-electron wavefunctions to the UPF file.
The GIPAW reconstruction to be performed using the PAW data and
projectors for the valence wavefunctions.
In the default case, the GIPAW valence wavefunction and projectors
are independent from the PAW ones and must be then specified as
explained above in lgipaw_reconstruction.
Setting this to .true. always implies @ref lgipaw_reconstruction = .true.
}
}
var author -type CHARACTER {
info { Name of the author. }
default { 'anonymous' }
}
var file_chi -type CHARACTER {
info { file containing output PP chi functions }
default { ' ' }
}
var file_beta -type CHARACTER {
info { file containing output PP beta functions }
default { ' ' }
}
var file_qvan -type CHARACTER {
info { file containing output PP qvan functions }
default { ' ' }
}
var file_screen -type CHARACTER {
info { file containing output screening potential }
default { ' ' }
}
var file_core -type CHARACTER {
info { file containing output total and core charge }
default { ' ' }
}
var file_wfcaegen -type CHARACTER {
info { file with the all-electron wfc for generation }
default { ' ' }
}
var file_wfcncgen -type CHARACTER {
info { file with the norm-conserving wfc for generation }
default { ' ' }
}
var file_wfcusgen -type CHARACTER {
info { file with the ultra-soft wfc for generation }
default { ' ' }
}
}
card PseudoPotentialGenerationCards -nameless 1 {
choose {
when -test "rel=0 OR rel=2" {
syntax {
line {
var nwfs -type INTEGER {
info { number of wavefunctions to be pseudized }
}
}
table PP_wfs {
rows -start 1 -end nwfs {
col nls -type CHARACTER {
info {
Wavefunction label (same as in the all-electron configuration).
}
}
col nns -type INTEGER {
info {
Principal quantum number (referred to the PSEUDOPOTENTIAL case;
nns=1 for lowest s, nns=2 for lowest p, and so on).
}
}
col lls -type INTEGER { info { Angular momentum quantum number. } }
col ocs -type REAL { info { Occupation number (same as in the all-electron configuration). } }
col ener -type REAL {
info {
Energy (Ry) used to pseudize the corresponding state.
If 0.d0, use the one-electron energy of the all-electron state.
Do not use 0.d0 for unbound states!
}
}
col rcut -type REAL { info { Matching radius (a.u.) for norm conserving PP. } }
col rcutus -type REAL {
info { Matching radius (a.u.) for ultrasoft PP - only for pseudotype=3. }
}
col jjs -type REAL { info { The total angular momentum (0.0 is allowed for complete shells). } }
}
}
}
message {
* if @ref lloc>-1 the state with @ref lls=@ref lloc must be the last
* if @ref lloc>0 in the relativistic case, both states with @ref jjs=@ref lloc-1/2
and @ref jjs=@ref lloc+1/2 must be the last two
}
}
otherwise {
syntax {
line {
var nwfs
}
table PP_wfs {
rows -start 1 -end nwfs {
col nls
col nns
col lls
col ocs
col ener
col rcut
col rcutus
}
}
}
}
}
}
namelist TEST {
label { needed only if iswitch=2 or iswitch=4, optional if iswitch=3 }
var nconf -type INTEGER {
info { the number of configurations to be tested. For iswitch = 4 nconf=2 }
default { 1 }
}
var file_pseudo -type CHARACTER {
status { ignored if iswitch=3 }
info {
File containing the PP.
* If the file name contains ".upf" or ".UPF",
the file is assumed to be in UPF format;
* else if the file name contains ".rrkj3" or ".RRKJ3",
the old RRKJ format is first tried;
* otherwise, the old NC format is read.
IMPORTANT: in the latter case, all calculations are done
using the SEMILOCAL form, not the separable nonlocal form.
Use the UPF format if you want to test the separable form!
}
default { ' ' }
}
vargroup -type REAL {
var ecutmin
var ecutmax
var decut
info {
Parameters (Ry) used for test with a basis set of spherical
Bessel functions j_l(qr) . The hamiltonian at fixed scf
potential is diagonalized for various values of ecut:
@ref ecutmin, @ref ecutmin+@ref decut, @ref ecutmin+2*@ref decut ... up to @ref ecutmax.
This yields an indication of convergence with the
corresponding plane-wave cutoff in solids, and shows
in an unambiguous way if there are "ghost" states
}
default {
decut=5.0 Ry; ecutmin=ecutmax=0Ry
}
status { specify @ref ecutmin and @ref ecutmax if you want to perform this test }
}
var rm -type REAL {
info { Radius of the box used with spherical Bessel functions. }
default { 30 a.u. }
}
dimension configts -start 1 -end nconf -type CHARACTER {
info {
A string containing the test valence electronic
configuration nc, nc=1,nconf. Same syntax as for "config".
If configts(i) is not set, the electron configuration
is read from the cards following the namelist.
}
}
dimension lsdts -start 1 -end nconf -type INTEGER {
default { 1 }
see { lsd }
info {
0 or 1. It is the value of lsd used in the i-th test.
Allows to make simultaneously spin-polarized and
spin-unpolarized tests.
}
}
var frozen_core -type LOGICAL {
default { .false. }
info {
If .true. only the core wavefunctions of the first
configuration are calculated. The eigenvalues, orbitals
and energies of the other configurations are calculated
with the core of the first configuration.
The first configuration must be spin-unpolarized.
}
}
var rcutv -type REAL {
info { Cutoff distance (CUT) for the inclusion of LDA-1/2 potential.
Needed (mandatory) only if iswitch = 4 }
default { -1.0 }
}
}
card PseudoPotentialTestCards -nameless 1 {
label { IMPORTANT: this card has to be specified for each missing configts(i) }
choose {
when -test "lsd=1" {
syntax {
line {
var nwfts -type INTEGER { info { number of wavefunctions } }
}
table test_wfs {
rows -start 1 -end nwfts {
col elts -type CHARACTER { see { nls } }
col nnts -type INTEGER { see { nns } }
col llts -type INTEGER { see { lls } }
col octs -type REAL { see { ocs } }
col enerts -type REAL { status { not used } }
col rcutts -type REAL { status { not used } }
col rcutusts -type REAL { status { not used } }
col iswts -type INTEGER { info { spin index (1 or 2, used in lsda case) } }
}
}
}
}
elsewhen -test "rel=2" {
syntax {
line { var nwfts }
table test_wfs {
rows -start 1 -end nwfts {
col elts
col nnts
col llts
col octs
col enerts
col rcutts
col rcutusts
col jjts -type REAL { info { total angular momentum of the state } }
}
}
}
}
otherwise {
syntax {
line { var nwfts }
table test_wfs {
rows -start 1 -end nwfts {
col elts
col nnts
col llts
col octs
col enerts
col rcutts
col rcutusts
}
}
}
}
}
}
section -title {Notes} {
text {
For PP generation you do not need to specify namelist &test, UNLESS:
1. you want to use a different configuration for unscreening wrt the
one used to generate the PP. This is useful for PP with semicore
states: use semicore states ONLY to produce the PP, use semicore
AND valence states (if occupied) to make the unscreening
2. you want to specify some more states for PAW style reconstruction of
all-electron orbitals from pseudo-orbitals
}
subsection -title {Output files written} {
text {
* file_tests "prefix".test results of transferability test
for each testing configuration N:
* file_wavefunctions "prefix"N.wfc all-electron KS orbitals
* file_wavefunctionsps "prefix"Nps.wfc pseudo KS orbitals
if lsd=1:
* file_wavefunctions "prefix"N.wfc.up all-electron KS up orbitals
* file_wavefunctions "prefix"N.wfc.dw all-electron KS down orbitals
if rel=2 and lsmall=.true.:
* file_wavefunctions "prefix".wfc.small all-electron KS small component
if parameters for logarithmic derivatives are specified:
* file_logder "prefix"Nps.dlog all-electron logarithmic derivatives
* file_logderps "prefix"Nps.dlog pseudo logarithmic derivatives
"N" is not present if there is just one testing configuration.
}
}
subsection -title {Recipes to reproduce old all-electron atomic results with the ld1 program} {
text {
* The Hartree results in Phys. Rev. 59, 299 (1940) or in
Phys. Rev. 59, 306 (1940) can be reproduced with:
rel=0,
isic=1,
dft='NOX-NOC'
* The Herman-Skillman tables can be reproduced with:
rel=0,
isic=0,
latt=1,
dft='SL1-NOC'
* Data on the paper Liberman, Waber, Cromer Phys. Rev. 137, A27 (1965) can be
reproduced with:
rel=2,
isic=0,
latt=1,
dft='SL1-NOC'
* Data on the paper S. Cohen Phys. Rev. 118, 489 (1960) can be reproduced with:
rel=2,
isic=1,
latt=0,
dft='NOX-NOC'
* The revised PBE described in PRL 80, 890 (1998) can be obtained with:
isic=0
latt=0
dft='SLA-PW-RPB-PBC' or 'dft='revPBE'
* The relativistic energies of closed shell atoms reported in PRB 64 235126 (2001)
can be reproduced with:
isic=0
latt=0
cau_fact=137.0359895
dft='sla-vwn' for the LDA case
dft='PBE' for the PBE case
* The NIST results in PRA 55, 191 (1997):
LDA:
rel=0
dft='sla-vwn'
LSD:
rel=0
lsd=1
dft='sla-vwn'
RLDA
rel=2
rel_dist='average'
dft='rxc-vwn'
ScRLDA:
rel=1
dft='rxc-vwn'
}
}
}
}