The Fleur input file

Fleur expects all input (execpt some computational settings specified on the command line) in the inp.xml file in the current working directory. This file is usually produced by the input-generator 'inpgen' but Fleur workflows require adjustments to this file.

Example inp.xml file

<?xml version="1.0" encoding="UTF-8" standalone="no"?>
<fleurInput fleurInputVersion="0.31">
   <comment>
      Si bulk                                                                         
   </comment>
   <calculationSetup>
      <cutoffs Kmax="3.70" Gmax="11.10" GmaxXC="9.20" numbands="0"/>
      <scfLoop itmax="15" minDistance=".00001" maxIterBroyd="99"
               imix="Anderson" alpha=".05" precondParam="0.0" spinf="2.0"/>
      <coreElectrons ctail="T" frcor="F" kcrel="0" coretail_lmax="0"/>
      <magnetism jspins="1" l_noco="F" swsp="F" lflip="F"/>
      <soc theta=".00000000" phi=".00000000" l_soc="F" spav="F"/>
      <expertModes gw="0" secvar="F"/>
      <geometryOptimization l_f="F" forcealpha="1.00" forcemix="BFGS"
                            epsdisp=".00001" epsforce=".00001"/>
      <ldaU l_linMix="F" mixParam=".050000" spinf="1.000000"/>
      <bzIntegration valenceElectrons="20.00000000" mode="hist" 
                     fermiSmearingEnergy=".00100000">
         <kPointCount count="17" gamma="F"/>
         <altKPointSet purpose="bands">
            <kPointCount count="   240" gamma="F"/>
         </altKPointSet>
      </bzIntegration>
      <energyParameterLimits ellow="-.80000000" elup="1.00000000"/>
   </calculationSetup>
   <cell>
      <symmetryFile filename="sym.out"/>
      <bulkLattice scale="1.0000000000" latnam="any">
         <bravaisMatrix>
            <row-1>.0000000000 5.1306085335 5.1306085335</row-1>
            <row-2>5.1306085335 .0000000000 5.1306085335</row-2>
            <row-3>5.1306085335 5.1306085335 .0000000000</row-3>
         </bravaisMatrix>
      </bulkLattice>
   </cell>
   <xcFunctional name="pbe" relativisticCorrections="F"/>
   <atomSpecies>
      <species name="Si-1" element="Si" atomicNumber="14" coreStates="2"
               magMom=".00000000" flipSpin="T">
         <mtSphere radius="2.17000000" gridPoints="717"
                   logIncrement=".01600000"/>
         <atomicCutoffs lmax="8" lnonsphr="6"/>
         <energyParameters s="3" p="3" d="3" f="4"/>
         <lo type="SCLO" l="1" n="2" eDeriv="0"/>
      </species>
   </atomSpecies>
   <atomGroups>
      <atomGroup species="Si-1">
         <relPos label="   1">1.00/8.00 1.00/8.00 1.00/8.00</relPos>
         <relPos label="   2">-1.00/8.00 -1.00/8.00 -1.00/8.00</relPos>
         <force calculate="T" relaxXYZ="TTT"/>
      </atomGroup>
   </atomGroups>
   <output dos="F" band="F" vacdos="F" slice="F" mcd="F">
      <checks vchk="F" cdinf="F"/>
      <densityOfStates ndir="0" minEnergy="-.50" maxEnergy=".50" 
                       sigma=".0150"/>
      <vacuumDOS layers="0" integ="F" star="F" nstars="0" locx1=".00000"
                 locy1=".00000" locx2=".00000" locy2=".00000" nstm="0"
                 tworkf=".00000"/>
      <unfoldingBand unfoldBand="F" supercellX="1" 
                     supercellY="1" supercellZ="1"/>
      <plotting iplot="0" score="F" plplot="F"/>
      <chargeDensitySlicing numkpt="0" minEigenval=".00000" 
                            maxEigenval=".00000" nnne="0" pallst="F"/>
      <specialOutput eonly="F" bmt="F"/>
      <magneticCircularDichroism energyLo="-10.0000" energyUp=".0000"/>
   </output>
   <!-- We include the file relax.inp here to enable relaxations 
        (see documentation) -->
   <xi:include xmlns:xi="http://www.w3.org/2001/XInclude" 
               href="relax.xml">
      <xi:fallback/>
   </xi:include>
</fleurInput>

Being an XML file inp.xml starts with some general XML information in line 1. The rest of the file is enclosed in the <fleurInput> element that carries as an attribute the version number of the input file format. Within <fleurInput> the input file consists of several sections to be discussed in detail below.

Comment and Constants sections

<comment>
   Si bulk
</comment>

The comment section is optional. It encloses a simple line of text that is written out as part of the inp.xml into the out.xml file.

<constants>
   <constant name="myConst" value="5.1673552752"/>
</constants>

The constants section is optional and not part of the example inp.xml file shown above. It can be used to define constants that can then be used in other parts of the XML input file, e.g., the lattice setup or the declaration of the atom positions. The constants element may enclose multiple constant definitions. Each one has to provide the name and value of the respective constant.

Calculation Setup section

<calculationSetup>
   <cutoffs Kmax="3.70" Gmax="11.10" GmaxXC="9.20" numbands="0"/>
   <scfLoop itmax="15" minDistance=".00001" maxIterBroyd="99"
            imix="Anderson" alpha=".05" precondParam="0.0" spinf="2.0"/>
   <coreElectrons ctail="T" frcor="F" kcrel="0" coretail_lmax="0"/>
   <magnetism jspins="1" l_noco="F" swsp="F" lflip="F"/>
   <soc theta=".00000000" phi=".00000000" l_soc="F" spav="F"/>
   <expertModes gw="0" secvar="F"/>
   <geometryOptimization l_f="F" forcealpha="1.00" forcemix="BFGS"
                         epsdisp=".00001" epsforce=".00001"/>
   <ldaU l_linMix="F" mixParam=".050000" spinf="1.000000"/>
   <bzIntegration valenceElectrons="20.00000000" mode="hist" 
                  fermiSmearingEnergy=".00100000">
      <kPointCount count="17" gamma="F"/>
      <altKPointSet purpose="bands">
         <kPointCount count="   240" gamma="F"/>
      </altKPointSet>
   </bzIntegration>
   <energyParameterLimits ellow="-.80000000" elup="1.00000000"/>
</calculationSetup>

The calculation setup section covers the input of general numerical parameters controlling the Fleur calculation.

Tag Attribute Description
cutoffs The main cutoff parameters for the calculation
Kmax The cutoff for the basis functions
Gmax The cutoff for the density
GmaxXC The cutoff used for the potential when the XC functional is calculated
numbands The number of eigenvalues to be calculated at each k point. A value of 0 is converted to a default value.
scfLoop Parameters controlling the number of SCF loop iterations and the mixing scheme
itmax The number of SCF loop iterations to be performed by Fleur
maxIterBroyd The number of iterations to be taken into account by Broyden-like mixing schemes
imix The mixing scheme. This can be one of "straight", "Broyden1", "Broyden2", and "Anderson"
alpha The mixing factor
precondParam The preconditioning parameter for bulk metals. Typical value: 0.8. Choose higher mixing parameter, e.g. 0.25.
spinf The spin mixing factor enhancement
coreElectrons Parameters for the treatment of the core electrons
ctail Use core tail corrections.
frcor The frozen core approximation can be activated here.
kcrel If true fully relativistic core routines are used, otherwise only spin-polarized routines.
coretail_lmax Cutoff for the expansion of the core-tail into other MT spheres. Also relevant for initial charge generation.
magnetism Parameters for controlling the degree of magnetism considered in the calculation
jspins The number of spins to be considered: 1 for nonmagnetic calculations and 2 for calculations incorporating magnetism.
l_noco Set this to true to consider not only collinear but also non-collinear magnetism
l_J Set this to true to calculate Heisenberg parameters.
swsp If true generate spin-polarized from unpolarized density.
lflip If true flip spin directions for each atom with set flipSpin flag.
soc Parameters needed for calculations with spin-orbit coupling
theta The spin quantization axis is given by the theta and phi angles.
phi The spin quantization axis is given by the theta and phi angles.
l_soc This switch is used to toggle the consideration of spin-orbit coupling in the respective calculation .
spav Construct spin-orbit operator from spin-averaged potential.
nocoParams See the section on the non-collinear magnetism setup for details.
expertModes Parameters for the control of certain advanced Fleur calculation modes
gw The different output modes for GW approximation calculations are set here.
secvar Treat the nonspherical part of the Hamiltonian in second variation.
geometryOptimization Parameters required for force calculations and the optimization of atom positions
l_f l_f is used to switch on the calculation of forces.
forcealpha The mixing parameter for the forces when performing structural optimizations.
forcemix The mixing procedure for the forces. Choose one of "BFGS" and "straight"
epsdisp The convergence criterion for atom displacements.
epsforce The convergence criterion for the forces.
For further options on relaxation..
ldaU Optional parameters for the density matrix mixing in LDA+U calculations. See the LDA+U setup for details.
bzIntegration Parameters required for the Brillouin zone integration
valenceElectrons The number of electrons to be represented within the valence electron framework.
mode The Brillouin zone integration mode. It can be one of
hist - Use the histogram mode. This is the default.
gauss - Use Gaussian smearing.
tria - Use the tetrahedron method.
fermiSmearingEnergy The Fermi smearing can be parametrized by this energy in Hartree.
fermiSmearingTemp As an alternative to fermiSmearingEnergy a Fermi smearing temperature can be set in Kelvin.
kPointCount See the section on the k point set setup for details.
kPointMesh See the section on the k point set setup for details.
kPointList See the section on the k point set setup for details.
altKPointSet See the section on the k point set setup for details.
energyParameterLimits Boundaries for energy parameters
ellow
elup

Cell section

<cell>
   <symmetryFile filename="sym.out"/>
   <bulkLattice scale="1.0000000000" latnam="any">
      <bravaisMatrix>
         <row-1>.0000000000 5.1306085335 5.1306085335</row-1>
         <row-2>5.1306085335 .0000000000 5.1306085335</row-2>
         <row-3>5.1306085335 5.1306085335 .0000000000</row-3>
      </bravaisMatrix>
   </bulkLattice>
</cell>

The cell section covers the declaration of the symmetry operations available in the unit cell and the definition of the lattice vectors. For details please see the sections on the Bravais lattice setup and the unit cell symmetry setup.

Exchange correlation functional section

<xcFunctional name="pbe" relativisticCorrections="F"/>

The exchange correlation functional section consists of a single xml element with the two attributes name and relativisticCorrections. The XC functional is specified by the name attribute:

LDA functionals
x-a
wign
mjw
pz The functional by Perdew and Zunger
vwn The functional by Vosko, Wilk, and Nusair
bh The functional by Barth and Hedin
GGA functionals
pw91 The functional by Perdew and Wang
pbe The functional by Perdew, Burke, and Ernzerhof
rpbe The revPBE functional by Zhang and Yang
Rpbe The RPBE functional by Hammer, Hansen, and Nørskov
wc The functional by Wu and Cohen
Tag Attribute Description
relativisticCorrections A boolean switch used to activate optional relativistic corrections according to MacDonnald-Vosko.

Atom species section

   <atomSpecies>
      <species name="Si-1" element="Si" atomicNumber="14" coreStates="2"
               magMom=".00000000" flipSpin="T">
         <mtSphere radius="2.17000000" gridPoints="717"
                   logIncrement=".01600000"/>
         <atomicCutoffs lmax="8" lnonsphr="6"/>
         <energyParameters s="3" p="3" d="3" f="4"/>
         <lo type="SCLO" l="1" n="2" eDeriv="0"/>
      </species>
   </atomSpecies>

The atomSpecies section is a tool to set identical numerical parameters for the atoms of different atom groups without introducing redundancy. For this several species with a unique name can be defined in the section. In the following atomGroups section each atom group is associated to one of the species.

Tag Attribute Description
species The element defining a species. There can be multiple of these elements in this section.
name A name for the species. This should be unique within the set of species.
element The abbreviation of the chemical element.
atomicNumber The atomic number of the chemical element.
coreStates The number of states whose electrons are to be treated as core electrons.
magMom The magnetic moment of each atom in the unit cell belonging to this species.
flipSpin A boolean switch to flip the spin direction for the associated atoms.

A species element contains other elements to determine its numerical parameters. Please note that the order of these elements in the input file is predefined:

Tag Attribute Description
mtSphere This element is used to define the properties of the muffin-tin spheres.
radius The radius of the MT sphere.
gridPoints The number of grid points on the radial mesh for this MT sphere.
logIncrement The logarithmic increment of the radial mesh.
atomicCutoffs This element covers the l-cutoffs.
lmax The general l-cutoff for all atoms of the species.
lnonsphr The reduced l-cutoff for the calculation of contributions originating from non-spherical part of the potential
lmaxAPW If present the APW+lo approach will be used. This is the cutoff defining the lmax up to which LAPWs are used if no APW+lo local orbital is defined for the respective l. See also the respective article by G.K.H. Madsen.
energyParameters This element sets the energy parameters.
s The main quantum number for the valence s electrons.
p The main quantum number for the valence p electrons.
d The main quantum number for the valence d electrons.
f The main quantum number for the valence f electrons.
special
socscale A float in range from 0.0 to 1.0. Scales the magnitude of SOC at the species. socscale="0.0" switches SOC off.
vca_charge A float to specify an extra charge for calculations in the virtual crystal approximation for this species.
lda Logical switch to use LDA for this atom.
b_field_mt Zeeman field applied to this atom.
electronConfig See the description of the electron configuration setup for details.
ldaU Up to one ldaU element can be present to define a U parameter for this atom and a certain l channel.
l This is the l channel the U is supposed to affect.
U This is the U parameter in eV.
J This is the J parameter in eV.
l_amf If true the around mean field limit is employed, otherwise the atomic limit.
For a more detailed description have a look at [[xmlLDAUSetup
lo This element is used to introduce local orbitals to each of the associated atoms. It can be present multiple times. See the section on the local orbital setup for a more profound discussion.
type The type of the LO. This can be SCLO for semicore LOs or HELO for LOs at higher energies
l The angular momentum quantum number belonging to this LO
n The main quantum number for this LO
eDeriv The energy derivative of the additional radial function introduced with this LO. This is by default 0 to obtain conventional LOs. For HDLOs it has to be set to a finite positive integer value.

Atom groups section

<atomGroups>
   <atomGroup species="Si-1">
      <relPos label="   1">1.00/8.00 1.00/8.00 1.00/8.00</relPos>
      <relPos label="   2">-1.00/8.00 -1.00/8.00 -1.00/8.00</relPos>
      <force calculate="T" relaxXYZ="TTT"/>
   </atomGroup>
</atomGroups>

The atom groups section covers parameters relevant for each group of symmetry equivalent atoms.

Tag Attribute Description
atomGroup There is at least one atom group. Each atom in the unit cell has to be in one.
species The name of the species of this group's atoms.

Each atom group element encloses certain other elements:

Tag Attribute Description
relPos See below
filmPos See below
force Some switches associated to force calculations.
calculate This boolean switch determines whether forces are calculated for the atoms of this group.
relaxXYZ Three boolean switches used declare in which directions the atom position may be optimized in force relaxation calculations.
nocoParams See the description of the non-collinear magnetism setup for details.

The atom positions are defined within each atomGroup of symmetry equivalent atoms in the atom groups section of the input file. They can be provided as relative or film positions:

Relative positions

<atomGroup species="W-1">
   <relPos label="   1">.000000 1.0/2.0 .060000</relPos>
   <relPos label="   2">1.0/2.0 .000000 -.060000</relPos>
   <force calculate="T" relaxXYZ="TTT"/>
</atomGroup>

Typically for bulk materials the atom positions are provided in relative coordinates as fractions of the three lattice vectors. For this the relPos tag is used. In the example the atom group consists of two atoms at two different positions. The first one is the representative atom. As shown the relative coordinates are provided as three numbers within the relPos element. Note that each coordinate can also be provided by a short mathematical expression that does not contain any spaces, e.g., 1.0/4.0.

Each atom can be associated to a label that can be specified in the inpgen input. If there is no explicit labeling in the inpgen input the labels just associate a unique number with each atom. This labeling is only meant as a help to the user to keep an overview for setups of very large unit cells.

Film positions

<atomGroup species="W-2">
   <filmPos label="   1">1.0/2.0 1.0/2.0 -12.0314467594</filmPos>
   <filmPos label="   2">1.0/2.0 1.0/2.0 12.0314467594</filmPos>
   <force calculate="T" relaxXYZ="TTT"/>
</atomGroup>

For calculations on films the atom positions are provided within the filmPos element. Here, the first two of the coordinates are relative, while the third coordinate in the direction normal to the film plane is absolute and in atomic units (Bohr radii).

Output section

<output dos="F" band="F" vacdos="F" slice="F" mcd="F">
   <checks vchk="F" cdinf="F"/>
   <densityOfStates ndir="0" minEnergy="-.50" maxEnergy=".50" 
                    sigma=".0150"/>
   <vacuumDOS layers="0" integ="F" star="F" nstars="0" locx1=".00000"
              locy1=".00000" locx2=".00000" locy2=".00000" nstm="0"
              tworkf=".00000"/>
   <unfoldingBand unfoldBand="F" supercellX="1" 
                  supercellY="1" supercellZ="1"/>
   <plotting iplot="0" score="F" plplot="F"/>
   <chargeDensitySlicing numkpt="0" minEigenval=".00000" 
                         maxEigenval=".00000" nnne="0" pallst="F"/>
   <specialOutput eonly="F" bmt="F"/>
   <magneticCircularDichroism energyLo="-10.0000" energyUp=".0000"/>
</output>

The output section is optional. It covers parameters relevant for the generation of special output.

Tag Attribute Description
output
dos A boolean switch that determines whether a density of states has to be calculated.
band A boolean switch that determines whether a band structure calculation should be performed.
vacdos A boolean switch that determines whether a vacuum DOS has to be calculated.
slice A boolean switch controlling whether a slice has to be calculated. The associated parameters are set in the chargeDensitySlicing element.
coreSpec A boolean switch controlling whether a core spectrum has to be calculated. The associated parameters are set in the coreSpectrum element.
mcd A boolean switch controlling whether a magnetic circular dichroism spectrum has to be calculated. The associated parameters are set in the magneticCircularDichroism element.

If an attribute of the output element is set to true the associated enclosed element has to be present:

Tag Attribute Description
checks The checks element covers several switches to perform and write out certain tests.
vchk Continuity checks of the potential at the MT and vacuum boundaries.
cdinf Calculation of partial charges and the continuity of the density.
disp Calculation of the distance between the in- and output potential.
densityOfStates Parameters for DOS calculations are set in this element.
ndir Select the DOS calculation mode. For details have a look at the respective section.
minEnergy The lower energy boundary of the window for the DOS plot in Htr.
maxEnergy The upper energy boundary of the window for the DOS plot in Htr.
sigma The Gaussian smearing factor used in the plot whenever the tetrahedron method is not used.
vacuumDOS
layers The number of layers in which the vacuum DOS is integrated. The value can be between 1 and 99.
integ If integ is true the vacuum DOS is also integrated in the direction normal to the film.
star if star is true, star coefficients are calculated at values of izlay for 0 (=q) to nstars-1.
nstars The number of star functions to be used (0th star is given by value of q=charge integrated in 2D)
locx1, locy1, locx2, locy2 The four attributes loc[xy]1, loc[xy]2 are used to calculate a local DOS at a certain vertical z position (or integrated in z). The local DOS is restricted to an area in the 2D unit cell which is defined by the corner points given by loc[xy]1 and loc[xy]2. The two corners have to be provided in internal coordinates.
nstm For the consideration of STM: 0: s-Tip, 1: p_z-Tip, 2: d_z^2-Tip (following Chen's derivative rule)
tworkf For the consideration of STM: Workfunction of Tip (in Hartree units) is needed for d_z^2-Orbital.
unfoldingBand
unfoldBand A boolean switch that defines if unfolding is used and additional weights are written.
supercellX The size of the supercell (in units of simple unit cells) (iteger value) in X direction.
supercellY The size of the supercell (in units of simple unit cells) (iteger value) in Y direction.
supercellZ The size of the supercell (in units of simple unit cells) (iteger value) in Z direction.
plotting The plotting element groups several switches to plot the density and the potential.
iplot Calculate a charge density plot. For details see the respective section.
score If true the core charge is excluded from the plot.
plplot This switch allows the plotting of the potential from its respective files.
chargeDensitySlicing
numkpt This is the number of the k-point which is used for the slice (k=0 : all k-points are used)
minEigenval This is the lower boundary for eigenvalues in the slice.
maxEigenval This is the upper boundary for eigenvalues in the slice.
nnne The number of eigenvalues used for the slice (nnne=0 : all eigenvalues between boundaries are taken into account)
pallst Set this to true if states above the Fermi level are plotted.
specialOutput
form66 Use this switch to write out a formatted eigenvector file.
eonly This switch can be used to prevent prevent the output of eigenvectors into associated file.
bmt This switch is used to generate a charge density file 'cdnbmt' with suppressed magnetization in the interstitial and vacuum regions.
coreSpectrum See the section on core spectrum calculations for details.
magneticCircularDichroism See the section on magnetic circular dichroism calculations for details.

Force therorem section

<forceTheorem>
   <MAE theta="0.0 0.1*Pi" phi="0.0 0.2*Pi" />
</forceTheorem>

It is possible to apply the magnetic force theorem to perform efficient approximative calculations on the magnetocrystalline anisotropy energy (MAE), spin-spiral dispersions, the Dzyaloshinskii Moriya interaction (DMI), and the Heisenberg exchange interaction. To do so a forceTheorem section has to be inserted directly below the output section. For details see the section on the respective calculations.

Inclusion of the relax.xml file

<!-- We include the file relax.inp here to enable relaxations 
   (see documentation) -->
<xi:include xmlns:xi="http://www.w3.org/2001/XInclude" 
            href="relax.xml">
   <xi:fallback/>
</xi:include>

It is possible to automatically include other xml files into the Fleur input file by using the respective XInclude feature. By default this is done for the relax.xml file. The syntax above also provides a fallback that is used when the relax.xml file is not present: In that case it is just not included, but a message about this is written out to the terminal.

When starting new Fleur calculations the relax.xml file is typically not present. It is automatically created when geometry optimizations for the atom positions are performed. Its contents are the atom displacements calculated so far. These automatically modify the atom positions provided in the Fleur input file. Furthermore a history of atom positions and related forces for the previous optimization steps is found in the relax.xml file. This is needed for the determination of the next displacements for ongoing atom position optimization calculations.

Details on how geometry optimizations are performed are found in the respective section.