2 ! Copyright (C) 2002-2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
6 ! main routine for the EXCITING code
14 ! perform the appropriate task
20 write(*,'("EXCITING version ",I1.1,".",I2.2,".",I3.3)') version
43 case(72,73,82,83,142,143,152,153)
73 write(*,'("Error(main): task not defined : ",I8)') task
82 ! !TITLE: The EXCITING Code Manual\\ Version 0.9.150
83 ! !AUTHORS: J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl
85 ! !INTRODUCTION: Introduction
86 ! Welcome to the {\sf EXCITING} Code Manual!
87 ! The {\sf EXCITING} code is a state-of-the-art full-potential linearised
88 ! augmented-plane-wave (FP-LAPW) code for determining the properties of
89 ! crystalline solids. It was developed mainly at the
90 ! Karl-Franzens-Universit\"{a}t Graz as part of the {\sf EXCITING} EU Research
91 ! and Training Network project \cite{exciting}. The guiding philosophy during
92 ! the implementation of the code was to keep it as simple as possible for both
93 ! users and developers without compromising on its capabilities. All the
94 ! routines are released under either the GNU General Public License (GPL) or
95 ! the GNU Lesser General Public License (LGPL) in the hope that they may
96 ! inspire other scientists to implement new and, well, exciting developments
97 ! in the field of density functional theory and beyond.
99 ! \section{Acknowledgements}
100 ! Lots of people contributed to the {\sf EXCITING} code with ideas, checking
101 ! and testing, writing code or documentation and general encouragement. They
102 ! include Lars Nordstr\"{o}m, Clas Persson, Christian Brouder, Rickard
103 ! Armiento, Andrew Chizmeshya, Per Anderson, Igor Nekrasov, Fredrik Bultmark,
104 ! Sushil Auluck, Frank Wagner, Fateh Kalarasse, J\"{u}rgen Spitaler, Stefano
105 ! Pittalis, Nektarios Lathiotakis, Tobias Burnus, Stephan Sagmeister,
106 ! Christian Meisenbichler, Francesco Cricchio, S\'{e}bastien Leb\`{e}gue,
107 ! Yigang Zhang, Fritz K\"{o}rmann and Alexey Baranov. Special mention of David
108 ! Singh's very useful book {\it Planewaves, Pseudopotentials and the LAPW
109 ! Method} \cite{singh} must also be made. Finally we would like to acknowledge
110 ! the generous support of Karl-Franzens-Universit\"{a}t Graz, as well as the
111 ! EU Marie-Curie Research Training Networks initiative.
114 ! Kay Dewhurst, Sangeeta Sharma and Claudia Ambrosch-Draxl
117 ! Edinburgh, Berlin and Leoben, February 2008
121 ! Unless explicitly stated otherwise, {\sf EXCITING} uses atomic units. In
122 ! this system $\hbar=1$, the electron mass $m=1$, the Bohr radius $a_0=1$ and
123 ! the electron charge $e=1$ (note that the electron charge is positive, so
124 ! that the atomic numbers $Z$ are negative). Thus, the atomic unit of length
125 ! is 0.52917720859(36) \AA, and the atomic unit of energy is the Hartree which
126 ! equals 27.21138386(68) eV. The unit of the external magnetic fields is
127 ! defined such that one unit of magnetic field in {\tt exciting.in} equals
128 ! 1717.2445320376 Tesla.
130 ! \section{Compiling and running {\sf EXCITING}}
131 ! \subsection{Compiling the code}
132 ! Unpack the code from the archive file. Run the command
136 ! in the {\tt exciting} directory and select the appropriate system and
137 ! compiler. We highly recommend that you edit the file {\tt make.inc} and tune
138 ! the compiler options for your particular system. You can also make use of
139 ! machine-optimised {\tt BLAS/LAPACK} libraries if they are available, but
140 ! make sure they are version $3.x$. Setting the {\sf OpenMP} options of your
141 ! compiler will enable {\sf EXCITING} to run in parallel mode on
142 ! multiprocessor systems. Following this, run
146 ! This will hopefully compile the entire code and all the libraries into one
147 ! executable, {\tt exciting}, located in the {\tt src} directory. It will also
148 ! compile a few useful auxilliary programs, namely {\tt spacegroup} for
149 ! producing crystal geometries from spacegroup data, {\tt species} for
150 ! generating species files, and {\tt eos} for fitting equations of state to
151 ! energy-volume data. If you want to compile everything all over again, then
152 ! run {\tt make clean} from the {\tt exciting} directory, followed by
155 ! \subsection{Running the code}
156 ! As a rule, all input files for the code are in lower case and end with the
157 ! extension {\tt .in}. All output files are uppercase and have the extension
158 ! {\tt .OUT}. For most cases, the user will only need to modify the file
159 ! {\tt exciting.in}. In this file input parameters are arranged in blocks.
160 ! Each block consists of a block name on one line and the block variables on
161 ! subsequent lines. Almost all blocks are optional: the code uses reasonable
162 ! default values in cases where they are absent. Blocks can appear in any
163 ! order, if a block is repeated then the second instance is used. Comment
164 ! lines can be included in the input file and begin with the {\tt !}
167 ! The only other input files are those describing the atomic species which go
168 ! into the crystal. These files are found in the {\tt species} directory and
169 ! are named with the element symbol and the extension {\tt .in}, for example
170 ! {\tt Sb.in}. They contain parameters like the atomic charge, mass,
171 ! muffin-tin radius, occupied atomic states and the type of linearisation
172 ! required. Users should not have to modify these files in the majority of
175 ! The best way to learn to use {\sf EXCITING} is to run the examples included
176 ! with the package. These can be found in the {\tt examples} directory and use
177 ! many of the code's capabilities. The following section which describes all
178 ! the input parameters will be of invaluable assistance.
180 ! \section{Input blocks}
181 ! This section lists all the input blocks available. It is arranged with the
182 ! name of the block followed by a table which lists each parameter name, what
183 ! the parameter does, its type and default value. A horizontal line in the
184 ! table indicates a new line in {\tt exciting.in}. Below the table is a brief
185 ! overview of the block's function.
187 ! \subsection{{\tt atoms}}
188 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
190 ! {\tt nspecies} & number of species & integer & 0 \\
192 ! {\tt spfname(i)} & species filename for species {\tt i} & string & - \\
194 ! {\tt natoms(i)} & number of atoms for species {\tt i} & integer & - \\
196 ! {\tt atposl(j,i)} & atomic position in lattice coordinates for atom {\tt j}
198 ! {\tt bfcmt(j,i)} & muffin-tin external magnetic field in Cartesian
199 ! coordinates for atom {\tt j} & real(3) & - \\
201 ! \end{tabularx}\newline\newline
202 ! Defines the atomic species as well as their positions in the unit cell and
203 ! the external magnetic field applied throughout the muffin-tin. These fields
204 ! are used to break spin symmetry and should be considered infinitesimal as
205 ! they do not contribute directly to the total energy. Collinear calculations
206 ! are more efficient if the field is applied in the $z$-direction. One could,
207 ! for example, set up an anti-ferromagnetic crystal by pointing the field on
208 ! one atom in the positive $z$-direction and in the opposite direction on
209 ! another atom. If {\tt molecule} is {\tt .true.} then the atomic positions
210 ! are assumed to be in Cartesian coordinates. See also {\tt sppath},
211 ! {\tt bfieldc} and {\tt molecule}.
213 ! \subsection{{\tt autokpt}}
214 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
216 ! {\tt autokpt} & {\tt .true.} if the {\bf k}-point set is to be determined
217 ! automatically & logical & {\tt .false.} \\
219 ! \end{tabularx}\newline\newline
220 ! See {\tt rlambda} for details.
222 ! \subsection{{\tt autormt}}
223 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
225 ! {\tt autormt} & {\tt .true.} if muffin-tin radii should be determined
226 ! automatically & logical & {\tt .false.} \\
228 ! \end{tabularx}\newline\newline
229 ! See {\tt rmtapm} for details.
231 ! \subsection{{\tt avec}}
232 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
234 ! {\tt avec(1) } & first lattice vector & real(3) & $(1.0,0.0,0.0)$ \\
236 ! {\tt avec(2) } & second lattice vector & real(3) & $(0.0,1.0,0.0)$ \\
238 ! {\tt avec(3) } & third lattice vector & real(3) & $(0.0,0.0,1.0)$ \\
240 ! \end{tabularx}\newline\newline
241 ! Lattice vectors of the crystal in atomic units (Bohr). If {\tt molecule} is
242 ! {\tt .true.} then these vectors are not used.
244 ! \subsection{{\tt beta0}}
245 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
247 ! {\tt beta0 } & initial mixing parameter and increment & real & $0.1$ \\
249 ! \end{tabularx}\newline\newline
250 ! This sets the initial parameter used for mixing the old and new potentials
251 ! during the self-consistent cycle. For some materials, such as magnetic
252 ! metals, this should be made smaller to avoid instability. The code
253 ! automatically adjusts the mixing parameter to the optimial size. Making
254 ! {\tt beta0} too large can result in instability and poor convergence. See
255 ! {\tt betamax} as well as the routine {\tt mixer}.
257 ! \subsection{{\tt betamax}}
258 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
260 ! {\tt betamax } & maximum mixing parameter & real & $1.0$ \\
262 ! \end{tabularx}\newline\newline
263 ! The mixing parameter is adjusted in increments of {\tt beta0} to optimise
264 ! that rate of convergece. {\tt betamax} sets the upper limit to this
265 ! parameter. Making this too large can result in poor convergence due to
268 ! \subsection{{\tt bfieldc}}
269 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
271 ! {\tt bfieldc} & global external magnetic field in Cartesian coordinates &
272 ! real(3) & $(0.0,0.0,0.0)$ \\
274 ! \end{tabularx}\newline\newline
275 ! This is a constant magnetic field applied thoughout the entire unit cell
276 ! and enters the second-variational Hamiltonian as
277 ! $$ \frac{g_e\alpha}{4}\,\vec{\sigma}\cdot{\bf B}_{\rm ext}, $$
278 ! where $g_e$ is the electron $g$-factor (2.0023193043718). This field is
279 ! normally used to break spin symmetry for spin-polarised calculations and
280 ! considered to be infinitesimal with no direct contribution to the total
281 ! energy. In cases where the magnetic field is finite (for example when
282 ! computing magnetic response) the external ${\bf B}$-field energy reported in
283 ! {\tt INFO.OUT} should be added to the total by hand. This field is applied
284 ! throughout the entire unit cell. To apply magnetic fields in particular
285 ! muffin-tins use the {\tt bfcmt} vectors in the {\tt atoms} block. Collinear
286 ! calculations are more efficient if the field is applied in the
289 ! \subsection{{\tt chgexs}}
290 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
292 ! {\tt chgexs } & excess electronic charge & real & $0.0$ \\
294 ! \end{tabularx}\newline\newline
295 ! This controls the amount of charge in the unit cell beyond that required to
296 ! maintain neutrality. It can be set positive or negative depending on whether
297 ! electron or hole doping is required.
299 ! \subsection{{\tt deband}}
300 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
302 ! {\tt deband} & initial band energy step size & real & $0.0025$ \\
304 ! \end{tabularx}\newline\newline
305 ! The initial step length used when searching for the band energy, which is
306 ! used as the APW linearisation energy. This is done by first searching
307 ! upwards in energy until the radial wavefunction at the muffin-tin radius is
308 ! zero. This is the energy at the top of the band, denoted $E_{\rm t}$. A
309 ! downward search is now performed from $E_{\rm t}$ until the slope of the
310 ! radial wavefunction at the muffin-tin radius is zero. This energy,
311 ! $E_{\rm b}$, is at the bottom of the band. The band energy is taken as
312 ! $(E_{\rm t}+E_{\rm b})/2$. If either $E_{\rm t}$ or $E_{\rm b}$ cannot be
313 ! found then the band energy is set to the default value.
315 ! \subsection{{\tt deltaem}}
316 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
318 ! {\tt deltaem} & the size of the ${\bf k}$-vector displacement used when
319 ! calculating numerical derivatives for the effective mass tensor & real &
322 ! \end{tabularx}\newline\newline
323 ! See {\tt ndspem} and {\tt vklem}.
325 ! \subsection{{\tt deltaph}}
326 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
328 ! {\tt deltaph} & the size of the atomic displacement used for calculating
329 ! dynamical matrices & real & $0.03$ \\
331 ! \end{tabularx}\newline\newline
332 ! Phonon calculations are performed by constructing a supercell corresponding
333 ! to a particular ${\bf q}$-vector and making a small periodic displacement of
334 ! the atoms. The magnitude of this displacement is given by {\tt deltaph}.
335 ! This should not be made too large, as anharmonic terms could then become
336 ! significant, neither should it be too small as this can introduce numerical
339 ! \subsection{{\tt dos}}
340 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
342 ! {\tt nwdos} & number of frequency/energy points in the DOS or optics plot &
344 ! {\tt ngrdos} & effective {\bf k}-point mesh size to be used for Brillouin
345 ! zone integration & integer & $100$ \\
346 ! {\tt nsmdos} & level of smoothing applied to DOS/optics output & integer &
349 ! {\tt wintdos} & frequency/energy window for the DOS or optics plot &
350 ! real(2) & $(-0.5,0.5)$ \\
352 ! \end{tabularx}\newline\newline
353 ! DOS and optics plots require integrals of the kind
354 ! $$ g(\omega_i)=\frac{\Omega}{(2\pi)^3}\int_{\rm BZ} f({\bf k})
355 ! \delta(\omega_i-e({\bf k}))d{\bf k}. $$
356 ! These are calculated by first interpolating the functions $e({\bf k})$ and
357 ! $f({\bf k})$ with the trilinear method on a much finer mesh whose size is
358 ! determined by {\tt ngrdos}. Then the $\omega$-dependent histogram of the
359 ! integrand is accumulated over the fine mesh. If the output function is noisy
360 ! then either {\tt ngrdos} should be increased or {\tt nwdos} decreased.
361 ! Alternatively, the output function can be artificially smoothed up to a
362 ! level given by {\tt nsmdos}. This is the number of successive 3-point
363 ! averages to be applied to the function $g$.
365 ! \subsection{{\tt dtauoep}}
366 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
368 ! {\tt dtauoep} & step length increment for the exact exchange iterative
369 ! solver & real & $0.5$ \\
371 ! \end{tabularx}\newline\newline
372 ! See {\tt maxitoep} and {\tt tau0oep}
374 ! \subsection{{\tt epschg}}
375 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
377 ! {\tt epschg} & maximum allowed error in the calculated total charge beyond
378 ! which a warning message will be issued & real & $1\times 10^{-3}$ \\
382 ! \subsection{{\tt epsforce}}
383 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
385 ! {\tt epsforce} & convergence tolerance for the forces during a structural
386 ! optimisation run & real & $5\times 10^{-4}$ \\
388 ! \end{tabularx}\newline\newline
389 ! If the mean absolute value of the atomic forces is less than {\tt epsforce}
390 ! then the structural optimisation run is ended. See {\tt tasks}.
392 ! \subsection{{\tt epslat}}
393 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
395 ! {\tt epslat } & vectors with lengths less than this are considered zero &
396 ! real & $10^{-6}$ \\
398 ! \end{tabularx}\newline\newline
399 ! Sets the tolerance for determining if a vector or its components are zero.
400 ! This is to account for any numerical error in real or reciprocal space
403 ! \subsection{{\tt epsocc}}
404 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
406 ! {\tt epsocc} & smallest occupancy for which a state will contribute to the
407 ! density & real & $1\times 10^{-8}$ \\
411 ! \subsection{{\tt epspot}}
412 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
414 ! {\tt epspot} & convergence criterion for the effective potential and field &
415 ! real & $1\times 10^{-6}$ \\
417 ! \end{tabularx}\newline\newline
418 ! If the RMS change in the effective potential and magnetic field is smaller
419 ! than {\tt epspot}, then the self-consistent loop is considered converged
420 ! and exited. For structural optimisation runs this results in the forces
421 ! being calculated, the atomic positions updated and the loop restarted. See
424 ! \subsection{{\tt evalmin}}
425 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
427 ! {\tt evalmin} & valence eigenvalue minimum & real & $-4.5$ \\
429 ! \end{tabularx}\newline\newline
430 ! Any valence states with eigenvalues below {\tt evalmin} are not occupied and
431 ! a warning message is issued.
433 ! \subsection{{\tt fixspin}}
434 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
436 ! {\tt fixspin} & 0 for no fixed spin moment (FSM), 1 for total FSM, 2 for
437 ! local muffin-tin FSM, and 3 for both total and local FSM & integer & 0 \\
439 ! \end{tabularx}\newline\newline
440 ! Set to 1, 2 or 3 for fixed spin moment calculations. See also
441 ! {\tt momfix}, {\tt mommtfix}, {\tt taufsm} and {\tt spinpol}.
443 ! \subsection{{\tt fracinr}}
444 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
446 ! {\tt fracinr} & fraction of the muffin-tin radius up to which {\tt lmaxinr}
447 ! is used as the angular momentum cut-off & real & $0.25$ \\
449 ! \end{tabularx}\newline\newline
452 ! \subsection{{\tt gmaxvr}}
453 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
455 ! {\tt gmaxvr} & maximum length of $|{\bf G}|$ for expanding the interstitial
456 ! density and potential & real & $12.0$ \\
458 ! \end{tabularx}\newline\newline
459 ! See also {\tt rgkmax}.
461 ! \subsection{{\tt intraband}}
462 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
464 ! {\tt intraband} & {\tt .true.} if the intraband (Drude-like) contribution is
465 ! to be added to the dieletric tensor & logical & {\tt .false.} \\
467 ! \end{tabularx}\newline\newline
469 ! \subsection{{\tt kstlist}}
470 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
472 ! {\tt kstlist(i)} & $i$th ${\bf k}$-point and state pair & integer(2) &
475 ! \end{tabularx}\newline\newline
476 ! This is a user-defined list of ${\bf k}$-point and state index pairs which
477 ! are those used for plotting wavefunctions and writing ${\bf L}$, ${\bf S}$
478 ! and ${\bf J}$ expectation values. Only the first pair is used by the
479 ! aforementioned tasks. The list should be terminated by a blank line.
481 ! \subsection{{\tt lda+u}}
482 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
484 ! {\tt ldapu} & type of LDA+$U$ calculation & integer & 0 \\
486 ! {\tt is} & species number & integer & - \\
487 ! {\tt l} & angular momentum value & integer & -1 \\
488 ! {\tt u} & the desired $U$ value & real & $0.0$ \\
489 ! {\tt j} & the desired $J$ value & real & $0.0$ \\
491 ! \end{tabularx}\newline\newline
492 ! This block contains the parameters required for an LDA+$U$ calculation, with
493 ! the list of parameters for each species terminated with a blank line. The
494 ! type of calculation required is set with the parameter {\tt ldapu}.
495 ! Currently implemented are:\newline\newline
496 ! \begin{tabularx}{\textwidth}[h]{lX}
497 ! 0 & No LDA+$U$ calculation \\
498 ! 1 & Fully localised limit (FLL) \\
499 ! 2 & Around mean field (AFM) \\
500 ! 3 & An interpolation between FLL and AFM \\
501 ! \end{tabularx}\newline\newline
502 ! See (amongst others) Phys. Rev. B {\bf 67}, 153106 (2003), Phys. Rev. B
503 ! {\bf 52}, R5467 (1995), and Phys. Rev. B {\bf 60}, 10673 (1999).
505 ! \subsection{{\tt lmaxapw}}
506 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
508 ! {\tt lmaxapw} & angular momentum cut-off for the APW functions & integer &
513 ! \subsection{{\tt lmaxinr}}
514 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
516 ! {\tt lmaxinr} & angular momentum cut-off for themuffin-tin density and
517 ! potential on the inner part of the muffin-tin & integer & 2 \\
519 ! \end{tabularx}\newline\newline
520 ! Close to the nucleus, the density and potential is almost spherical and
521 ! therefore the spherical harmonic expansion can be truncated a low angular
522 ! momentum. See also {\tt fracinr}.
524 ! \subsection{{\tt lmaxmat}}
525 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
527 ! {\tt lmaxmat} & angular momentum cut-off for the outer-most loop in the
528 ! hamiltonian and overlap matrix setup & integer & 5 \\
532 ! \subsection{{\tt lmaxvr}}
533 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
535 ! {\tt lmaxvr} & angular momentum cut-off for the muffin-tin density and
536 ! potential & integer & 7 \\
540 ! \subsection{{\tt lradstp}}
541 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
543 ! {\tt lradstp} & radial step length for determining coarse radial mesh &
546 ! \end{tabularx}\newline\newline
547 ! Some muffin-tin functions (such as the density) are calculated on a coarse
548 ! radial mesh and then interpolated onto a fine mesh. This is done for the
549 ! sake of efficiency. {\tt lradstp} defines the step size in going from the
550 ! fine to the coarse radial mesh. If it is too large, loss of precision may
553 ! \subsection{{\tt maxitoep}}
554 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
556 ! {\tt maxitoep} & maximum number of iterations when solving the exact
557 ! exchange integral equations & integer & 120 \\
559 ! \end{tabularx}\newline\newline
560 ! See {\tt tau0oep} and {\tt dtauoep}.
562 ! \subsection{{\tt maxscl}}
563 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
565 ! {\tt maxscl } & maximum number of self-consistent loops allowed & integer &
568 ! \end{tabularx}\newline\newline
569 ! This determines after how many loops the self-consistent cycle will
570 ! terminate if the convergence criterion is not met. If {\tt maxscl} is $1$
571 ! then the density and potential file, {\tt STATE.OUT}, will {\bf not} be
572 ! written to disk at the end of the loop. See {\tt epspot}.
574 ! \subsection{{\tt molecule}}
575 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
577 ! {\tt molecule} & {\tt .true.} if the system is an isolated molecule &
578 ! logical & {\tt .false.} \\
580 ! \end{tabularx}\newline\newline
581 ! If {\tt molecule} is {\tt .true.}, then the atomic positions, ${\bf a}$,
582 ! given in the {\tt atoms} block are assumed to be in Cartesian coordinates.
583 ! The lattice vectors are also set up automatically with the $i$th lattice
585 ! $$ {\bf A}^i=A_i\hat{\bf e}^i, $$
587 ! $$ A_i=\max_{\alpha,\beta}\left|{\bf a}^{\alpha}_i-{\bf a}^{\beta}_i\right|
589 ! with $\alpha$ and $\beta$ labeling atoms, and $d_{\rm vac}$ determines the
590 ! size of the vacuum around the molecule. The last variable is set by the
591 ! input parameter {\tt vacuum}.
593 ! \subsection{{\tt momfix}}
594 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
596 ! {\tt momfix} & the desired total moment for a FSM calculation &
597 ! real(3) & $(0.0,0.0,0.0)$ \\
599 ! \end{tabularx}\newline\newline
600 ! Note that all three components must be specified (even for collinear
601 ! calculations). See {\tt fixspin}, {\tt taufsm} and {\tt spinpol}.
603 ! \subsection{{\tt mommtfix}}
604 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
606 ! {\tt is} & species number & integer & 0 \\
607 ! {\tt ia} & atom number & integer & 0 \\
608 ! {\tt mommtfix} & the desired muffin-tin moment for a FSM calculation &
609 ! real(3) & $(0.0,0.0,0.0)$ \\
611 ! \end{tabularx}\newline\newline
612 ! The local muffin-tin moments are specified for a subset of atoms, with the
613 ! list terminated with a blank line. Note that all three components must be
614 ! specified (even for collinear calculations). See {\tt fixspin}, {\tt taufsm}
617 ! \subsection{{\tt ndspem}}
618 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
620 ! {\tt ndspem} & the number of {\bf k}-vector displacements in each direction
621 ! around {\tt vklem} when computing the numerical derivatives for the
622 ! effective mass tensor & integer & 1 \\
624 ! \end{tabularx}\newline\newline
625 ! See {\tt deltaem} and {\tt vklem}.
627 ! \subsection{{\tt nempty}}
628 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
630 ! {\tt nempty} & the number of empty states & integer & 5 \\
632 ! \end{tabularx}\newline\newline
633 ! Defines the number of eigenstates beyond that required for charge
634 ! neutrality. When running metals it is not known {\it a priori} how many
635 ! states will be below the Fermi energy for each {\bf k}-point. Setting
636 ! {\tt nempty} greater than zero allows the additional states to act as a
637 ! buffer in such cases. Furthermore, magnetic calculations use the
638 ! first-variational eigenstates as a basis for setting up the
639 ! second-variational Hamiltonian, and thus {\tt nempty} will determine the
640 ! size of this basis set. Convergence with respect to this quantity should be
643 ! \subsection{{\tt ngridk}}
644 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
646 ! {\tt ngridk } & the {\bf k}-point mesh sizes & integer(3) & $(1,1,1)$ \\
648 ! \end{tabularx}\newline\newline
650 ! The ${\bf k}$-vectors are generated using
651 ! $$ {\bf k}=(\frac{i_1}{n_1},\frac{i_2}{n_2},\frac{i_3}{n_3})
652 ! +{\bf v}_{\rm off}, $$
653 ! where $i_j$ runs from 0 to $n_j-1$ and $0\le{\bf v}_{{\rm off};j}<1$ for
654 ! $j=1,2,3$. See also {\tt reducek} and {\tt vkloff}.
656 ! \subsection{{\tt ngridq}}
657 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
659 ! {\tt ngridq } & the phonon {\bf q}-point mesh sizes & integer(3) &
662 ! \end{tabularx}\newline\newline
663 ! Same as {\tt ngridk}, except that this mesh is for the phonon
664 ! {\bf q}-points. See also {\tt reduceq}.
666 ! \subsection{{\tt nosource}}
667 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
669 ! {\tt nosource} & when set to {\tt .true.}, source fields are projected out
670 ! of the exchange-correlation magnetic field & logical & {\tt .false.} \\
672 ! \end{tabularx}\newline\newline
673 ! Experimental feature.
675 ! \subsection{{\tt nosym}}
676 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
678 ! {\tt nosym} & when set to {\tt .true.} no symmetries, apart from the
679 ! identity, are used anywhere in the code & logical & {\tt .false.} \\
683 ! \subsection{{\tt notes}}
684 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
686 ! {\tt notes(i)} & the $i$th line of the notes & string & - \\
688 ! \end{tabularx}\newline\newline
689 ! This block allows users to add their own notes to the file {\tt INFO.OUT}.
690 ! The block should be terminated with a blank line, and no line should exceed
693 ! \subsection{{\tt nprad}}
694 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
696 ! {\tt nprad} & radial polynomial order & integer & 4 \\
698 ! \end{tabularx}\newline\newline
699 ! This sets the polynomial order for the predictor-corrector method when
700 ! solving the radial Dirac and Schr\"odinger equations, as well as for
701 ! performing radial interpolation in the plotting routines.
703 ! \subsection{{\tt nstfsp}}
704 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
706 ! {\tt nstfsp} & number of states to be included in the Fermi surface plot
707 ! file & integer & 6 \\
711 ! \subsection{{\tt nwrite}}
712 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
714 ! {\tt nwrite} & number of iterations after which {\tt STATE.OUT} is to be
715 ! written & integer & 0 \\
717 ! \end{tabularx}\newline\newline
718 ! Normally, the density and potentials are written to the file {\tt STATE.OUT}
719 ! only after completion of the self-consistent loop. By setting {\tt nwrite}
720 ! to a positive integer the file will be written during the loop every
721 ! {\tt nwrite} iterations.
723 ! \subsection{{\tt optcomp}}
724 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
726 ! {\tt optcomp} & the components of the first- or second-order optical tensor
727 ! to be calculated & integer(3) & $(1,1,1)$ \\
729 ! \end{tabularx}\newline\newline
730 ! This selects which components of the optical tensor you would like to plot.
731 ! Only the first two are used for the first-order tensor.
733 ! \subsection{{\tt phwrite}}
734 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
736 ! {\tt nphwrt} & number of {\bf q}-points for which phonon modes are to be
737 ! found & integer & 1 \\
739 ! {\tt vqlwrt(i)} & the $i$th {\bf q}-point in lattice coordinates & real(3) &
742 ! \end{tabularx}\newline\newline
743 ! This is used in conjunction with {\tt task}=230. The code will write the
744 ! phonon frequencies and eigenvectors to the file {\tt PHONON.OUT} for all the
745 ! {\bf q}-points in the list. The {\bf q}-points can be anywhere in the
746 ! Brillouin zone and do not have to lie on the mesh defined by {\tt ngridq}.
747 ! Obviously, all the dynamical matrices have to be computed first using
750 ! \subsection{{\tt plot1d}}
751 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
753 ! {\tt nvp1d} & number of vertices & integer & 2 \\
754 ! {\tt npp1d} & number of plotting points & integer & 200 \\
756 ! {\tt vvlp1d(i)} & lattice coordinates for vertex {\tt i} & real(3) &
757 ! $(0.0,0.0,0.0)\rightarrow(1.0,1.0,1.0)$ \\
759 ! \end{tabularx}\newline\newline
760 ! Defines the path in either real or reciprocal space along which the 1D plot
761 ! is to be produced. The user should provide {\tt nvp1d} vertices in lattice
764 ! \subsection{{\tt plot2d}}
765 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
767 ! {\tt vclp2d(1)} & first corner (origin) & real(3) & $(0.0,0.0,0.0)$ \\
769 ! {\tt vclp2d(2)} & second corner & real(3) & $(1.0,0.0,0.0)$ \\
771 ! {\tt vclp2d(3)} & third corner & real(3) & $(0.0,1.0,0.0)$ \\
773 ! {\tt np2d} & number of plotting points in both directions & integer(2) &
776 ! \end{tabularx}\newline\newline
777 ! Defines corners of the parallelogram and the mesh size used for producing 2D
780 ! \subsection{{\tt plot3d}}
781 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
783 ! {\tt nup3d} & number of unit cells to plot & integer(3) & $(1,1,1)$ \\
785 ! {\tt np3d} & number of plotting points each direction & integer(3) &
788 ! \end{tabularx}\newline\newline
789 ! Defines the number of unit cells in each direction to be plotted in 3D as
790 ! well as the size of the plotting mesh. The {\tt nup3d} parameter is also
791 ! used to define the number of reciprocal lattice unit cells to be plotted for
792 ! Fermi surface plots.
794 ! \subsection{{\tt primcell}}
795 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
797 ! {\tt primcell} & {\tt .true.} if the primitive unit cell should be found
798 ! & logical & {\tt .false.} \\
800 ! \end{tabularx}\newline\newline
801 ! Allows the primitive unit cell to be determined automatically from the
802 ! conventional cell. This is done by searching for lattice vectors among all
803 ! those which connect atomic sites, and using the three shortest which produce
804 ! a unit cell with non-zero volume.
806 ! \subsection{{\tt reducebf}}
807 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
809 ! {\tt reducebf} & reduction factor for the external magnetic fields & real &
812 ! \end{tabularx}\newline\newline
813 ! After each iteration the external magnetic fields are multiplied with
814 ! {\tt reducebf}. This allows for a large external magnetic field at the start
815 ! of the self-consistent loop to break spin symmetry, while at the end of the
816 ! loop the field will be effectively zero, i.e. infinitesimal. See
817 ! {\tt bfieldc} and {\tt atoms}.
819 ! \subsection{{\tt reducek}}
820 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
822 ! {\tt reducek} & set to {\tt .true.} if the ${\bf k}$-point set is to be
823 ! reduced with the crystal symmetries & logical & {\tt .true.} \\
825 ! \end{tabularx}\newline\newline
826 ! See also {\tt ngridk} and {\tt vkloff}.
828 ! \subsection{{\tt reduceq}}
829 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
831 ! {\tt reduceq} & set to {\tt .true.} if the ${\bf q}$-point set is to be
832 ! reduced with the crystal symmetries & logical & {\tt .true.} \\
834 ! \end{tabularx}\newline\newline
835 ! See also {\tt ngridq}.
837 ! \subsection{{\tt rgkmax}}
838 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
840 ! {\tt rgkmax} & $R^{\rm MT}_{\rm min}\times\max(|{\bf G}+{\bf k}|)$ & real &
843 ! \end{tabularx}\newline\newline
844 ! This sets the maximum length for the ${\bf G}+{\bf k}$ vectors, defined as
845 ! {\tt rgkmax} divided by the smallest muffin-tin radius.
847 ! \subsection{{\tt rlambda}}
848 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
850 ! {\tt rlambda } & maximum de Broglie wavelength of {\bf k}-vectors & real &
853 ! \end{tabularx}\newline\newline
854 ! Used for the automatic determination of the {\bf k}-point mesh. If
855 ! {\tt autokpt} is set to {\tt .true.} then the mesh sizes will be determined
856 ! by $n_i=\lambda/|{\bf A}_i|+1$.
858 ! \subsection{{\tt rmtapm}}
859 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
861 ! {\tt rmtapm } & parameters governing the automatic generation of the
862 ! muffin-tin radii & real(2) & $(0.25, 0.95)$ \\
864 ! \end{tabularx}\newline\newline
865 ! When {\tt autormt} is set to true, the muffin-tin radii are found
866 ! automatically from the formula
867 ! $$ R_i\propto 1+\zeta|Z_i|^{1/3}, $$
868 ! where $Z_i$ is the atomic number of the $i$th species, $\zeta$ is stored in
869 ! {\tt rmtapm(1)} and the value which governs the distance between the
870 ! muffin-tins is stored in {\tt rmtapm(2)}. When {\tt rmtapm(2)} $=1$, the
871 ! closest muffin-tins will touch.
873 ! \subsection{{\tt scale}}
874 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
876 ! {\tt scale } & lattice vector scaling factor & real & $1.0$ \\
878 ! \end{tabularx}\newline\newline
879 ! Scaling factor for all three lattice vectors. Applied in conjunction with
880 ! {\tt scale1}, {\tt scale2} and {\tt scale3}.
882 ! \subsection{{\tt scale1/2/3}}
883 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
885 ! {\tt scale1/2/3 } & separate scaling factors for each lattice vector &
890 ! \subsection{{\tt scissor}}
891 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
893 ! {\tt scissor} & the scissors correction & real & $0.0$ \\
895 ! \end{tabularx}\newline\newline
896 ! This is the scissors shift applied to states above the Fermi energy. Affects
897 ! DOS, optics and band structure plots.
899 ! \subsection{{\tt scrpath}}
900 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
902 ! {\tt scrpath} & scratch space path & string & {\tt ./} \\
904 ! \end{tabularx}\newline\newline
905 ! This is the path to scratch space where the eigenvector file
906 ! {\tt EIGVEC.OUT} will be written. If the local directory is accessed via a
907 ! network then {\tt scrpath} can be set to a directory on the local disk, for
908 ! example {\tt /tmp/}. Note that the forward slash {\tt /} at the end of the
909 ! string must be included.
911 ! \subsection{{\tt spinorb}}
912 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
914 ! {\tt spinorb} & set to {\tt .true.} if a spin-orbit coupling is required
915 ! & logical & {\tt .false.} \\
917 ! \end{tabularx}\newline\newline
918 ! If {\tt spinorb} is {\tt .true.}, then a $\boldsymbol\sigma\cdot{\bf L}$
919 ! term is added to the second-variational Hamiltonian. See {\tt spinpol}.
921 ! \subsection{{\tt spinpol}}
922 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
924 ! {\tt spinpol} & set to {\tt .true.} if a spin-polarised calculation is
925 ! required & logical & {\tt .false.} \\
927 ! \end{tabularx}\newline\newline
928 ! If {\tt spinpol} is {\tt .true.}, then the spin-polarised Hamiltonian is
929 ! solved as a second-variational step using two-component spinors in the
930 ! effective magnetic field. The first variational scalar wavefunctions are
931 ! used as a basis for setting this Hamiltonian.
933 ! \subsection{{\tt spinsprl}}
934 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
936 ! {\tt spinsprl} & set to {\tt .true.} if a spin-spiral calculation is
937 ! required & logical & {\tt .false.} \\
939 ! \end{tabularx}\newline\newline
940 ! Experimental feature for the calculation of spin-spiral states. See
941 ! {\tt vqlss} for details.
943 ! \subsection{{\tt sppath}}
944 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
946 ! {\tt sppath} & path where the species files can be found & string &
949 ! \end{tabularx}\newline\newline
950 ! Note that the forward slash {\tt /} at the end of the string must be
953 ! \subsection{{\tt stype}}
954 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
956 ! {\tt stype} & integer defining the type of smearing to be used & integer &
959 ! \end{tabularx}\newline\newline
960 ! A smooth approximation to the Dirac delta function is needed to compute the
961 ! occupancies of the Kohn-Sham states. The variable {\tt swidth} determines
962 ! the width of the approximate delta function. Currently implemented are
964 ! \begin{tabularx}{\textwidth}[h]{lX}
966 ! 1 & Methfessel-Paxton order 1, Phys. Rev. B {\bf 40}, 3616 (1989) \\
967 ! 2 & Methfessel-Paxton order 2 \\
971 ! \subsection{{\tt swidth}}
972 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
974 ! {\tt swidth} & width of the smooth approximation to the Dirac delta
975 ! function & real & $0.01$ \\
977 ! \end{tabularx}\newline\newline
978 ! See {\tt stype} for details.
980 ! \subsection{{\tt tasks}}
981 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
983 ! {\tt task(i) } & the $i$th task & integer & $-1$ \\
985 ! \end{tabularx}\newline\newline
986 ! A list of tasks for the code to perform sequentially. The list should be
987 ! terminated with a blank line. Each task has an associated integer as
988 ! follows:\newline\newline
989 ! \begin{tabularx}{\textwidth}[h]{lX}
990 ! -1 & Write out the version number of the code. \\
991 ! 0 & Ground state run starting from the atomic densities. \\
992 ! 1 & Resumption of ground state run using density in {\tt STATE.OUT}. \\
993 ! 2 & Structural optimisation run starting from the atomic densities, with
994 ! atomic positions written to {\tt GEOMETRY.OUT}. \\
995 ! 3 & Resumption of structural optimisation run using density in
996 ! {\tt STATE.OUT} but with positions from {\tt exciting.in}. \\
997 ! 5 & Ground state Hartree-Fock run (experimental feature). \\
998 ! 10 & Total, partial and interstitial density of states (DOS). \\
999 ! 15 & Output ${\bf L}$, ${\bf S}$ and ${\bf J}$ total expectation values. \\
1000 ! 16 & Output ${\bf L}$, ${\bf S}$ and ${\bf J}$ expectation values for each
1001 ! {\bf k}-point and state in {\tt kstlist}. \\
1002 ! 20 & Band structure plot. \\
1003 ! 21 & Band structure plot which includes angular momentum characters for
1005 ! 25 & Compute the effective mass tensor at the {\bf k}-point given by
1007 ! 31, 32, 33 & 1/2/3D charge density plot. \\
1008 ! 41, 42, 43 & 1/2/3D exchange-correlation and Coulomb potential plots. \\
1009 ! 51, 52, 53 & 1/2/3D electron localisation function (ELF) plot. \\
1010 ! 61, 62, 63 & 1/2/3D wavefunction plot:
1011 ! $\left|\Phi_{i{\bf k}}({\bf r})\right|^2$. \\
1012 ! 72, 73 & 2/3D plot of magnetisation vector field, ${\bf m}({\bf r})$. \\
1013 ! 82, 83 & 2/3D plot of exchange-correlation magnetic vector field,
1014 ! ${\bf B}_{\rm xc}({\bf r})$. \\
1015 ! 91, 92, 93 & 1/2/3D plot of $\nabla\cdot{\bf B}_{\rm xc}({\bf r})$. \\
1016 ! 100 & 3D Fermi surface plot using the scalar product
1017 ! $p({\bf k})=\Pi_i(\epsilon_{i{\bf k}}-\epsilon_{\rm F})$. \\
1018 ! 101 & 3D Fermi surface plot using separate bands (minus the Fermi
1020 ! 110 & Calculation of M\"{o}ssbauer contact charge densities and magnetic
1021 ! fields at the nuclear sites. \\
1022 ! 115 & Calculation of the electric field gradient (EFG) at the nuclear
1024 ! 120 & Output of the momentum matrix elements
1025 ! $\langle\Phi_{i{\bf k}}|-i\nabla|\Phi_{j{\bf k}}\rangle$. \\
1026 ! 121 & Linear optical response tensor. \\
1027 ! 122 & Magneto optical Kerr effect angle. \\
1028 ! 142, 143 & 2/3D plot of the electric field
1029 ! ${\bf E}({\bf r})\equiv\nabla V_{\rm C}({\bf r})$. \\
1030 ! 152, 153 & 2/3D plot of
1031 ! ${\bf m}({\bf r})\times{\bf B}_{\rm xc}({\bf r})$. \\
1032 ! 162 & Scanning-tunneling microscopy (STM) image. \\
1033 ! 200 & Calculation of dynamical matrices on a {\bf q}-point set defined by
1035 ! 210 & Phonon density of states. \\
1036 ! 220 & Phonon dispersion plot. \\
1037 ! 230 & Phonon frequencies and eigenvectors for an arbitrary
1038 ! ${\bf q}$-point. \\
1039 ! 250 & Write the atomic geometry to file for plotting with {\sf XCrySDen}
1043 ! \subsection{{\tt tau0atm}}
1044 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1046 ! {\tt tau0atm} & the step size to be used for structural optimisation &
1049 ! \end{tabularx}\newline\newline
1050 ! The position of atom $\alpha$ is updated on step $m$ of a structural
1051 ! optimisation run using
1052 ! $$ {\bf r}_{\alpha}^{m+1}={\bf r}_{\alpha}^m+\tau_{\alpha}^m
1053 ! \left({\bf F}_{\alpha}^m+{\bf F}_{\alpha}^{m-1}\right), $$
1054 ! where $\tau_{\alpha}$ is set to {\tt tau0atm} for $m=0$, and incremented by
1055 ! the same amount if the atom is moving in the same direction between steps.
1056 ! If the direction changes then $\tau_{\alpha}$ is reset to {\tt tau0atm}.
1058 ! \subsection{{\tt tau0oep}}
1059 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1061 ! {\tt tau0oep} & initial step length for the exact exchange iterative
1062 ! solver & real & $0.5$ \\
1064 ! \end{tabularx}\newline\newline
1065 ! See {\tt maxitoep} and {\tt dtauoep}.
1067 ! \subsection{{\tt taufsm}}
1068 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1070 ! {\tt taufsm} & the step size to be used when finding the effective magnetic
1071 ! field in fixed spin moment calculations & real & $0.01$ \\
1073 ! \end{tabularx}\newline\newline
1074 ! An effective magnetic field, ${\bf B}_{\rm FSM}$, is required for fixing the
1075 ! spin moment to a given value, $\boldsymbol{\mu}_{\rm FSM}$. This is found by
1076 ! adding a vector to the field which is proportional to the difference between
1077 ! the moment calculated in the $i$th self-consistent loop and the required
1079 ! $$ {\bf B}_{\rm FSM}^{i+1}={\bf B}_{\rm FSM}^i+\lambda\left(
1080 ! \boldsymbol{\mu}^i-\boldsymbol{\mu}_{\rm FSM}\right), $$
1081 ! where $\lambda$ is proportional to {\tt taufsm}. See also {\tt fixspin},
1082 ! {\tt momfix} and {\tt spinpol}.
1084 ! \subsection{{\tt tfibs}}
1085 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1087 ! {\tt tfibs} & set to {\tt .true.} if the IBS correction to the force should
1088 ! be calculated & logical & {\tt .true.} \\
1090 ! \end{tabularx}\newline\newline
1091 ! Because calculation of the incomplete basis set (IBS) correction to the
1092 ! force is fairly time-consuming, it can be switched off by setting
1093 ! {\tt tfibs} to {\tt .false.} This correction can then be included only when
1094 ! necessary, i.e. when the atoms are close to equilibrium in a structural
1097 ! \subsection{{\tt tforce}}
1098 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1100 ! {\tt tforce} & set to {\tt .true.} if the force should be calculated at the
1101 ! end of the self-consistent cycle & logical & {\tt .false.} \\
1103 ! \end{tabularx}\newline\newline
1104 ! This variable is automatically set to {\tt .true.} when performing
1105 ! structural optimisation.
1107 ! \subsection{{\tt tshift}}
1108 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1110 ! {\tt tshift} & set to {\tt .true.} if the crystal can be shifted so that the
1111 ! atom closest to the origin is exactly at the origin &
1112 ! logical & {\tt .true.} \\
1116 ! \subsection{{\tt usegdft}}
1117 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1119 ! {\tt usegdft} & set to {\tt .true.} if the generalised DFT correction of
1120 ! L. Fritsche and Y. M. Gu, Phys. Rev. {\bf B} 48, 4250 (1993), is to be
1121 ! used & logical & {\tt .false.} \\
1123 ! \end{tabularx}\newline\newline
1124 ! Experimental feature.
1126 ! \subsection{{\tt vacuum}}
1127 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1129 ! {\tt vacuum} & the size of the vacuum region around a molecule & real &
1132 ! \end{tabularx}\newline\newline
1133 ! See {\tt molecule}.
1135 ! \subsection{{\tt vklem}}
1136 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1138 ! {\tt vklem} & the ${\bf k}$-point in lattice coordinates at which to compute
1139 ! the effective mass tensors & real(3) & $(0.0,0.0,0.0)$ \\
1141 ! \end{tabularx}\newline\newline
1142 ! See {\tt deltaem} and {\tt ndspem}.
1144 ! \subsection{{\tt vqlss}}
1145 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1147 ! {\tt vqlss} & the ${\bf q}$-vector of the spin-spiral state in lattice
1148 ! coordinates & real(3) & $(0.0,0.0,0.0)$ \\
1150 ! \end{tabularx}\newline\newline
1151 ! Spin-spirals arise from spinor states assumed to be of the form
1152 ! $$ \Psi^{\bf q}_{\bf k}({\bf r})=
1153 ! \left( \begin{array}{c}
1154 ! U^{{\bf q}\uparrow}_{\bf k}({\bf r})e^{i({\bf k+q/2})\cdot{\bf r}} \\
1155 ! U^{{\bf q}\downarrow}_{\bf k}({\bf r})e^{i({\bf k-q/2})\cdot{\bf r}} \\
1156 ! \end{array} \right). $$
1157 ! These are determined using a second-variational approach, and give rise to a
1158 ! magnetisation density of the form
1159 ! $$ {\bf m}^{\bf q}({\bf r})=(m_x({\bf r})\cos({\bf q \cdot r}),
1160 ! m_y({\bf r})\sin({\bf q \cdot r}),m_z({\bf r})), $$
1161 ! where $m_x$, $m_y$ and $m_z$ are lattice periodic. See also {\tt spinprl}.
1163 ! \subsection{{\tt vkloff}}
1164 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1166 ! {\tt vkloff } & the ${\bf k}$-point offset vector in lattice coordinates &
1167 ! real(3) & $(0.0,0.0,0.0)$ \\
1169 ! \end{tabularx}\newline\newline
1172 ! \subsection{{\tt xctype}}
1173 ! \begin{tabularx}{\textwidth}[h]{|l|X|c|c|}
1175 ! {\tt xctype} & integer defining the type of exchange-correlation functional
1176 ! to be used & integer & 3 \\
1178 ! \end{tabularx}\newline\newline
1179 ! Currently implemented are:\newline\newline
1180 ! \begin{tabularx}{\textwidth}[h]{lX}
1181 ! 1 & No exchange-correlation funtional ($E_{\rm xc}\equiv 0$) \\
1182 ! 2 & LDA, Perdew-Zunger/Ceperley-Alder, {\it Phys. Rev. B} {\bf 23}, 5048
1184 ! 3 & LSDA, Perdew-Wang/Ceperley-Alder, {\it Phys. Rev. B} {\bf 45}, 13244
1186 ! 4 & LDA, X-alpha approximation, J. C. Slater, {\it Phys. Rev.} {\bf 81}, 385
1188 ! 20 & GGA, Perdew-Burke-Ernzerhof, {\it Phys. Rev. Lett.} {\bf 77}, 3865
1190 ! 21 & GGA, Revised PBE, Zhang-Yang, {\it Phys. Rev. Lett.} {\bf 80}, 890
1192 ! 22 & GGA, PBEsol, arXiv:0707.2088v1 (2007) \\
1193 ! 26 & GGA, Wu-Cohen exchange (WC06) with PBE correlation, {\it Phys. Rev. B}
1194 ! {\bf 73}, 235116 (2006) \\
1195 ! 30 & GGA, Armiento-Mattsson (AM05) spin-unpolarised functional,
1196 ! {\it Phys. Rev. B} {\bf 72}, 085108 (2005) \\
1199 ! \section{Contributing to {\sf EXCITING}}
1200 ! Please bear in mind when writing code for the {\sf EXCITING} project that
1201 ! it should be an exercise in physics and not software engineering. All code
1202 ! should therefore be kept as simple and concise as possible, and above all it
1203 ! should be easy for anyone to locate and follow the {\sf Fortran}
1204 ! representation of the original mathematics. We would also appreciate the
1205 ! following conventions being adhered to:
1207 ! \item Strict {\sf Fortran} 90/95 should be used. Features which are marked
1208 ! as obsolescent in F90/95 should be avoided. These include assigned format
1209 ! specifiers, labeled do-loops, computed goto statements and statement
1211 ! \item Modules should be used in place of common blocks for declaring
1212 ! global variables. Use the existing modules to declare new global variables.
1213 ! \item Any code should be written in lower-case free form style, starting
1214 ! from column one. Try and keep the length of each line to fewer than 80
1215 ! characters using the \& character for line continuation.
1216 ! \item Every function or subroutine, no matter how small, should be in its
1217 ! own file named {\tt routine.f90}, where {\tt routine} is the function or
1218 ! subroutine name. It is recommended that the routines are named so as to
1219 ! make their purpose apparent from the name alone.
1220 ! \item Use of {\tt implicit none} is mandatory. Remember also to define the
1221 ! {\tt intent} of any passed arguments.
1222 ! \item Local allocatable arrays must be deallocated on exit of the routine to
1223 ! prevent memory leakage. Use of automatic arrays should be limited to arrays
1225 ! \item Every function or subroutine must be documented with the {\sf Protex}
1226 ! source code documentation system. This should include a short \LaTeX\
1227 ! description of the algorithms and methods involved. Equations which need to
1228 ! be referenced should be labeled with {\tt routine\_1}, {\tt routine\_2}
1229 ! etc. The authorship of each new piece of code or modification should be
1230 ! indicated in the {\tt REVISION HISTORY} part of the header. See the
1231 ! {\sf Protex} documentation for details.
1232 ! \item Ensure as much as possible that a routine will terminate the program
1233 ! when given improper input instead of continuing with erroneous results.
1234 ! Specifically, functions should have a well-defined domain for which they
1235 ! return accurate results. Input outside that domain should result in an
1236 ! error message and termination.
1237 ! \item Report errors prior to termination with a short description, for
1241 ! write(*,'("Error(readinput): invalid spnst : ",I8)') spnst(is)
1242 ! write(*,'(" for species ",I4)') is
1245 ! \item Wherever possible, real numbers outputted as ASCII data should be
1246 ! formatted with the {\tt G18.10} specifier.
1247 ! \item Avoid redundant or repeated code: check to see if the routine you need
1248 ! already exists, before writing a new one.
1249 ! \item All reading in of ASCII data should be done in the subroutine
1250 ! {\tt readinput}. For binary data, separate routines for reading and writing
1251 ! should be used (for example {\tt writestate} and {\tt readstate}).
1252 ! \item Input file names should be in lowercase and have the extension
1253 ! {\tt .in} . All output file names should be in uppercase with the extension
1255 ! \item All internal units should be atomic. Input and output units should be
1256 ! atomic by default and clearly stated otherwise. Rydbergs should not be used
1257 ! under any circumstances.
1259 ! \subsection{Licensing}
1260 ! Routines which constitute the main part of the code are released under the
1261 ! GNU General Public License (GPL). Library routines are released under the
1262 ! less restrictive GNU Lesser General Public License (LGPL). Both licenses
1263 ! are contained in the file {\tt COPYING}. Any contribution to the code must
1264 ! be licensed at the authors' discretion under either the GPL or LGPL.
1265 ! Author(s) of the code retain the copyrights. Copyright and (L)GPL
1266 ! information must be included at the beginning of every file, and no code
1267 ! will be accepted without this.
1269 ! \bibliographystyle{unsrt}
1270 ! \bibliography{exciting}