The Elk Code
elk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2011 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.
5 
6 ! main routine for the Elk code
7 program elk
8 use modmain
9 use modmpi
10 use modomp
11 use modvars
12 use modramdisk
13 use moddelf
14 implicit none
15 ! local variables
16 logical exist
17 character(64) str
18 ! initialise MPI execution environment
19 call mpi_init(ierror)
20 ! duplicate mpi_comm_world
21 call mpi_comm_dup(mpi_comm_world,mpicom,ierror)
22 ! determine the number of MPI processes
23 call mpi_comm_size(mpicom,np_mpi,ierror)
24 ! determine the local MPI process number
25 call mpi_comm_rank(mpicom,lp_mpi,ierror)
26 ! determine if the local process is the master
27 if (lp_mpi == 0) then
28  mp_mpi=.true.
29  write(str,'("Elk code version ",I0,".",I0,".",I0," started")') version
30  call writebox(6,trim(str))
31 else
32  mp_mpi=.false.
33 end if
34 10 continue
35 ! read input files
36 call readinput
37 ! initialise OpenMP variables
38 call omp_init
39 ! initialise the MKL library
40 call mkl_init
41 if (mp_mpi) then
42  write(*,*)
43  write(*,'("Number of MPI processes : ",I6)') np_mpi
44  write(*,'("Number of OpenMP threads per MPI process : ",I4)') maxthd
45  write(*,'("Total number of threads (MPI x OpenMP) : ",I6)') np_mpi*maxthd
46  write(*,'("Maximum OpenMP nesting level : ",I4)') maxlvl
47  write(*,'("Number of OpenMP threads at first nesting level : ",I4)') maxthd1
48  write(*,'("Number of MKL threads : ",I4)') maxthdmkl
49 ! check if Elk is already running in this directory
50  inquire(file='RUNNING',exist=exist)
51  if (exist) then
52  write(*,*)
53  write(*,'("Info(elk): several copies of Elk may be running in this path")')
54  write(*,'("(this could be intentional, or result from a previous crash,")')
55  write(*,'(" or arise from an incorrect MPI compilation)")')
56  else
57  open(50,file='RUNNING')
58  close(50)
59  end if
60  if (batch) then
61  wrtvars=.true.
62  write(*,*)
63  write(*,'("Batch mode enabled")')
64  end if
65  if (wrtvars) then
66 ! delete the VARIABLES.OUT file
67  call delfile('VARIABLES.OUT')
68 ! write version number to VARIABLES.OUT
69  call writevars('version',nv=3,iva=version)
70  end if
71 else
72  wrtvars=.false.
73 end if
74 ! initialise the RAM disk if required
75 if (ramdisk) then
76  call initrd
77  if (mp_mpi) then
78  write(*,*)
79  write(*,'("RAM disk enabled")')
80  if (.not.wrtdsk) then
81  write(*,*)
82  write(*,'("Info(elk): some direct access files are not written to disk")')
83  end if
84  end if
85 else
86  wrtdsk=.true.
87 end if
88 ! perform the tasks
89 do itask=1,ntasks
91  if (mp_mpi) then
92  write(str,'("Current task : ",I0)') task
93  call writebox(6,trim(str))
94  end if
95 ! increment the batch variables if required
96  if (batch) call batchdv
97 ! check if task can be run with MPI
98  if (lp_mpi > 0) then
99  if (all(task /= [0,1,2,3,5,15,16,28,29,61,62,63,110,120,121,125,135,136, &
100  162,170,180,185,200,201,202,205,208,209,240,241,270,271,300,320,330,331, &
101  350,351,352,371,372,373,380,390,420,421,440,460,461,462,463,471,478,600, &
102  601,610,620,630,640,680,700,701,720,725])) then
103  write(*,'("Info(elk): MPI process ",I6," idle for task ",I6)') lp_mpi,task
104  goto 20
105  end if
106  end if
107 ! write task to VARIABLES.OUT
108  if (wrtvars) call writevars('task',iv=task)
109  select case(task)
110  case(0,1)
111  call gndstate
112  case(2,3)
113  call geomopt
114  case(5)
115  call hartfock
116  case(10)
117  call writedos
118  case(14)
119  call writesf
120  case(15,16)
121  call writelsj
122  case(20,21,22,23,24)
123  call bandstr
124  case(25)
125  call effmass
126  case(28,29)
127  call mae
128  case(31,32,33)
129  call rhoplot
130  case(41,42,43)
131  call potplot
132  case(51,52,53)
133  call elfplot
134  case(61,62,63,162)
135  call wfplot
136  case(65)
137  call wfcrplot
138  case(68)
139  call rdstatus
140  case(71,72,73,81,82,83,141,142,143,151,152,153)
141  call vecplot
142  case(91,92,93)
143  call dbxcplot
144  case(100,101,103,104)
145  call fermisurf
146  case(102)
147  call fermisurfbxsf
148  case(105)
149  call nesting
150  case(110)
151  call mossbauer
152  case(115)
153  call writeefg
154  case(120)
155  call writepmat
156  case(121)
157  call dielectric
158  case(122)
159  call moke
160  case(125)
161  call nonlinopt
162  case(130)
163  call writeexpmat
164  case(135)
165  call writewfpw
166  case(140)
167  call elnes
168  case(150)
169  call writeevsp
170  case(160)
171  call torque
172  case(170)
173  call writeemd
174  case(171,172,173)
175  call emdplot
176  case(180)
177  call writeepsinv
178  case(185)
179  call writehmlbse
180  case(186)
181  call writeevbse
182  case(187)
183  call dielectric_bse
184  case(190)
185  call geomplot
186  case(195)
187  call sfacrho
188  case(196)
189  call sfacmag
190  case(200,201,202)
191  call phononsc
192  case(205)
193  call phonon
194  case(208,209)
195  call bornechg
196  case(210)
197  call phdos
198  case(220)
199  call phdisp
200  case(230)
201  call writephn
202  case(240,241)
203  call ephcouple
204  case(245)
205  call phlwidth
206  case(250)
207  call alpha2f
208  case(260)
209  call eliashberg
210  case(270,271)
211  call gndsteph
212  case(280)
213  call ephdos
214  case(285)
215  call aceplot
216  case(300)
217  call rdmft
218  case(320)
219  call tddftlr
220  case(330,331)
221  call tddftsplr
222  case(341,342,343)
223  call wxcplot
224  case(350,351,352)
225  call spiralsc
226  case(371,372,373)
227  call jprplot
228  case(380)
229  call piezoelt
230  case(390)
231  call magnetoelt
232  case(400)
233  call writetm
234  case(420,421)
235  call moldyn
236  case(430)
237  call writestrain
238  case(440)
239  call writestress
240  case(450)
241  call genafieldt
242  case(455)
243  call writeafpdt
244  case(456)
245  call writeefieldw
246  case(460,461,462,463)
247  call tddft
248  case(471)
249  call rhosplot
250  case(478)
251  call bornecdyn
252  case(480,481)
253  call dielectric_tdrt
254  case(500)
255  call testcheck
256  case(550)
257  call writew90
258  case(600,601)
259  call gwsefm
260  case(610)
261  call gwspecf
262  case(620)
263  call gwbandstr
264  case(630)
265  call writegwefm
266  case(640)
267  call gwdmat
268  case(680)
269  call tdhfc
270  case(700,701)
271  call gndstulr
272  case(710)
273  call writedosu
274  case(720,725)
275  call bandstrulr
276  case(731,732,733)
277  call rhouplot
278  case(741,742,743)
279  call potuplot
280  case(771,772,773)
281  call maguplot
282  case default
283  write(*,*)
284  write(*,'("Error(elk): task not defined : ",I8)') task
285  write(*,*)
286  stop
287  end select
288 ! reset the OpenMP thread variables
289  call omp_reset
290 20 continue
291 ! synchronise MPI processes
292  call mpi_barrier(mpicom,ierror)
293 ! check if a restart is requested
294  call checkrst
295  if (trestart) then
296  if (mp_mpi) then
297  write(*,*)
298  write(*,'("Restarting Elk")')
299  end if
300  goto 10
301  end if
302 end do
303 if (mp_mpi) then
304  call delfile('RUNNING')
305  call writebox(6,"Elk code stopped")
306 end if
307 ! terminate MPI execution environment
308 call mpi_finalize(ierror)
309 end program
310 
311 !BOI
312 ! !TITLE: {\huge{\sc The Elk Code Manual}}\\ \Large{\sc Version 10.7.8}\\ \vskip 20pt \includegraphics[height=1cm]{elk_silhouette.pdf}
313 ! !AUTHORS: {\sc J. K. Dewhurst, S. Sharma} \\ {\sc L. Nordstr\"{o}m, F. Cricchio, O. Gr\aa n\"{a}s} \\ {\sc E. K. U. Gross}
314 ! !AFFILIATION:
315 ! !INTRODUCTION: Introduction
316 ! Welcome to the Elk Code! Elk is an all-electron full-potential linearised
317 ! augmented-plane-wave (FP-LAPW) code for determining the properties of
318 ! crystalline solids. It was developed originally at the
319 ! Karl-Franzens-Universit\"{a}t Graz as part of the EXCITING EU Research and
320 ! Training Network project\footnote{EXCITING code developed under the Research
321 ! and Training Network EXCITING funded by the EU, contract No.
322 ! HPRN-CT-2002-00317}. The guiding philosophy during the implementation of the
323 ! code was to keep it as simple as possible for both users and developers
324 ! without compromising on its capabilities. All the routines are released
325 ! under either the GNU General Public License (GPL) or the GNU Lesser General
326 ! Public License (LGPL) in the hope that they may inspire other scientists to
327 ! implement new developments in the field of density functional theory and
328 ! beyond.
329 !
330 ! \section{Acknowledgments}
331 ! Lots of people contributed to the Elk code with ideas, checking and testing,
332 ! writing code or documentation and general encouragement. They include
333 ! Claudia Ambrosch-Draxl, Clas Persson, Fredrik Bultmark, Christian Brouder,
334 ! Rickard Armiento, Andrew Chizmeshya, Per Anderson, Igor Nekrasov, Sushil
335 ! Auluck, Frank Wagner, Fateh Kalarasse, J\"{u}rgen Spitaler, Stefano
336 ! Pittalis, Nektarios Lathiotakis, Tobias Burnus, Stephan Sagmeister,
337 ! Christian Meisenbichler, S\'{e}bastien Leb\`{e}gue, Yigang Zhang, Fritz
338 ! K\"{o}rmann, Alexey Baranov, Anton Kozhevnikov, Shigeru Suehara, Frank
339 ! Essenberger, Antonio Sanna, Tyrel McQueen, Tim Baldsiefen, Marty Blaber,
340 ! Anton Filanovich, Torbj\"{o}rn Bj\"{o}rkman, Martin Stankovski, Jerzy
341 ! Goraus, Markus Meinert, Daniel Rohr, Vladimir Nazarov, Kevin Krieger, Pink
342 ! Floyd, Arkardy Davydov, Florian Eich, Aldo Romero Castro, Koichi Kitahara,
343 ! James Glasbrenner, Konrad Bussmann, Igor Mazin, Matthieu Verstraete, David
344 ! Ernsting, Stephen Dugdale, Peter Elliott, Marcin Dulak, Jos\'{e} A. Flores
345 ! Livas, Stefaan Cottenier, Yasushi Shinohara, Michael Fechner, Yaroslav
346 ! Kvashnin, Tristan M\"uller, Arsenii Gerasimov, Manh Duc Le, Jon Lafuente
347 ! Bartolom\'{e}, Ren\'{e} Wirnata, Jagdish Kumar, Andrew Shyichuk, Nisha
348 ! Singh, Pietro Bonfa, Ronald Cohen, Alyn James, Chung-Yu Wang, Leon Kerber,
349 ! Yunfan Liang, Xavier Gonze, Mike Bruckhoff, Eddie Harris-Lee, Andreas
350 ! Fischer, Wenhan Chen and Jyoti Krishna. Special mention of David Singh's
351 ! very useful book on the LAPW method\footnote{D. J. Singh, {\it Planewaves,
352 ! Pseudopotentials and the LAPW Method} (Kluwer Academic Publishers, Boston,
353 ! 1994).} must also be made. Finally we would like to acknowledge the generous
354 ! support of Karl-Franzens-Universit\"{a}t Graz, the EU Marie-Curie Research
355 ! Training Networks initiative, the Max Born Institute and the Max Planck
356 ! Society.
357 !
358 ! \vspace{24pt}
359 ! Kay Dewhurst, Sangeeta Sharma \\
360 ! Lars Nordstr\"{o}m, Francesco Cricchio, Oscar Gr\aa n\"{a}s \\
361 ! Hardy Gross
362 !
363 ! \vspace{12pt}
364 ! Halle, Berlin, Uppsala and Jerusalem
365 ! \newpage
366 !
367 ! \section{Units}
368 ! Unless explicitly stated otherwise, Elk uses atomic units. In this system
369 ! $\hbar=1$, the electron mass $m=1$, the Bohr radius $a_0=1$ and the electron
370 ! charge $e=1$ (note that the electron charge is positive, so that the atomic
371 ! numbers $Z$ are negative). Thus the atomic unit of length is
372 ! 0.529177210903(80) \AA, and the atomic unit of energy is the Hartree which
373 ! equals 27.211386245988(53) eV. The unit of the external magnetic fields is
374 ! defined such that one unit of magnetic field in {\tt elk.in} equals
375 ! 1715.255541 Tesla.
376 !
377 ! \section{Compiling and running Elk}
378 ! \subsection{Compiling the code}
379 ! Unpack the code from the archive file. Edit the file {\tt make.inc} in the
380 ! {\tt elk} directory and adjust the compiler options for your computer
381 ! system. Use of machine-optimised BLAS/LAPACK and FFT libraries will result
382 ! in significant increase in performance. Following this, run
383 ! \begin{verbatim}
384 ! make
385 ! \end{verbatim}
386 ! This will hopefully compile the entire code and all the libraries into one
387 ! executable, {\tt elk}, located in the {\tt elk/src} directory. It will also
388 ! compile two useful auxiliary programs, namely {\tt spacegroup} for producing
389 ! crystal geometries from spacegroup data and {\tt eos} for fitting equations
390 ! of state to energy-volume data. If you want to compile everything all over
391 ! again, then run {\tt make clean} from the {\tt elk} directory, followed by
392 ! {\tt make}.
393 ! \subsubsection{Parallelism in Elk}
394 ! Three forms of parallelism are implemented in Elk, and all can be used in
395 ! combination with each other, with efficiency depending on the particular
396 ! task, crystal structure and computer system. You may need to contact your
397 ! system administrator for assistance with running Elk in parallel.
398 ! \begin{enumerate}
399 ! \item
400 ! OpenMP works for symmetric multiprocessors, i.e. computers that have many
401 ! cores with the same unified memory accessible to each. It is enabled by
402 ! setting the appropriate command-line options (e.g. {\tt -qopenmp} for the
403 ! Intel compiler) before compiling, and also at runtime by the environment
404 ! variable
405 ! \begin{verbatim}
406 ! export OMP_NUM_THREADS=n
407 ! \end{verbatim}
408 ! where n is the number of cores available on a particular node. The same can
409 ! be accomplished in {\tt elk.in} with
410 ! \begin{verbatim}
411 ! maxthd
412 ! n
413 ! \end{verbatim}
414 ! In addition, some vendor-supplied BLAS/LAPACK libraries use OpenMP
415 ! internally. The maximum number of threads used for LAPACK operations by
416 ! Intel's MKL can be set with
417 ! \begin{verbatim}
418 ! maxthdmkl
419 ! n
420 ! \end{verbatim}
421 ! \item
422 ! The message passing interface (MPI) is particularly suitable for running
423 ! Elk across multiple nodes of a cluster, with scaling to hundreds of
424 ! processors possible. To enable MPI, comment out the lines indicated in
425 ! {\tt elk/make.inc}. Then run {\tt make clean} followed by {\tt make}. If
426 ! $y$ is the number of nodes and $x$ is the number of cores per node, then at
427 ! runtime envoke
428 ! \begin{verbatim}
429 ! mpirun -np z ./elk
430 ! \end{verbatim}
431 ! where $z=x y$ is the total number of cores available on the machine.
432 ! Highest efficiency is obtained by using hybrid parallelism with OpenMP on
433 ! each node and MPI across nodes. This can be done by compiling the code
434 ! using the MPI Fortran compiler in combination with the OpenMP command-line
435 ! option. At runtime set {\tt export OMP\_NUM\_THREADS=x} and start the MPI
436 ! run with {\em one process per node} as follows
437 ! \begin{verbatim}
438 ! mpirun -pernode -np y ./elk
439 ! \end{verbatim}
440 ! The number of MPI processes is reported in the file {\tt INFO.OUT} which
441 ! serves as a check that MPI is running correctly. Note that version 2 of the
442 ! MPI libraries is required to run Elk.
443 ! \item
444 ! Phonon calculations use a simple form of parallelism by just examining the
445 ! run directory for dynamical matrix files. These files are of the form
446 ! \begin{verbatim}
447 ! DYN_Qqqqq_qqqq_qqqq_Sss_Aaa_Pp.OUT
448 ! \end{verbatim}
449 ! and contain a single row of a particular dynamical matrix. Elk simply finds
450 ! which {\tt DYN} files do not exist, chooses one and runs it. This way many
451 ! independent runs of Elk can be started in the same directory on a networked
452 ! file system (NFS), and will run until all the dynamical matrices files are
453 ! completed. Should a particular run crash, then delete the associated empty
454 ! {\tt DYN} file and rerun Elk.
455 ! \end{enumerate}
456 !
457 ! \subsection{Memory requirements}
458 ! Elk is a memory-bound code and runs best on processors with large caches and
459 ! a large number of memory channels per core. Some tasks in Elk require a
460 ! considerable amount of memory which can exceed the physical memory of the
461 ! computer. In such cases, the number of threads at the first nesting level
462 ! can be reduced with (for example)
463 ! \begin{verbatim}
464 ! maxthd1
465 ! -4
466 ! \end{verbatim}
467 ! which restricts the number of threads at the first nesting level to
468 ! maxthd/4. Deeper nesting levels, which generally require less memory, will
469 ! still utilise the full compliment of available threads.
470 ! \subsubsection{Stack space}
471 ! The latest versions of Elk use stack space aggressively. This is because
472 ! accessing variables is faster on the stack than on the heap. This can,
473 ! however, result in the code crashing as threads run out of their stack
474 ! space. To avoid this, increase the stack size for each OpenMP thread with
475 ! (for example)
476 ! \begin{verbatim}
477 ! export OMP_STACKSIZE=512M
478 ! \end{verbatim}
479 ! before running the code.
480 !
481 ! \subsection{Linking with the Libxc functional library}
482 ! Libxc is a library of exchange-correlation functionals. Elk can use the
483 ! complete set of LDA and GGA functionals available in Libxc as well as some
484 ! meta-GGAs. In order to enable this, first download and compile Libxc. This
485 ! should have produced the files {\tt libxc.a} and {\tt libxcf03.a} in the
486 ! Libxc directory {\tt src/.libs}. Copy these to the {\tt elk/src} directory
487 ! and then uncomment the lines indicated for Libxc in {\tt elk/make.inc}. Once
488 ! this is done, run {\tt make clean} followed by {\tt make}. To select a
489 ! particular functional of Libxc, use the block
490 ! \begin{verbatim}
491 ! xctype
492 ! 100 nx nc
493 ! \end{verbatim}
494 ! where {\tt nx} and {\tt nc} are, respectively, the numbers of the exchange
495 ! and correlation functionals in the Libxc library. See the file
496 ! {\tt elk/src/libxcf03.f90} for a list of the functionals and their
497 ! associated numbers.
498 !
499 ! \subsection{Running the code}
500 ! As a rule, all input files for the code are in lower case and end with the
501 ! extension {\tt .in}. All output files are uppercase and have the extension
502 ! {\tt .OUT}. For most cases, the user will only need to modify the file
503 ! {\tt elk.in}. In this file input parameters are arranged in blocks.
504 ! Each block consists of a block name on one line and the block variables on
505 ! subsequent lines. Almost all blocks are optional: the code uses reasonable
506 ! default values in cases where they are absent. Blocks can appear in any
507 ! order, if a block is repeated then the second instance is used. Comment
508 ! lines can be included in the input file and begin with the {\tt !}
509 ! character.
510 !
511 ! \subsubsection{Species files}
512 ! The only other input files are those describing the atomic species which go
513 ! into the crystal. These files are found in the {\tt species} directory and
514 ! are named with the element symbol and the extension {\tt .in}, for example
515 ! {\tt Sb.in}. They contain parameters like the atomic charge, mass,
516 ! muffin-tin radius, occupied atomic states and the type of linearisation
517 ! required. Here as an example is the copper species file {\tt Cu.in}:
518 ! \begin{verbatim}
519 ! 'Cu' : spsymb
520 ! 'copper' : spname
521 ! -29.0000 : spzn
522 ! 115837.2716 : spmass
523 ! 0.371391E-06 2.0000 34.8965 500 : rminsp, rmt, rmaxsp, nrmt
524 ! 10 : nstsp
525 ! 1 0 1 2.00000 T : nsp, lsp, ksp, occsp, spcore
526 ! 2 0 1 2.00000 T
527 ! 2 1 1 2.00000 T
528 ! 2 1 2 4.00000 T
529 ! 3 0 1 2.00000 T
530 ! 3 1 1 2.00000 F
531 ! 3 1 2 4.00000 F
532 ! 3 2 2 4.00000 F
533 ! 3 2 3 6.00000 F
534 ! 4 0 1 1.00000 F
535 ! 1 : apword
536 ! 0.1500 0 F : apwe0, apwdm, apwve
537 ! 1 : nlx
538 ! 2 2 : lx, apword
539 ! 0.1500 0 T : apwe0, apwdm, apwve
540 ! 0.1500 1 T
541 ! 4 : nlorb
542 ! 0 2 : lorbl, lorbord
543 ! 0.1500 0 F : lorbe0, lorbdm, lorbve
544 ! 0.1500 1 F
545 ! 1 2
546 ! 0.1500 0 F
547 ! 0.1500 1 F
548 ! 2 2
549 ! 0.1500 0 F
550 ! 0.1500 1 F
551 ! 1 3
552 ! 0.1500 0 F
553 ! 0.1500 1 F
554 ! -2.8652 0 T
555 ! \end{verbatim}
556 ! The input parameters are defined as follows:
557 ! \vskip 6pt
558 ! {\tt spsymb} \\
559 ! The symbol of the element.
560 ! \vskip 6pt
561 ! {\tt spname} \\
562 ! The name of the element.
563 ! \vskip 6pt
564 ! {\tt spzn} \\
565 ! Nuclear charge: should be negative since the electron charge is taken to be
566 ! postive in the code; it can also be fractional for purposes of doping.
567 ! \vskip 6pt
568 ! {\tt spmass} \\
569 ! Nuclear mass in atomic units.
570 ! \vskip 6pt
571 ! {\tt rminsp}, {\tt rmt}, {\tt rmaxsp}, {\tt nrmt} \\
572 ! Respectively, the minimum radius on logarithmic radial mesh; muffin-tin
573 ! radius; effective infinity for atomic radial mesh; and number of radial mesh
574 ! points to muffin-tin radius.
575 ! \vskip 6pt
576 ! {\tt nstsp} \\
577 ! Number of atomic states.
578 ! \vskip 6pt
579 ! {\tt nsp}, {\tt lsp}, {\tt ksp}, {\tt occsp}, {\tt spcore} \\
580 ! Respectively, the principal quantum number of the radial Dirac equation;
581 ! quantum number $l$; quantum number $k$ ($l$ or $l+1$); occupancy of atomic
582 ! state (can be fractional); {\tt .T.} if state is in the core and therefore
583 ! treated with the Dirac equation in the spherical part of the muffin-tin
584 ! Kohn-Sham potential.
585 ! \vskip 6pt
586 ! {\tt apword} \\
587 ! Default APW function order, i.e. the number of radial functions and
588 ! therefore the order of the radial derivative matching at the muffin-tin
589 ! surface.
590 ! \vskip 6pt
591 ! {\tt apwe0}, {\tt apwdm}, {\tt apwve} \\
592 ! Respectively, the default APW linearisation energy; the order of the energy
593 ! derivative of the APW radial function $\partial^m u(r)/\partial E^m$; and
594 ! {\tt .T.} if the linearisation energy is allowed to vary.
595 ! \vskip 6pt
596 ! {\tt nlx} \\
597 ! The number of exceptions to the default APW configuration. These should be
598 ! listed on subsequent lines for particular angular momenta. In this example,
599 ! the fixed energy APW with angular momentum $d$ ({\tt lx} $=2$) is replaced
600 ! with a LAPW, which has variable linearisation energy.
601 ! \vskip 6pt
602 ! {\tt nlorb} \\
603 ! Number of local-orbitals.
604 ! \vskip 6pt
605 ! {\tt lorbl}, {\tt lorbord} \\
606 ! Respectively, the angular momentum $l$ of the local-orbital; and the order
607 ! of the radial derivative which goes to zero at the muffin-tin surface.
608 ! \vskip 6pt
609 ! {\tt lorbe0}, {\tt lorbdm}, {\tt lorbve} \\
610 ! Respectively, the default local-orbital linearisation energy; the order of
611 ! the energy derivative of the local-orbital radial function; and {\tt .T.} if
612 ! the linearisation energy is allowed to vary.
613 !
614 ! \subsubsection{Examples}
615 ! The best way to learn to use Elk is to run the examples included with the
616 ! package. These can be found in the {\tt examples} directory and use many of
617 ! the code's capabilities. The following section which describes all the input
618 ! parameters will be of invaluable assistance.
619 !
620 ! \section{Input blocks}
621 ! This section lists all the input blocks available. It is arranged with the
622 ! name of the block followed by a table which lists each parameter name, what
623 ! the parameter does, its type and default value. A horizontal line in the
624 ! table indicates a new line in {\tt elk.in}. Below the table is a brief
625 ! overview of the block's function.
626 !
627 ! \block{atoms}{
628 ! {\tt nspecies} & number of species & integer & 0 \\
629 ! \hline
630 ! {\tt spfname(i)} & species filename for species $i$ & string & - \\
631 ! \hline
632 ! {\tt natoms(i)} & number of atoms for species $i$ & integer & - \\
633 ! \hline
634 ! {\tt atposl(j,i)} & atomic position in lattice coordinates for atom $j$
635 ! & real(3) & - \\
636 ! {\tt bfcmt(j,i)} & muffin-tin external magnetic field in Cartesian
637 ! coordinates for atom $j$ & real(3) & -}
638 ! Defines the atomic species as well as their positions in the unit cell and
639 ! the external magnetic field applied throughout the muffin-tin. These fields
640 ! are used to break spin symmetry and should be considered infinitesimal as
641 ! they do not contribute directly to the total energy. Collinear calculations
642 ! are more efficient if the field is applied in the $z$-direction. One could,
643 ! for example, set up an antiferromagnetic crystal by pointing the field on
644 ! one atom in the positive $z$-direction and in the opposite direction on
645 ! another atom. If {\tt molecule} is {\tt .true.} then the atomic positions
646 ! are assumed to be in Cartesian coordinates. See also {\tt sppath},
647 ! {\tt bfieldc} and {\tt molecule}.
648 !
649 ! \block{autokpt}{
650 ! {\tt autokpt} & {\tt .true.} if the $k$-point set is to be determined
651 ! automatically & logical & {\tt .false.}}
652 ! See {\tt radkpt} for details.
653 !
654 ! \block{autolinengy}{
655 ! {\tt autolinengy} & {\tt .true.} if the fixed linearisation energies are
656 ! to be determined automatically & logical & {\tt .false.}}
657 ! See {\tt dlefe} for details.
658 !
659 ! \block{autoswidth}{
660 ! {\tt autoswidth} & {\tt .true.} if the smearing parameter {\tt swidth}
661 ! should be determined automatically & logical & {\tt .false.}}
662 ! Calculates the smearing width from the $k$-point density, $V_{\rm BZ}/n_k$;
663 ! the valence band width, $W$; and an effective mass parameter, $m^{*}$;
664 ! according to
665 ! $$ \sigma=\frac{\sqrt{2W}}{m^{*}}\left(\frac{3}{4\pi}
666 ! \frac{V_{\rm BZ}}{n_k}\right)^{1/3}. $$
667 ! The variable {\tt mstar} then replaces {\tt swidth} as the control parameter
668 ! of the smearing width. A large value of $m^{*}$ gives a narrower smearing
669 ! function. Since {\tt swidth} is adjusted according to the fineness of the
670 ! ${\bf k}$-mesh, the smearing parameter can then be eliminated. It is not
671 ! recommended that {\tt autoswidth} be used in conjunction with the
672 ! Fermi-Dirac smearing function, since the electronic temperature will then be
673 ! a function of the $k$-point mesh. See T. Bj\"orkman and O. Gr\aa n\"as,
674 ! {\it Int. J. Quant. Chem.} DOI: 10.1002/qua.22476 (2010) for details. See
675 ! also {\tt stype} and {\tt swidth}.
676 !
677 ! \block{avec}{
678 ! {\tt avec(1)} & first lattice vector & real(3) & $(1.0,0.0,0.0)$ \\
679 ! \hline
680 ! {\tt avec(2)} & second lattice vector & real(3) & $(0.0,1.0,0.0)$ \\
681 ! \hline
682 ! {\tt avec(3)} & third lattice vector & real(3) & $(0.0,0.0,1.0)$}
683 ! Lattice vectors of the crystal in atomic units (Bohr).
684 !
685 ! \block{avecref}{
686 ! {\tt avecref(1)} & first reference lattice vector, etc. & real(3) &
687 ! $(0.0,0.0,0.0)$}
688 ! Reference lattice vectors for calculating the ${\bf G}$-vector grid and
689 ! derived quantities. If any of these elements are non-zero then the code
690 ! computes the corresponding reciprocal lattice vectors and set of
691 ! ${\bf G}$-vectors. These are then transformed to be identical to those
692 ! calculated with {\tt avec}. The purpose of this is to enable accurate
693 ! energy-volume curves, where the number of grid points should remain fixed
694 ! for all volumes.
695 !
696 ! \block{beta0}{
697 ! {\tt beta0} & adaptive mixing parameter & real & $0.05$}
698 ! This determines how much of the potential from the previous self-consistent
699 ! loop is mixed with the potential from the current loop. It should be made
700 ! smaller if the calculation is unstable. See {\tt betamax} and also the
701 ! routine {\tt mixadapt}.
702 !
703 ! \block{betamax}{
704 ! {\tt betamax} & maximum adaptive mixing parameter & real & $0.5$}
705 ! Maximum allowed mixing parameter used in routine {\tt mixadapt}.
706 !
707 ! \block{bfdmag}{
708 ! {\tt bfdmag} & {\tt .true.} if the external ${\bf B}$-field diamagnetic
709 ! coupling term should be included & logical & {\tt .false.}}
710 ! This causes the diamagnetic coupling term
711 ! $$ H_{\rm dia}=\frac{B^2 r^2}{8c^2}\big(1-(\hat{\bf B}\cdot
712 ! \hat{\bf r})^2\big) $$
713 ! to be included in the first-variational Hamiltonian. Note that because this
714 ! is a scalar potential, spin-polarisation does not have to be enabled for
715 ! this term to have an effect.
716 !
717 ! \block{bfieldc}{
718 ! {\tt bfieldc} & global external magnetic field in Cartesian coordinates &
719 ! real(3) & $(0.0,0.0,0.0)$}
720 ! This is a constant magnetic field applied throughout the entire unit cell
721 ! and enters the second-variational Hamiltonian as
722 ! $$ \frac{g_e}{4c}\,\vec{\sigma}\cdot{\bf B}, $$
723 ! where $g_e$ is the electron $g$-factor. This field is normally used to break
724 ! spin symmetry for spin-polarised calculations and considered to be
725 ! infinitesimal with no direct contribution to the total energy. In cases
726 ! where the magnetic field is finite (for example when computing magnetic
727 ! response) the external ${\bf B}$-field energy reported in {\tt INFO.OUT}
728 ! should be added to the total by hand. This field is applied throughout the
729 ! entire unit cell. To apply magnetic fields in particular muffin-tins use the
730 ! {\tt bfcmt} vectors in the {\tt atoms} block. Collinear calculations are
731 ! more efficient if the field is applied in the $z$-direction.
732 !
733 ! \block{bfieldcu}{
734 ! {\tt bfieldcu} & global external magnetic field in Cartesian coordinates &
735 ! real(3) & $(0.0,0.0,0.0)$}
736 ! As with {\tt bfieldc} but applied only in the ultracell during an ultra
737 ! long-range calculation.
738 !
739 ! \block{bforb}{
740 ! {\tt bforb} & {\tt .true.} if the external ${\bf B}$-field-orbit coupling
741 ! term should be included & logical & {\tt .false.}}
742 ! This causes the term corresponding to the coupling between an external
743 ! magnetic field and the orbit of an electron (orbital paramagnetism)
744 ! $$ \hat{H}_{{\bf B}{\rm o}}=\frac{1}{2c}{\bf B}\cdot\hat{\bf L} $$
745 ! to be added to the second-variational Hamiltonian. Spin-polarisation is
746 ! automatically enabled.
747 !
748 ! \block{broydpm}{
749 ! {\tt broydpm} & Broyden mixing parameters $\alpha$ and $w_0$ & real &
750 ! $(0.4,0.15)$}
751 ! See {\tt mixtype} and {\tt mixsdb}.
752 !
753 ! \block{c\_tb09}{
754 ! {\tt c\_tb09} & Tran-Blaha constant $c$ & real & -}
755 ! Sets the constant $c$ in the Tran-Blaha '09 functional. Normally this is
756 ! calculated from the density, but there may be situations where this needs to
757 ! be adjusted by hand. See {\it Phys. Rev. Lett.} {\bf 102}, 226401 (2009).
758 !
759 ! \block{chgexs}{
760 ! {\tt chgexs} & excess electronic charge & real & $0.0$}
761 ! This controls the amount of charge in the unit cell beyond that required to
762 ! maintain neutrality. It can be set positive or negative depending on whether
763 ! electron or hole doping is required.
764 !
765 ! \block{cmagz}{
766 ! {\tt cmagz} & .true. if $z$-axis collinear magnetism is to be enforced &
767 ! logical & {\tt .false.}}
768 ! This variable can be set to .true. in cases where the magnetism is
769 ! predominantly collinear in the $z$-direction, for example a ferromagnet with
770 ! spin-orbit coupling. This will make the calculation considerably faster at
771 ! the slight expense of precision.
772 !
773 ! \block{deltaem}{
774 ! {\tt deltaem} & the size of the ${\bf k}$-vector displacement used when
775 ! calculating numerical derivatives for the effective mass tensor & real &
776 ! $0.025$}
777 ! See {\tt ndspem} and {\tt vklem}.
778 !
779 ! \block{deltaph}{
780 ! {\tt deltaph} & size of the atomic displacement used for calculating
781 ! dynamical matrices & real & $0.01$}
782 ! Phonon calculations are performed by constructing a supercell corresponding
783 ! to a particular ${\bf q}$-vector and making a small periodic displacement of
784 ! the atoms. The magnitude of this displacement is given by {\tt deltaph}.
785 ! This should not be made too large, as anharmonic terms could then become
786 ! significant, neither should it be too small as this can introduce numerical
787 ! error.
788 !
789 ! \block{deltast}{
790 ! {\tt deltast} & size of the change in lattice vectors used for calculating
791 ! the stress tensor & real & $0.005$}
792 ! The stress tensor is computed by changing the lattice vector matrix $A$ by
793 ! $$ A\rightarrow (1+\delta t\,e_k)A, $$
794 ! where $\delta t$ is an infinitesimal equal in practice to {\tt deltast} and
795 ! $e_k$ is the $k^{\rm th}$ strain tensor. Numerical finite differences are
796 ! used to compute the stress tensor as the derivative of the total energy
797 ! $dE_k/dt$.
798 !
799 ! \block{dft+u}{
800 ! {\tt dftu} & type of DFT+$U$ calculation & integer & 0 \\
801 ! {\tt inpdftu} & type of input for DFT+U calculation & integer & 1 \\
802 ! \hline
803 ! {\tt is} & species number & integer & - \\
804 ! {\tt l} & angular momentum value & integer & -1 \\
805 ! {\tt u} & the desired $U$ value & real & $0.0$ \\
806 ! {\tt j} & the desired $J$ value & real & $0.0$}
807 ! This block contains the parameters required for an DFT+$U$ calculation, with
808 ! the list of parameters for each species terminated with a blank line. The
809 ! type of double counting required is set with the parameter {\tt dftu}.
810 ! Currently implemented are:
811 ! \vskip 6pt
812 ! \begin{tabularx}{\textwidth}[h]{lX}
813 ! 0 & No DFT+$U$ calculation \\
814 ! 1 & Fully localised limit (FLL) \\
815 ! 2 & Around mean field (AFM) \\
816 ! 3 & An interpolation between FLL and AFM
817 ! \end{tabularx}
818 ! \vskip 6pt
819 ! The type of input parameters is set with the parameter {\tt inpdftu}.
820 ! The current possibilities are:
821 ! \vskip 6pt
822 ! \begin{tabularx}{\textwidth}[h]{lX}
823 ! 1 & U and J \\
824 ! 2 & Slater parameters \\
825 ! 3 & Racah parameters \\
826 ! 4 & Yukawa screening length \\
827 ! 5 & U and determination of corresponding Yukawa screening length
828 ! \end{tabularx}
829 ! \vskip 6pt
830 ! See (amongst others) {\it Phys. Rev. B} {\bf 67}, 153106 (2003),
831 ! {\it Phys. Rev. B} {\bf 52}, R5467 (1995), {\it Phys. Rev. B} {\bf 60},
832 ! 10763 (1999), and {\it Phys. Rev. B} {\bf 80}, 035121 (2009).
833 !
834 ! \block{dlefe}{
835 ! {\tt dlefe} & difference between the fixed linearisation energy and the
836 ! Fermi energy & real & $-0.1$}
837 ! When {\tt autolinengy} is {\tt .true.} then the fixed linearisation energies
838 ! are set to the Fermi energy plus {\tt dlefe}.
839 !
840 ! \block{dncgga}{
841 ! {\tt dncgga} & small constant used to stabilise non-collinear GGA &
842 ! real & $1\times 10^{-8}$}
843 ! This small constant, $d$, is required in order to remove the infinite
844 ! gradients obtained when using `Kubler's trick' in conjunction with GGA and
845 ! non-collinear magnetism. It is applied by calculating the up and down
846 ! densities as
847 ! $$ \rho^{\uparrow}({\bf r})=\rho({\bf r})+\widetilde{m}({\bf r})
848 ! \qquad \rho^{\downarrow}({\bf r})=\rho({\bf r})-\widetilde{m}({\bf r}), $$
849 ! where $\widetilde{m}({\bf r})=\sqrt{{\bf m}^2({\bf r})+d}$,
850 ! and should be taken as the smallest value for which the exchange-correlation
851 ! magnetic field ${\bf B}_{\rm xc}$ is smooth.
852 !
853 ! \block{dosmsum}{
854 ! {\tt dosmsum} & {\tt .true.} if the partial DOS is to be summed over $m$ &
855 ! logical & {\tt .false.}}
856 ! By default, the partial density of states is resolved over $(l,m)$ quantum
857 ! numbers. If {\tt dosmsum} is set to {\tt .true.} then the partial DOS is
858 ! summed over $m$, and thus depends only on $l$.
859 !
860 ! \block{dosssum}{
861 ! {\tt dosssum} & {\tt .true.} if the partial DOS is to be summed over spin &
862 ! logical & {\tt .false.}}
863 ! By default, the partial density of states for spin-polarised systems is spin
864 ! resolved.
865 !
866 ! \block{dtimes}{
867 ! {\tt dtimes} & time step used in time evolution run & real & $0.1$}
868 ! See also {\tt tstime}.
869 !
870 ! \block{epsband}{
871 ! {\tt epsband} & convergence tolerance for determining band energies & real &
872 ! $1\times 10^{-12}$}
873 ! APW and local-orbital linearisation energies are determined from the band
874 ! energies. This is done by first searching upwards in energy until the radial
875 ! wavefunction at the muffin-tin radius is zero. This is the energy at the top
876 ! of the band, denoted $E_{\rm t}$. A downward search is now performed from
877 ! $E_{\rm t}$ until the slope of the radial wavefunction at the muffin-tin
878 ! radius is zero. This energy, $E_{\rm b}$, is at the bottom of the band. The
879 ! band energy is taken as $(E_{\rm t}+E_{\rm b})/2$. If either $E_{\rm t}$ or
880 ! $E_{\rm b}$ is not found, then the band energy is set to the default value.
881 !
882 ! \block{epschg}{
883 ! {\tt epschg} & maximum allowed error in the calculated total charge beyond
884 ! which a warning message will be issued & real & $1\times 10^{-3}$}
885 !
886 ! \block{epsengy}{
887 ! {\tt epsengy} & convergence criterion for the total energy & real &
888 ! $1\times 10^{-4}$}
889 ! See {\tt epspot}.
890 !
891 ! \block{epsforce}{
892 ! {\tt epsforce} & convergence tolerance for the forces during a geometry
893 ! optimisation run & real & $5\times 10^{-3}$}
894 ! If the mean absolute value of the atomic forces is less than {\tt epsforce}
895 ! then the geometry optimisation run is ended. See also {\tt tasks} and
896 ! {\tt latvopt}.
897 !
898 ! \block{epslat}{
899 ! {\tt epslat } & vectors with lengths less than this are considered zero &
900 ! real & $10^{-6}$}
901 ! Sets the tolerance for determining if a vector or its components are zero.
902 ! This is to account for any numerical error in real or reciprocal space
903 ! vectors.
904 !
905 ! \block{epsocc}{
906 ! {\tt epsocc} & smallest occupancy for which a state will contribute to the
907 ! density & real & $1\times 10^{-8}$}
908 !
909 ! \block{epspot}{
910 ! {\tt epspot} & convergence criterion for the Kohn-Sham potential and field &
911 ! real & $1\times 10^{-6}$}
912 ! If the RMS change in the Kohn-Sham potential and magnetic field is smaller
913 ! than {\tt epspot} and the absolute change in the total energy is less than
914 ! {\tt epsengy}, then the self-consistent loop is considered converged
915 ! and exited. For geometry optimisation runs this results in the forces being
916 ! calculated, the atomic positions updated and the loop restarted. See also
917 ! {\tt epsengy} and {\tt maxscl}.
918 !
919 ! \block{epsstress}{
920 ! {\tt epsstress} & convergence tolerance for the stress tensor during a
921 ! geometry optimisation run with lattice vector relaxation & real &
922 ! $2\times 10^{-3}$}
923 ! See also {\tt epsforce} and {\tt latvopt}.
924 !
925 ! \block{emaxelnes}{
926 ! {\tt emaxelnes} & maximum allowed initial-state eigenvalue for ELNES
927 ! calculations & real & $-1.2$}
928 !
929 ! \block{emaxrf}{
930 ! {\tt emaxrf} & energy cut-off used when calculating Kohn-Sham response
931 ! functions & real & $10^6$}
932 ! A typical Kohn-Sham response function is of the form
933 ! \begin{align*}
934 ! \chi_s({\bf r},{\bf r}',\omega)
935 ! \equiv\frac{\delta\rho({\bf r},\omega)}{\delta v_s({\bf r}',\omega)}
936 ! =\frac{1}{N_k}\sum_{i{\bf k},j{\bf k}'}(f_{i{\bf k}}-f_{j{\bf k}'})
937 ! \frac{\langle i{\bf k}|\hat{\rho}({\bf r})|j{\bf k}'\rangle
938 ! \langle j{\bf k}'|\hat{\rho}({\bf r}')|i{\bf k}\rangle}
939 ! {w+(\varepsilon_{i{\bf k}}-\varepsilon_{j{\bf k}'})+i\eta},
940 ! \end{align*}
941 ! where $\hat{\rho}$ is the density operator; $N_k$ is the number of
942 ! $k$-points; $\varepsilon_{i{\bf k}}$ and $f_{i{\bf k}}$ are the eigenvalues
943 ! and occupation numbers, respectively. The variable {\tt emaxrf} is an energy
944 ! window which limits the summation over states in the formula above so that
945 ! $|\varepsilon_{i{\bf k}}-\varepsilon_{\rm Fermi}|<{\tt emaxrf}$. Reducing
946 ! this can result in a faster calculation at the expense of accuracy.
947 !
948 ! \block{fracinr}{
949 ! {\tt fracinr} & fraction of the muffin-tin radius up to which {\tt lmaxi}
950 ! is used as the angular momentum cut-off & real & $0.01$}
951 ! If {\tt fracinr} is negative then the fraction is determined from
952 ! $f=\sqrt{({\tt lmaxi}+1)^2/({\tt lmaxo}+1)^2}$ in order to
953 ! maintain a minimum density of points throughout the muffin-tin. See
954 ! {\tt lmaxi} and {\tt lmaxo}.
955 !
956 ! \block{fsmtype}{
957 ! {\tt fsmtype} & fixed spin moment (FSM) type & integer & 0}
958 ! The magnetic moment, its direction or magnitude can be fixed both globally
959 ! as well as in each muffin-tin individually. The options are as follows:
960 ! \vskip 6pt
961 ! \begin{tabularx}{\textwidth}[h]{lX}
962 ! 0 & no FSM \\
963 ! 1 (-1) & total moment (direction) \\
964 ! 2 (-2) & individual muffin-tin moments (direction) \\
965 ! 3 (-3) & total and muffin-tin moments (direction) \\
966 ! 4 & total moment magnitude \\
967 ! 5 & individual muffin-tin moment magnitudes \\
968 ! 6 & total and muffin-tin moment magnitudes
969 ! \end{tabularx}
970 ! \vskip 6pt
971 ! See also {\tt momfix}, {\tt momfixm}, {\tt mommtfix}, {\tt mommtfixm},
972 ! {\tt taufsm} and {\tt spinpol}.
973 !
974 ! \block{ftmtype}{
975 ! {\tt ftmtype} & 1 to enable a fixed tensor moment (FTM) calculation,
976 ! 0 otherwise & integer & 0}
977 ! If {\tt ftmtype} is $-1$ then the symmetry corresponding to the tensor
978 ! moment is broken but no FTM calculation is performed. See also {\tt tm3fix}.
979 !
980 ! \block{fxclrc}{
981 ! {\tt fxclrc} & parameters for the dynamical long-range contribution (LRC) to
982 ! the TDDFT exchange-correlation kernel & real(2) & $(0.0,0.0)$}
983 ! These are the parameters $\alpha$ and $\beta$ for the kernel proposed in
984 ! {\it Phys. Rev. B} {\bf 72}, 125203 (2005), namely
985 ! $$ f_{xc}({\bf G},{\bf G}',{\bf q},\omega)=-\frac{\alpha+\beta\omega^2}{q^2}
986 ! \delta_{{\bf G},{\bf G}'}\delta_{{\bf G},{\bf 0}}. $$
987 !
988 ! \block{fxctype}{
989 ! {\tt fxctype} & integer defining the type of exchange-correlation kernel
990 ! $f_{\rm xc}$ & integer & $-1$}
991 ! The acceptable values are:
992 ! \vskip 6pt
993 ! \begin{tabularx}{\textwidth}[h]{lX}
994 ! $-1$ & $f_{\rm xc}$ defined by {\tt xctype} \\
995 ! 0,1 & RPA ($f_{\rm xc}=0$) \\
996 ! 200 & Long-range contribution (LRC) kernel, S. Botti {\it et al.},
997 ! {\it Phys. Rev. B} {\bf 72}, 125203 (2005); see {\tt fxclrc} \\
998 ! 210 & `Bootstrap' kernel, S. Sharma, J. K. Dewhurst, A. Sanna and
999 ! E. K. U. Gross, {\it Phys. Rev. Lett.} {\bf 107}, 186401 (2011) \\
1000 ! 211 & Single iteration bootstrap
1001 ! \end{tabularx}
1002 !
1003 ! \block{gmaxrf}{
1004 ! {\tt gmaxrf} & maximum length of $|{\bf G}|$ for computing response
1005 ! functions & real & $3.0$}
1006 !
1007 ! \block{gmaxvr}{
1008 ! {\tt gmaxvr} & maximum length of $|{\bf G}|$ for expanding the interstitial
1009 ! density and potential & real & $12.0$}
1010 ! This variable has a lower bound which is enforced by the code as follows:
1011 ! $$ {\rm gmaxvr}\rightarrow\max\,({\rm gmaxvr},2\times{\rm gkmax}
1012 ! +{\rm epslat}) $$
1013 ! See {\tt rgkmax}.
1014 !
1015 ! \block{hdbse}{
1016 ! {\tt hdbse} & {\tt .true.} if the direct term is to be included in the BSE
1017 ! Hamiltonian & logical & {\tt .true.}}
1018 !
1019 ! \block{highq}{
1020 ! {\tt highq} & {\tt .true.} if a high-quality parameter set should be used &
1021 ! logical & {\tt .false.}}
1022 ! Setting this to {\tt .true.} results in some default parameters being
1023 ! changed to ensure good convergence in most situations. These changes can be
1024 ! overruled by subsequent blocks in the input file. See also {\tt vhighq}.
1025 !
1026 ! \block{hmaxvr}{
1027 ! {\tt hmaxvr} & maximum length of ${\bf H}$-vectors & real & $6.0$}
1028 ! The ${\bf H}$-vectors are used for calculating X-ray and magnetic structure
1029 ! factors. They are also used in linear response phonon calculations for
1030 ! expanding the density and potential in plane waves. See also {\tt gmaxvr},
1031 ! {\tt vhmat}, {\tt reduceh}, {\tt wsfac} and {\tt hkmax}.
1032 !
1033 ! \block{hxbse}{
1034 ! {\tt hxbse} & {\tt .true.} if the exchange term is to be included in the BSE
1035 ! Hamiltonian & {\tt .true.}}
1036 !
1037 ! \block{hybrid}{
1038 ! {\tt hybrid} & {\tt .true} if a hybrid functional is to be used when running
1039 ! a Hartree-Fock calculation & logical & {\tt .false}}
1040 ! See also {\tt hybridc} and {\tt xctype}.
1041 !
1042 ! \block{hybridc}{
1043 ! {\tt hybridc} & hybrid functional mixing coefficient & real & $1.0$}
1044 !
1045 ! \block{intraband}{
1046 ! {\tt intraband} & {\tt .true.} if the intraband (Drude-like) contribution is
1047 ! to be added to the dieletric tensor & logical & {\tt .false.}}
1048 !
1049 ! \block{isgkmax}{
1050 ! {\tt isgkmax} & species for which the muffin-tin radius will be used for
1051 ! calculating {\tt gkmax} & integer & $-1$}
1052 ! The APW cut-off is determined from ${\tt gkmax}={\tt rgkmax}/R$. The
1053 ! variable {\tt isgkmax} determines which muffin-tin radius is to be used for
1054 ! $R$. These are the options:
1055 ! \vskip 6pt
1056 ! \begin{tabularx}{\textwidth}[h]{lX}
1057 ! -4 & Use the largest radius \\
1058 ! -3 & Use the smallest radius \\
1059 ! -2 & Use the fixed value $R=2.0$ \\
1060 ! -1 & Use the average of the muffin-tin radii \\
1061 ! $n\ge 1$ & Use the radius of species $n$
1062 ! \end{tabularx}
1063 !
1064 ! \block{kstlist}{
1065 ! {\tt kstlist(i)} & $i$th $k$-point and state pair & integer(2) & $(1,1)$}
1066 ! This is a user-defined list of $k$-point and state index pairs which are
1067 ! those used for plotting wavefunctions and writing ${\bf L}$, ${\bf S}$ and
1068 ! ${\bf J}$ expectation values. Only the first pair is used by the
1069 ! aforementioned tasks. The list should be terminated by a blank line.
1070 !
1071 ! \block{latvopt}{
1072 ! {\tt latvopt} & type of lattice vector optimisation to be performed during
1073 ! structural relaxation & integer & 0}
1074 ! Optimisation of the lattice vectors will be performed with ${\tt task}=2,3$
1075 ! when ${\tt latvopt}\ne 0$. When ${\tt latvopt}=1$ the lattice vector
1076 ! optimisation will be constrained only by symmetry. Optimisation over all
1077 ! symmetry-preserving strains except isotropic scaling is performed when
1078 ! ${\tt latvopt}=2$. If ${\tt latvopt}<0$ then the optimisation will be over
1079 ! strain number $|{\tt latvopt}|$. The list of symmetric strain tensors can be
1080 ! produced with ${\tt task}=430$. By default (${\tt latvopt}=0$) no lattice
1081 ! vector optimisation is performed during structural relaxation. See also
1082 ! {\tt tau0latv} and {\tt atpopt}.
1083 !
1084 ! \block{lmaxapw}{
1085 ! {\tt lmaxapw} & angular momentum cut-off for the APW functions & integer &
1086 ! $8$}
1087 !
1088 ! \block{lmaxdos}{
1089 ! {\tt lmaxdos} & angular momentum cut-off for the partial DOS plot &
1090 ! integer & $3$}
1091 !
1092 ! \block{lmaxi}{
1093 ! {\tt lmaxi} & angular momentum cut-off for the muffin-tin density and
1094 ! potential on the inner part of the muffin-tin & integer & 2}
1095 ! Close to the nucleus, the density and potential is almost spherical and
1096 ! therefore the spherical harmonic expansion can be truncated a low angular
1097 ! momentum. See also {\tt fracinr}.
1098 !
1099 ! \block{lmaxo}{
1100 ! {\tt lmaxo} & angular momentum cut-off for the muffin-tin density and
1101 ! potential & integer & 6}
1102 !
1103 ! \block{lmirep}{
1104 ! {\tt lmirep} & {\tt .true.} if the $Y_{lm}$ basis is to be transformed
1105 ! into the basis of irreducible representations of the site symmetries for
1106 ! DOS plotting & logical & {\tt .true.}}
1107 ! When lmirep is set to .true., the spherical harmonic basis is transformed
1108 ! into one in which the site symmetries are block diagonal. Band characters
1109 ! determined from the density matrix expressed in this basis correspond to
1110 ! irreducible representations, and allow the partial DOS to be resolved into
1111 ! physically relevant contributions, for example $e_g$ and $t_{2g}$.
1112 !
1113 ! \block{lorbcnd}{
1114 ! {\tt lorbcnd} & {\tt .true.} if conduction state local-orbitals are to be
1115 ! automatically added to the basis & logical & {\tt .false.}}
1116 ! Adding these higher energy local-orbitals can improve calculations which
1117 ! rely on accurate unoccupied states, such as the response function. See also
1118 ! {\tt lorbordc}.
1119 !
1120 ! \block{lorbordc}{
1121 ! {\tt lorbordc} & the order of the conduction state local-orbitals &
1122 ! integer & 2}
1123 ! See {\tt lorbcnd}.
1124 !
1125 ! \block{lradstp}{
1126 ! {\tt lradstp} & radial step length for determining coarse radial mesh &
1127 ! integer & 4}
1128 ! Some muffin-tin functions (such as the density) are calculated on a coarse
1129 ! radial mesh and then interpolated onto a fine mesh. This is done for the
1130 ! sake of efficiency. {\tt lradstp} defines the step size in going from the
1131 ! fine to the coarse radial mesh. If it is too large, loss of precision may
1132 ! occur.
1133 !
1134 ! \block{maxitoep}{
1135 ! {\tt maxitoep} & maximum number of iterations when solving the exact
1136 ! exchange integral equations & integer & 300}
1137 ! See {\tt tau0oep}.
1138 !
1139 ! \block{maxscl}{
1140 ! {\tt maxscl} & maximum number of self-consistent loops allowed & integer &
1141 ! 200}
1142 ! This determines after how many loops the self-consistent cycle will
1143 ! terminate if the convergence criterion is not met. If {\tt maxscl} is $1$
1144 ! then the density and potential file, {\tt STATE.OUT}, will {\bf not} be
1145 ! written to disk at the end of the loop. See {\tt epspot}.
1146 !
1147 ! \block{mbwgrf}{
1148 ! {\tt mbwgrf} & matrix bandwidth of response functions in the
1149 ! ${\bf G}$-vector basis & integer & $-1$}
1150 ! Setting this to a positive integer results in response functions and the
1151 ! screened interaction $W({\bf G},{\bf G}',{\bf q},\omega)$ being treated as a
1152 ! banded matrix in ${\bf G}$ and ${\bf G}'$ with bandwidth {\tt mbwgrf}.
1153 ! This can be used to speed up $GW$ calculations.
1154 !
1155 ! \block{mixsave}{
1156 ! {\tt mixsave} & {\tt .true.} if the mixer work array is to be saved during a
1157 ! ground-state run & logical & {\tt .false.}}
1158 ! If {\tt mixsave} is {\tt .true.}, then the mixer work array is saved to
1159 ! {\tt MIXWORK.OUT} every {\tt nwrite} iterations and at the end of the
1160 ! self-consistent loop. This array is subsequently read in at the beginning of
1161 ! a restarted calculation in order to improve convergence.
1162 !
1163 ! \block{mixtype}{
1164 ! {\tt mixtype } & type of mixing required for the potential & integer & 3}
1165 ! Currently implemented are:
1166 ! \vskip 6pt
1167 ! \begin{tabularx}{\textwidth}[h]{lX}
1168 ! 0 & Linear mixing \\
1169 ! 1 & Adaptive linear mixing \\
1170 ! 3 & Broyden mixing, {\it J. Phys. A: Math. Gen.} {\bf 17}, L317 (1984)
1171 ! \end{tabularx}
1172 !
1173 ! \block{mixsdb}{
1174 ! {\tt mixsdb} & subspace dimension for Broyden mixing & integer & 5}
1175 ! This is the number of mixing vectors which define the subspace in which the
1176 ! Hessian matrix is calculated. See {\tt mixtype} and {\tt broydpm}.
1177 !
1178 ! \block{molecule}{
1179 ! {\tt molecule} & {\tt .true.} if the system is an isolated molecule &
1180 ! logical & {\tt .false.}}
1181 ! If {\tt molecule} is {\tt .true.}, then the atomic positions given in the
1182 ! {\tt atoms} block are assumed to be in Cartesian coordinates.
1183 !
1184 ! \block{momfix}{
1185 ! {\tt momfix} & the desired total moment for a FSM calculation &
1186 ! real(3) & $(0.0,0.0,0.0)$}
1187 ! Note that all three components must be specified (even for collinear
1188 ! calculations). Applies when {\tt fsmtype} is 1(-1) or 3(-3). See
1189 ! {\tt fsmtype}, {\tt taufsm} and {\tt spinpol}.
1190 !
1191 ! \block{momfixm}{
1192 ! {\tt momfixm} & the desired total moment magnitude for a FSM calculation &
1193 ! real & $0.0$}
1194 ! This applies when {\tt fsmtype} is 4 or 6.
1195 !
1196 ! \block{mommtfix}{
1197 ! {\tt is} & species number & integer & 0 \\
1198 ! {\tt ia} & atom number & integer & 0 \\
1199 ! {\tt mommtfix} & desired muffin-tin moment for a FSM calculation &
1200 ! real(3) & $(0.0,0.0,0.0)$}
1201 ! The local muffin-tin moments are specified for a subset of atoms, with the
1202 ! list terminated with a blank line. Note that all three components must be
1203 ! specified (even for collinear calculations). Applies when {\tt fsmtype} is
1204 ! 2(-2) or 3(-3). Note that the moment is not fixed in a particular muffin-tin
1205 ! when the magnitude of any component of the corresponding {\tt mommtfix} is
1206 ! $\ge 1000$. See {\tt fsmtype}, {\tt taufsm} and {\tt spinpol}.
1207 !
1208 ! \block{mommtfixm}{
1209 ! {\tt is} & species number & integer & 0 \\
1210 ! {\tt ia} & atom number & integer & 0 \\
1211 ! {\tt mommtfixm} & desired muffin-tin moment magnitude for a FSM
1212 ! calculation & real & $-1$}
1213 ! This applies when {\tt fsmtype} is 5 or 6. Note that the moment magnitude is
1214 ! not fixed in a particular muffin-tin when the corresponding {\tt mommtfixm}
1215 ! is negative.
1216 !
1217 ! \block{mrmtav}{
1218 ! {\tt mrmtav} & order of averaging applied to the muffin-tin radii &
1219 ! integer & 0}
1220 ! Crystal structures with muffin-tin radii which are widely varying in size
1221 ! can cause calculations to become unstable. Applying a simple averaging
1222 ! procedure to the radii reduces this variation and can improve stability.
1223 ! The larger {\tt mrmtav}, the more equal the muffin-tin radii will become.
1224 ! See the routine {\tt rmtavrg}.
1225 !
1226 ! \block{mstar}{
1227 ! {\tt mstar} & value of the effective mass parameter used for adaptive
1228 ! determination of {\tt swidth} & real & $10.0$}
1229 ! See {\tt autoswidth}.
1230 !
1231 ! \block{mustar}{
1232 ! {\tt mustar} & Coulomb pseudopotential, $\mu^*$, used in the
1233 ! McMillan-Allen-Dynes equation & real & $0.15$}
1234 ! This is used when calculating the superconducting critical temperature with
1235 ! the formula {\it Phys. Rev. B 12, 905 (1975)}
1236 ! $$ T_c=\frac{\omega_{\rm log}}{1.2 k_B}\exp\left[\frac{-1.04(1+\lambda)}
1237 ! {\lambda-\mu^*(1+0.62\lambda)}\right], $$
1238 ! where $\omega_{\rm log}$ is the logarithmic average frequency and $\lambda$
1239 ! is the electron-phonon coupling constant.
1240 !
1241 ! \block{ncbse}{
1242 ! {\tt ncbse} & number of conduction states to be used for BSE calculations &
1243 ! integer & 3}
1244 ! See also {\tt nvbse}.
1245 !
1246 ! \block{ndspem}{
1247 ! {\tt ndspem} & the number of {\bf k}-vector displacements in each direction
1248 ! around {\tt vklem} when computing the numerical derivatives for the
1249 ! effective mass tensor & integer & 1}
1250 ! See {\tt deltaem} and {\tt vklem}.
1251 !
1252 ! \block{nempty}{
1253 ! {\tt nempty} & the number of empty states per atom and spin & real & $4.0$ }
1254 ! Defines the number of eigenstates beyond that required for charge
1255 ! neutrality. When running metals it is not known {\it a priori} how many
1256 ! states will be below the Fermi energy for each $k$-point. Setting
1257 ! {\tt nempty} greater than zero allows the additional states to act as a
1258 ! buffer in such cases. Furthermore, magnetic calculations use the
1259 ! first-variational eigenstates as a basis for setting up the
1260 ! second-variational Hamiltonian, and thus {\tt nempty} will determine the
1261 ! size of this basis set. Convergence with respect to this quantity should be
1262 ! checked.
1263 !
1264 ! \block{ngridk}{
1265 ! {\tt ngridk } & the $k$-point mesh sizes & integer(3) & $(1,1,1)$}
1266 ! The ${\bf k}$-vectors are generated using
1267 ! $$ {\bf k}=(\frac{i_1+v_1}{n_1},\frac{i_2+v_2}{n_2},\frac{i_3+v_3}{n_3}), $$
1268 ! where $i_j$ runs from 0 to $n_j-1$ and $0\le v_j<1$ for $j=1,2,3$. The
1269 ! vector ${\bf v}$ is given by the variable {\tt vkloff}. See also
1270 ! {\tt reducek}.
1271 !
1272 ! \block{ngridq}{
1273 ! {\tt ngridq } & the phonon $q$-point mesh sizes & integer(3) & $(1,1,1)$}
1274 ! Same as {\tt ngridk}, except that this mesh is for the phonon $q$-points
1275 ! and other tasks. See also {\tt reduceq}.
1276 !
1277 ! \block{nosource}{
1278 ! {\tt nosource} & when set to {\tt .true.}, source fields are projected out
1279 ! of the exchange-correlation magnetic field & logical & {\tt .false.}}
1280 ! Experimental feature.
1281 !
1282 ! \block{notes}{
1283 ! {\tt notes(i)} & the $i$th line of the notes & string & -}
1284 ! This block allows users to add their own notes to the file {\tt INFO.OUT}.
1285 ! The block should be terminated with a blank line, and no line should exceed
1286 ! 80 characters.
1287 !
1288 ! \block{npmae}{
1289 ! {\tt npmae } & number or distribution of directions for MAE calculations &
1290 ! integer & $-1$}
1291 ! Automatic determination of the magnetic anisotropy energy (MAE) requires
1292 ! that the total energy is determined for a set of directions of the total
1293 ! magnetic moment. This variable controls the number or distribution of these
1294 ! directions. The convention is:
1295 ! \vskip 6pt
1296 ! \begin{tabularx}{\textwidth}[h]{lX}
1297 ! $-4,-3,-2,-1$ & Cardinal directions given by the primitive translation
1298 ! vectors $n_1{\bf A}_1+n_2{\bf A}_2+n_3{\bf A}_3$, where
1299 ! $1\le n_i\le|{\tt npmae}|$ \\
1300 ! 2 & Cartesian $x$ and $z$ directions \\
1301 ! 3 & Cartesian $x$, $y$ and $z$ directions \\
1302 ! $4,5,\ldots$ & Even distribution of {\tt npmae} directions
1303 ! \end{tabularx}
1304 !
1305 ! \block{ntemp}{
1306 ! {\tt ntemp} & number of temperature steps & integer & 40}
1307 ! This is the number of temperature steps to be used in the Eliashberg gap
1308 ! and thermodynamic properties calculations.
1309 !
1310 ! \block{num\_wann}{
1311 ! {\tt num\_wann} & number of Wannier90 wavefunctions & integer & 0}
1312 ! If ${\tt num\_wann}>0$ then this is the number of Wannier wavefunctions
1313 ! to be found by the Wannier90 package. If ${\tt num\_wann}\le 0$ then the
1314 ! number of wavefunctions is given by ${\tt num\_bands}+{\tt num\_wann}$.
1315 !
1316 ! \block{nvbse}{
1317 ! {\tt nvbse} & number of valence states to be used for BSE calculations &
1318 ! integer & 2}
1319 ! See also {\tt ncbse}.
1320 !
1321 ! \block{nwrite}{
1322 ! {\tt nwrite} & number of self-consistent loops after which {\tt STATE.OUT}
1323 ! is to be written & integer & 0}
1324 ! Normally, the density and potentials are written to the file {\tt STATE.OUT}
1325 ! only after completion of the self-consistent loop. By setting {\tt nwrite}
1326 ! to a positive integer the file will instead be written every {\tt nwrite}
1327 ! loops.
1328 !
1329 ! \block{nxoapwlo}{
1330 ! {\tt nxoapwlo} & extra order of radial functions to be added to the existing
1331 ! APW and local-orbital set & integer & 0}
1332 ! Setting this variable will result in the APWs and local-orbitals for all
1333 ! species becoming higher order with corresponding increase in derivative
1334 ! matching at the muffin-tin surface. For example, setting {\tt nxoapwlo}=1
1335 ! turns all APWs into LAPWs.
1336 !
1337 ! \block{optcomp}{
1338 ! {\tt optcomp} & the components of the first- or second-order optical tensor
1339 ! to be calculated & integer(3) & $(1,1,1)$}
1340 ! This selects which components of the optical tensor you would like to plot.
1341 ! Only the first two are used for the first-order tensor. Several components
1342 ! can be listed one after the other with a blank line terminating the list.
1343 !
1344 ! \block{phwrite}{
1345 ! {\tt nphwrt} & number of $q$-points for which phonon modes are to be found &
1346 ! integer & 1 \\
1347 ! \hline
1348 ! {\tt vqlwrt(i)} & the $i$th $q$-point in lattice coordinates & real(3) &
1349 ! $(0.0,0.0,0.0)$}
1350 ! This is used in conjunction with {\tt task}=230. The code will write the
1351 ! phonon frequencies and eigenvectors to the file {\tt PHONON.OUT} for all the
1352 ! $q$-points in the list. The $q$-points can be anywhere in the Brillouin zone
1353 ! and do not have to lie on the mesh defined by {\tt ngridq}. Obviously, all
1354 ! the dynamical matrices have to be computed first using {\tt task}=200.
1355 !
1356 ! \block{plot1d}{
1357 ! {\tt nvp1d} & number of vertices & integer & 2 \\
1358 ! {\tt npp1d} & number of plotting points & integer & 200 \\
1359 ! \hline
1360 ! {\tt vvlp1d(i)} & lattice coordinates for vertex $i$ & real(3) &
1361 ! $(0.0,0.0,0.0)\rightarrow(1.0,1.0,1.0)$}
1362 ! Defines the path in either real or reciprocal space along which the 1D plot
1363 ! is to be produced. The user should provide {\tt nvp1d} vertices in lattice
1364 ! coordinates.
1365 !
1366 ! \block{plot2d}{
1367 ! {\tt vclp2d(0)} & zeroth corner (origin) & real(3) & $(0.0,0.0,0.0)$ \\
1368 ! \hline
1369 ! {\tt vclp2d(1)} & first corner & real(3) & $(1.0,0.0,0.0)$ \\
1370 ! \hline
1371 ! {\tt vclp2d(2)} & second corner & real(3) & $(0.0,1.0,0.0)$ \\
1372 ! \hline
1373 ! {\tt np2d} & number of plotting points in both directions & integer(2) &
1374 ! $(40,40)$}
1375 ! Defines the corners of a parallelogram and the grid size used for producing
1376 ! 2D plots.
1377 !
1378 ! \block{plot3d}{
1379 ! {\tt vclp3d(0)} & zeroth corner (origin) & real(3) & $(0.0,0.0,0.0)$ \\
1380 ! \hline
1381 ! {\tt vclp3d(1)} & first corner & real(3) & $(1.0,0.0,0.0)$ \\
1382 ! \hline
1383 ! {\tt vclp3d(2)} & second corner & real(3) & $(0.0,1.0,0.0)$ \\
1384 ! \hline
1385 ! {\tt vclp3d(3)} & third corner & real(3) & $(0.0,0.0,1.0)$ \\
1386 ! \hline
1387 ! {\tt np3d} & number of plotting points each direction & integer(3) &
1388 ! $(20,20,20)$}
1389 ! Defines the corners of a box and the grid size used for producing 3D plots.
1390 !
1391 ! \block{primcell}{
1392 ! {\tt primcell} & {\tt .true.} if the primitive unit cell should be found
1393 ! & logical & {\tt .false.}}
1394 ! Allows the primitive unit cell to be determined automatically from the
1395 ! conventional cell. This is done by searching for lattice vectors among all
1396 ! those which connect atomic sites, and using the three shortest which produce
1397 ! a unit cell with non-zero volume.
1398 !
1399 ! \block{pulse}{
1400 ! {\tt n} & number of pulses & integer & - \\
1401 ! \hline
1402 ! {\tt a0(i)} & polarisation vector (including amplitude) & real(3) & - \\
1403 ! {\tt w(i)} & frequency & real & - \\
1404 ! {\tt phi(i)} & phase in degrees & real & - \\
1405 ! {\tt rc(i)} & chirp rate & real & - \\
1406 ! {\tt t0(i)} & peak time & real & - \\
1407 ! {\tt d(i)} & full-width at half-maximum & real & -}
1408 ! Parameters used to generate a time-dependent vector potential ${\bf A}(t)$
1409 ! representing a laser pulse. The total vector potential is the sum of
1410 ! individual pulses and is given by the formula
1411 ! $$ {\bf A}(t)=\sum_{i=1}^n {\bf A}_0^i\exp
1412 ! \left[-(t-t_0^i)^2/2\sigma_i^2\right]
1413 ! \sin\left[w_i(t-t_0^i)+\phi_i+r_{\rm c}^i t^2/2\right], $$
1414 ! where $\sigma=d/2\sqrt{2\ln 2}$. See also {\tt ramp}.
1415 !
1416 ! \block{q0cut}{
1417 ! {\tt q0cut} & Q-vector cut-off for the ultra long-range Coulomb Green's
1418 ! function & real & $0.0$}
1419 ! The ultra long-range Poisson equation is solved using the Green's function
1420 ! $4\pi/|{\bf G}+{\bf Q}|^2$. This is set to zero for all
1421 ! $|{\bf G}+{\bf Q}|<{\tt q0cut}$ when {\tt q0cut} is positive. For negative
1422 ! {\tt q0cut}, a Yukawa-type screening of the form
1423 ! $4\pi/(|{\bf G}+{\bf Q}|^2+{\tt q0cut}^2)$ is employed. Setting this
1424 ! variable to be small but finite can improve the stability of a
1425 ! self-consistent calculation.
1426 !
1427 ! \block{radkpt}{
1428 ! {\tt radkpt} & radius of sphere used to determine $k$-point density &
1429 ! real & $40.0$}
1430 ! Used for the automatic determination of the $k$-point mesh. If {\tt autokpt}
1431 ! is set to {\tt .true.} then the mesh sizes will be determined by
1432 ! $n_i=R_k|{\bf B}_i|+1$, where ${\bf B}_i$ are the primitive reciprocal
1433 ! lattice vectors.
1434 !
1435 ! \block{ramp}{
1436 ! {\tt n} & number of ramps & integer & - \\
1437 ! \hline
1438 ! {\tt a0(i)} & polarisation vector (including amplitude) & real(3) & - \\
1439 ! {\tt t0(i)} & ramp start time & real & - \\
1440 ! {\tt c1(i)} & linear coefficient of ${\bf A}(t)$ & real & - \\
1441 ! {\tt c2(i)} & quadratic coefficient & real & -}
1442 ! Parameters used to generate a time-dependent vector potential ${\bf A}(t)$
1443 ! representing a constant or linearly increasing electric field
1444 ! ${\bf E}(t)=-\partial{\bf A}(t)/\partial t$. The vector potential is given
1445 ! by
1446 ! $$ {\bf A}(t)=\sum_{i=1}^n {\bf A}_0^i
1447 ! \left[c_1(t-t_0)+c_2(t-t_0)^2\right]\Theta(t-t_0). $$
1448 !
1449 ! \block{readadu}{
1450 ! {\tt readadu} & set to {\tt .true.} if the interpolation constant for
1451 ! DFT+$U$ should be read from file rather than calculated & logical &
1452 ! {\tt .false.}}
1453 ! When {\tt dftu}=3, the DFT+$U$ energy and potential are interpolated
1454 ! between FLL and AFM. The interpolation constant, $\alpha$, is normally
1455 ! calculated from the density matrix, but can also be read in from the file
1456 ! {\tt ALPHADU.OUT}. This allows the user to fix $\alpha$, but is also
1457 ! necessary when calculating forces, since the contribution of the potential
1458 ! of the variation of $\alpha$ with respect to the density matrix is not
1459 ! computed. See {\tt dft+u}.
1460 !
1461 ! \block{reducebf}{
1462 ! {\tt reducebf} & reduction factor for the external magnetic fields & real &
1463 ! $1.0$}
1464 ! After each self-consistent loop, the external magnetic fields are multiplied
1465 ! with {\tt reducebf}. This allows for a large external magnetic field at the
1466 ! start of the self-consistent loop to break spin symmetry, while at the end
1467 ! of the loop the field will be effectively zero, i.e. infinitesimal. See
1468 ! {\tt bfieldc} and {\tt atoms}.
1469 !
1470 ! \block{reduceh}{
1471 ! {\tt reduceh} & set to {\tt .true.} if the reciprocal ${\bf H}$-vectors
1472 ! should be reduced by the symmorphic crystal symmetries & logical & .true.}
1473 ! See {\tt hmaxvr} and {\tt vmat}.
1474 !
1475 ! \block{reducek}{
1476 ! {\tt reducek} & type of reduction of the $k$-point set & integer & 1}
1477 ! Types of reduction are defined by the symmetry group used:
1478 ! \vskip 6pt
1479 ! \begin{tabularx}{\textwidth}[h]{lX}
1480 ! 0 & no reduction \\
1481 ! 1 & reduce with full crystal symmetry group (including non-symmorphic
1482 ! symmetries) \\
1483 ! 2 & reduce with symmorphic symmetries only
1484 ! \end{tabularx}
1485 ! \vskip 6pt
1486 ! See also {\tt ngridk} and {\tt vkloff}.
1487 !
1488 ! \block{reduceq}{
1489 ! {\tt reduceq} & type of reduction of the $q$-point set & integer & 1}
1490 ! See {\tt reducek} and {\tt ngridq}.
1491 !
1492 ! \block{rgkmax}{
1493 ! {\tt rgkmax} & $R^{\rm MT}_{\rm min}\times\max\{|{\bf G}+{\bf k}|\}$ &
1494 ! real & $7.0$}
1495 ! This sets the maximum length for the ${\bf G}+{\bf k}$ vectors, defined as
1496 ! {\tt rgkmax} divided by the average muffin-tin radius. See {\tt isgkmax}.
1497 !
1498 ! \block{rmtall}{
1499 ! {\tt rmtall} & muffin-tin radius for all species & real & $-1.0$}
1500 ! If {\tt rmtall} is positive then all muffin-tin radii are set to this value.
1501 !
1502 ! \block{rmtdelta}{
1503 ! {\tt rmtdelta} & minimum allowed distance between muffin-tin surfaces &
1504 ! real & $0.05$}
1505 !
1506 ! \block{rmtscf}{
1507 ! {\tt rmtscf} & muffin-tin radius scaling factor & real & $1.0$}
1508 ! All muffin-tin radii read from the species files are scaled by this factor.
1509 !
1510 ! \block{rndavec}{
1511 ! {\tt rndavec} & lattice vector randomisation amplitude & real & $0.0$}
1512 ! Setting this to a number larger than zero causes the code add random jitter
1513 ! to the lattice vectors proportional to this amplitude. This can be used to
1514 ! break lattice symmetry.
1515 !
1516 ! \block{rotavec}{
1517 ! {\tt axang} & axis-angle representation of lattice vector rotation &
1518 ! real(4) & $(0.0,0.0,0.0,0.0)$}
1519 ! This determines the rotation matrix which is applied to the lattice vectors
1520 ! prior to any calculation. The first three components specify the axis and
1521 ! the last component is the angle in degrees. The `right-hand rule' convention
1522 ! is followed.
1523 !
1524 ! \block{scale}{
1525 ! {\tt scale } & lattice vector scaling factor & real & $1.0$}
1526 ! Scaling factor for all three lattice vectors. Applied in conjunction with
1527 ! {\tt scale1}, {\tt scale2} and {\tt scale3}.
1528 !
1529 ! \block{scale1/2/3}{
1530 ! {\tt scale1/2/3 } & separate scaling factors for each lattice vector &
1531 ! real & $1.0$}
1532 !
1533 ! \block{scissor}{
1534 ! {\tt scissor} & the scissor correction & real & $0.0$}
1535 ! This is the scissor shift applied to states above the Fermi energy
1536 ! {\it Phys. Rev. B} {\bf 43}, 4187 (1991). Affects optics calculations only.
1537 !
1538 ! \block{scrpath}{
1539 ! {\tt scrpath} & scratch space path & string & null}
1540 ! This is the scratch space path where the eigenvector files {\tt EVALFV.OUT}
1541 ! and {\tt EVALSV.OUT} will be written. If the run directory is accessed via a
1542 ! network then {\tt scrpath} can be set to a directory on the local disk, for
1543 ! example {\tt /tmp/}. Note that the forward slash {\tt /} at the end of the
1544 ! path must be included.
1545 !
1546 ! \block{socscf}{
1547 ! {\tt socscf} & scaling factor for the spin-orbit coupling term in the
1548 ! Hamiltonian & real & $1.0$}
1549 ! This can be used to enhance the effect of spin-orbit coupling in order to
1550 ! accurately determine the magnetic anisotropy energy (MAE).
1551 !
1552 ! \block{spincore}{
1553 ! {\tt spincore} & set to {\tt .true.} if the core should be spin-polarised
1554 ! & logical & {\tt .false.}}
1555 !
1556 ! \block{spinorb}{
1557 ! {\tt spinorb} & set to {\tt .true.} if a spin-orbit coupling is required
1558 ! & logical & {\tt .false.}}
1559 ! If {\tt spinorb} is {\tt .true.}, then a $\boldsymbol\sigma\cdot{\bf L}$
1560 ! term is added to the second-variational Hamiltonian. See {\tt spinpol}.
1561 !
1562 ! \block{spinpol}{
1563 ! {\tt spinpol} & set to {\tt .true.} if a spin-polarised calculation is
1564 ! required & logical & {\tt .false.}}
1565 ! If {\tt spinpol} is {\tt .true.}, then the spin-polarised Hamiltonian is
1566 ! solved as a second-variational step using two-component spinors in the
1567 ! Kohn-Sham magnetic field. The first variational scalar wavefunctions are
1568 ! used as a basis for setting this Hamiltonian.
1569 !
1570 ! \block{spinsprl}{
1571 ! {\tt spinsprl} & set to {\tt .true.} if a spin-spiral calculation is
1572 ! required & logical & {\tt .false.}}
1573 ! Experimental feature for the calculation of spin-spiral states. See
1574 ! {\tt vqlss} for details.
1575 !
1576 ! \block{sppath}{
1577 ! {\tt sppath} & path where the species files can be found & string & null}
1578 ! Note that the forward slash {\tt /} at the end of the path must be included.
1579 !
1580 ! \block{ssdph}{
1581 ! {\tt ssdph} & set to {\tt .true.} if a complex de-phasing factor is to be
1582 ! used in spin-spiral calculations & logical & {\tt .true.}}
1583 ! If this is {\tt .true.} then spin-spiral wavefunctions in each muffin-tin at
1584 ! position ${\bf r}_{\alpha}$ are de-phased by the matrix
1585 ! $$ \begin{pmatrix} e^{-i{\bf q}\cdot{\bf r}_{\alpha}/2} & 0 \\
1586 ! 0 & e^{i{\bf q}\cdot{\bf r}_{\alpha}/2} \end{pmatrix}. $$
1587 ! In simple situations, this has the advantage of producing magnon dynamical
1588 ! matrices which are already in diagonal form. This option should be used with
1589 ! care, and a full understanding of the spin-spiral configuration is required.
1590 ! See {\tt spinsprl}.
1591 !
1592 ! \block{stype}{
1593 ! {\tt stype} & integer defining the type of smearing to be used & integer &
1594 ! $3$}
1595 ! A smooth approximation to the Dirac delta function is needed to compute the
1596 ! occupation numbers of the Kohn-Sham states. The variable {\tt swidth}
1597 ! determines the width of the approximate delta function. Currently
1598 ! implemented are
1599 ! \vskip 6pt
1600 ! \begin{tabularx}{\textwidth}[h]{lX}
1601 ! 0 & Gaussian \\
1602 ! 1 & Methfessel-Paxton order 1, Phys. Rev. B {\bf 40}, 3616 (1989) \\
1603 ! 2 & Methfessel-Paxton order 2 \\
1604 ! 3 & Fermi-Dirac
1605 ! \end{tabularx}
1606 ! \vskip 6pt
1607 ! See also {\tt autoswidth}, {\tt swidth} and {\tt tempk}.
1608 !
1609 ! \block{swidth}{
1610 ! {\tt swidth} & width of the smooth approximation to the Dirac delta
1611 ! function & real & $0.001$}
1612 ! See {\tt stype} for details and the variable {\tt tempk}.
1613 !
1614 ! \newpage
1615 ! \block{tasks}{
1616 ! {\tt task(i) } & the $i$th task & integer & $-1$}
1617 ! A list of tasks for the code to perform sequentially. The list should be
1618 ! terminated with a blank line. Each task has an associated integer as
1619 ! follows:
1620 ! \vskip 6pt
1621 ! \begin{tabularx}{\textwidth}[h]{lX}
1622 ! 0 & Ground-state run starting from the atomic densities. \\
1623 ! 1 & Resumption of ground-state run using density in {\tt STATE.OUT}. \\
1624 ! 2 & Geometry optimisation run starting from the atomic densities, with
1625 ! atomic positions written to {\tt GEOMETRY.OUT}. \\
1626 ! 3 & Resumption of geometry optimisation run using density in {\tt STATE.OUT}
1627 ! but with positions from {\tt elk.in}. \\
1628 ! 5 & Ground-state Hartree-Fock run. \\
1629 ! 10 & Total, partial and interstitial density of states (DOS). \\
1630 ! 14 & Plots the smooth Dirac delta and Heaviside step functions used by the
1631 ! code to calculate occupation numbers. \\
1632 ! 15 & Output ${\bf L}$, ${\bf S}$ and ${\bf J}$ total expectation values. \\
1633 ! 16 & Output ${\bf L}$, ${\bf S}$ and ${\bf J}$ expectation values for each
1634 ! $k$-point and state in {\tt kstlist}. \\
1635 ! 20 & Band structure plot. \\
1636 ! 21 & Band structure plot which includes total and angular momentum
1637 ! characters for every atom. \\
1638 ! 22 & Band structure plot which includes $(l,m)$ character for every atom. \\
1639 ! 23 & Band structure plot which includes spin character for every atom. \\
1640 ! 25 & Compute the effective mass tensor at the $k$-point given by
1641 ! {\tt vklem}. \\
1642 ! 31/2/3 & 1/2/3D charge density plot. \\
1643 ! 41/2/3 & 1/2/3D exchange-correlation and Coulomb potential plots. \\
1644 ! 51/2/3 & 1/2/3D electron localisation function (ELF) plot. \\
1645 ! 61/2/3 & 1/2/3D wavefunction plot:
1646 ! $\left|\Psi_{i{\bf k}}({\bf r})\right|^2$. \\
1647 ! 65 & Write the core wavefunctions to file for plotting. \\
1648 ! 68 & Output the status of the RAM disk. \\
1649 ! 71/2/3 & 1/2/3D plot of magnetisation vector field, ${\bf m}({\bf r})$. \\
1650 ! 81/2/3 & 1/2/3D plot of exchange-correlation magnetic vector field,
1651 ! ${\bf B}_{\rm xc}({\bf r})$. \\
1652 ! 91/2/3 & 1/2/3D plot of $\nabla\cdot{\bf B}_{\rm xc}({\bf r})$. \\
1653 ! 100 & 3D Fermi surface plot using the scalar product
1654 ! $p({\bf k})=\prod_i(\epsilon_{i{\bf k}}-\epsilon_{\rm F})$. \\
1655 ! 101 & 3D Fermi surface plot using separate bands (minus the Fermi
1656 ! energy). \\
1657 ! 102 & 3D Fermi surface which can be plotted with XCrysDen. \\
1658 ! 105 & 3D nesting function plot. \\
1659 ! 110 & Calculation of M\"{o}ssbauer contact charge densities and magnetic
1660 ! fields at the nuclear sites. \\
1661 ! 115 & Calculation of the electric field gradient (EFG) at the nuclear
1662 ! sites. \\
1663 ! 120 & Output of the momentum matrix elements
1664 ! $\langle\Psi_{i{\bf k}}|-i\nabla|\Psi_{j{\bf k}}\rangle$. \\
1665 ! 121 & Linear optical dielectric response tensor calculated within the random
1666 ! phase approximation (RPA) and in the $q\rightarrow 0$ limit, with no
1667 ! microscopic contributions. \\
1668 ! 122 & Magneto optical Kerr effect (MOKE) angle. \\
1669 ! 125 & Non-linear optical second harmonic generation. \\
1670 ! 130 & Output matrix elements of the type
1671 ! $\langle\Psi_{i{\bf k+q}}|e^{i{\bf q}\cdot{\bf r}}|
1672 ! \Psi_{j{\bf k}}\rangle$.
1673 ! \end{tabularx}
1674 !
1675 ! \begin{tabularx}{\textwidth}[h]{lX}
1676 ! 135 & Output all wavefunctions expanded in the plane wave basis up to a
1677 ! cut-off defined by {\tt hkmax}. \\
1678 ! 140 & Energy loss near edge structure (ELNES). \\
1679 ! 141/2/3 & 1/2/3D plot of the electric field
1680 ! ${\bf E}({\bf r})\equiv\nabla V_{\rm C}({\bf r})$. \\
1681 ! 150 & Write out the atomic eigenvalues for each species. \\
1682 ! 151/2/3 & 1/2/3D plot of
1683 ! ${\bf m}({\bf r})\times{\bf B}_{\rm xc}({\bf r})$. \\
1684 ! 160 & Calculates the total exchange-correlation spin-torque acting on the
1685 ! system: $\boldsymbol\tau=\int d^3r\, {\bf m}({\bf r})\times
1686 ! {\bf B}_{\rm xc}({\bf r}).$ \\
1687 ! 162 & 2D scanning-tunneling microscopy (STM) image. \\
1688 ! 170 & Writes the electron momentum density to {\tt EMD.OUT}. \\
1689 ! 171/2/3 & 1/2/3D plot of the electron momentum density. \\
1690 ! 180 & Generate the RPA inverse dielectric function with local contributions
1691 ! $\epsilon^{-1}({\bf G},{\bf G}',{\bf q},\omega)$ and write it to file. \\
1692 ! 185 & Write the Bethe-Salpeter equation (BSE) Hamiltonian to file. \\
1693 ! 186 & Diagonalise the BSE Hamiltonian and write the eigenvectors and
1694 ! eigenvalues to file. \\
1695 ! 187 & Output the BSE dielectric response function. \\
1696 ! 190 & Write the atomic geometry to file for plotting with XCrySDen and
1697 ! V\_Sim. \\
1698 ! 195 & Calculation of X-ray density structure factors. \\
1699 ! 196 & Calculation of magnetic structure factors. \\
1700 ! 200 & Calculation of phonon dynamical matrices on a $q$-point set defined by
1701 ! {\tt ngridq} using the supercell method. \\
1702 ! 202 & Phonon dry run: just produce a set of empty {\tt DYN} files. \\
1703 ! 205 & Calculation of phonon dynamical matrices using density functional
1704 ! perturbation theory (DFPT). \\
1705 ! 208 & Calculation of static Born effective charges. \\
1706 ! 209 & Born effective charge dry run: produce a set of empty {\tt BEC}
1707 ! files. \\
1708 ! 210 & Phonon density of states. \\
1709 ! 220 & Phonon dispersion plot. \\
1710 ! 230 & Phonon frequencies and eigenvectors for an arbitrary $q$-point. \\
1711 ! 240/1 & Generate the ${\bf q}$-dependent phonon linewidths and
1712 ! electron-phonon coupling constants and write them to file. Task 241 also
1713 ! writes the complete set of electron-phonon coupling matrix elements to
1714 ! {\tt EPHMAT.OUT}. \\
1715 ! 245 & Phonon linewidths plot. \\
1716 ! 250 & Eliashberg function $\alpha^2F(\omega)$, electron-phonon coupling
1717 ! constant $\lambda$, and the McMillan-Allen-Dynes critical temperature
1718 ! $T_c$. \\
1719 ! 260 & Solves the Eliashberg equations to find the superconducting gap. \\
1720 ! 270 & Electron-phonon Bogoliubov equation ground-state. \\
1721 ! 271 & Resumption of the Bogoliubov ground-state. \\
1722 ! 280 & Electron-phonon Bogoliubov density of states. \\
1723 ! 285 & Electron and phonon anomalous correlation entropy (ACE). \\
1724 ! 300 & Reduced density matrix functional theory (RDMFT) calculation. \\
1725 ! 320 & Time-dependent density functional theory (TDDFT) calculation of the
1726 ! dielectric response function including microscopic contributions. \\
1727 ! 330/1 & TDDFT calculation of the spin-polarised response function for
1728 ! arbitrary ${\bf q}$-vectors. Task 331 writes the entire response function
1729 ! $\overleftrightarrow{\chi}({\bf G},{\bf G}',{\bf q},\omega)$ to file.
1730 ! \end{tabularx}
1731 !
1732 ! \begin{tabularx}{\textwidth}[h]{lX}
1733 ! 341/2/3 & 1/2/3D plot of
1734 ! $w_{\rm xc}({\bf r})\equiv \left.\delta E_{\rm xc}[\rho,\tau]/
1735 ! \delta \tau({\bf r})\right|_{\rho}$ which is calculated for meta-GGA
1736 ! functionals. \\
1737 ! 350/1/2 & Spin-spiral supercell calculations (spin-orbit coupling can be
1738 ! included). Task 352 is a dry run: produce a set of empty {\tt SS} files. \\
1739 ! 371/2/3 & 1/2/3D plot of the paramagnetic current density
1740 ! ${\bf j}_{\rm p}({\bf r})$. \\
1741 ! 380 & Piezoelectric tensor. \\
1742 ! 390 & Magnetoelectric tensor. \\
1743 ! 400 & Calculation of tensor moments and corresponding DFT+$U$ Hartree-Fock
1744 ! energy contributions. \\
1745 ! 420/1 & Molecular dynamics (MD) calculation within the adiabatic
1746 ! approximation. Task 421 restarts an interrupted MD calculation. \\
1747 ! 430 & Write the strain tensors to {\tt STRAIN.OUT}. \\
1748 ! 440 & Write the stress tensor components corresponding to the strain tensors
1749 ! to {\tt STRESS.OUT}. \\
1750 ! 450 & Generates a laser pulse in the form of a time-dependent vector
1751 ! potential ${\bf A}(t)$ and writes it to {\tt AFIELDT.OUT}. \\
1752 ! 455 & Writes the time-dependent power density and total energy density of
1753 ! ${\bf A}(t)$. \\
1754 ! 456 & Writes the frequency-dependent electric field ${\bf E}(\omega)$
1755 ! corresponding to ${\bf A}(t)$. \\
1756 ! 460/1/2/3 & Time evolution run using TDDFT under the influence of
1757 ! ${\bf A}(t)$. Tasks 462 and 463 include nuclear Ehrenfest dynamics. Tasks
1758 ! 461 and 463 restart interrupted calculations. \\
1759 ! 471 & 1D plot of the static charge density. \\
1760 ! 478 & Calculation of the dynamical Born effective charges. \\
1761 ! 480/1 & Computes the dielectric function from the time-dependent current
1762 ! density in {\tt JTOT\_TD.OUT}. Task 481 assumes that the vector potential
1763 ! ${\bf A}(t)$ is a step function at $t=0$, while task 480 makes no such
1764 ! assumptions. \\
1765 ! 500 & Checks the test files generated when {\tt test} is {\tt .true.} \\
1766 ! 550 & Writes the files required by Wannier90. \\
1767 ! 600 & Output the $GW$ self-energy matrix elements. \\
1768 ! 610 & Generates the $GW$ spectral function. \\
1769 ! 620 & Generates the $GW$ band structure, i.e. the ${\bf k}$-dependent
1770 ! spectral function). \\
1771 ! 630 & Writes the $GW$ Fermi energy to {\tt GWEFERMI.OUT} \\
1772 ! 640 & Determines the $GW$ density matrix in terms of natural orbitals and
1773 ! occupation numbers. The files {\tt EVECSV.OUT} and {\tt OCCSV.OUT} are
1774 ! overwritten by these and can then be used for other calculations. \\
1775 ! 700/1 & Ultra long-range ground-state calculation. Task 701 is for
1776 ! restarting from the potential in {\tt STATE\_ULR.OUT}. \\
1777 ! 720/5 & Ultra long-range band structure and spectral function. Task 720 is
1778 ! for the central $k$-point (i.e. $\kappa=0$), while task 725 averages over
1779 ! all $\kappa$-points.\\
1780 ! 731/2/3 & 1/2/3D plot of the ultra long-range density. \\
1781 ! 741/2/3 & 1/2/3D plot of the ultra long-range Kohn-Sham potential. \\
1782 ! 771/2/3 & 1/2/3D plot of the ultra long-range magnetisation.
1783 ! \end{tabularx}
1784 !
1785 ! \block{tau0atp}{
1786 ! {\tt tau0atp} & the step size to be used for atomic position optimisation &
1787 ! real & $0.25$}
1788 ! The position of atom $\alpha$ is updated on step $m$ of a geometry
1789 ! optimisation run using
1790 ! $$ {\bf r}_{\alpha}^{m+1}={\bf r}_{\alpha}^m+\tau_{\alpha}^m
1791 ! \left({\bf F}_{\alpha}^m+{\bf F}_{\alpha}^{m-1}\right), $$
1792 ! where $\tau_{\alpha}$ is set to {\tt tau0atp} for $m=0$, and incremented by
1793 ! the same amount if the atom is moving in the same direction between steps.
1794 ! If the direction changes then $\tau_{\alpha}$ is reset to {\tt tau0atp}.
1795 !
1796 ! \block{tau0latv}{
1797 ! {\tt tau0latv} & the step size to be used for lattice vector optimisation &
1798 ! real & $0.25$}
1799 ! This parameter is used for lattice vector optimisation in a procedure
1800 ! identical to that for atomic position optimisation. See {\tt tau0atp} and
1801 ! {\tt latvopt}.
1802 !
1803 ! \block{tau0oep}{
1804 ! {\tt tau0oep} & initial step length for the OEP iterative solver & real &
1805 ! $0.5$}
1806 ! The optimised effective potential is determined using an interative method
1807 ! [Phys. Rev. Lett. 98, 196405 (2007)]. This variable sets the step length as
1808 ! described in the article. See {\tt maxitoep}.
1809 !
1810 ! \block{taufsm}{
1811 ! {\tt taufsm} & the step size to be used when finding the effective magnetic
1812 ! field in fixed spin moment calculations & real & $0.01$}
1813 ! An effective magnetic field, ${\bf B}_{\rm FSM}$, is required for fixing the
1814 ! spin moment to a given value, ${\bf M}_{\rm FSM}$. This is found by adding a
1815 ! vector to the field which is proportional to the difference between the
1816 ! moment calculated in the $i$th self-consistent loop and the required moment:
1817 ! $$ {\bf B}_{\rm FSM}^{i+1}={\bf B}_{\rm FSM}^i+\lambda\left({\bf M}^i
1818 ! -{\bf M}_{\rm FSM}\right), $$
1819 ! where $\lambda$ is proportional to {\tt taufsm}. See also {\tt fsmtype},
1820 ! {\tt momfix} and {\tt spinpol}.
1821 !
1822 ! \block{tempk}{
1823 ! {\tt tempk} & temperature $T$ of the electronic system in kelvin & real & -}
1824 ! Assigning a value to this variable sets {\tt stype} to 3 (Fermi-Dirac) and
1825 ! the smearing width to $k_{\rm B}T$.
1826 !
1827 ! \block{tforce}{
1828 ! {\tt tforce} & set to {\tt .true.} if the force should be calculated at the
1829 ! end of the self-consistent cycle & logical & {\tt .false.}}
1830 ! This variable is automatically set to {\tt .true.} when performing geometry
1831 ! optimisation.
1832 !
1833 ! \block{tefvit}{
1834 ! {\tt tefvit} & set to {\tt .true.} if the first-variational eigenvalue
1835 ! equation should be solved iteratively & logical & {\tt .false.}}
1836 !
1837 ! \block{tefvr}{
1838 ! {\tt tefvr} & set to {\tt .true.} if a real symmetric eigenvalue solver
1839 ! should be used for crystals which have inversion symmetry & logical &
1840 ! {\tt .true.}}
1841 ! For crystals with inversion symmetry, the first-variational Hamiltonian and
1842 ! overlap matrices can be made real by using appropriate transformations. In
1843 ! this case, a real symmetric (instead of complex Hermitian) eigenvalue solver
1844 ! can be used. This makes the calculation about three times faster.
1845 !
1846 ! \block{tm3fix}{
1847 ! {\tt ntmfix} & number of tensor moments (TM) to be fixed & integer & 0 \\
1848 ! \hline
1849 ! {\tt is(i)} & species number for entry $i$ & integer & - \\
1850 ! {\tt ia(i)} & atom number & integer & - \\
1851 ! {\tt l(i)} & $l$ of TM & integer & - \\
1852 ! \hline
1853 ! {\tt (k, p, r, t)(i)} & indices for the 3-index TM & integer & - \\
1854 ! \hline
1855 ! {\tt wkpr(t)(i)} & real TM value & real & - }
1856 ! This block sets up the fixed tensor moment (FTM). There should be as many
1857 ! TM entries as {\tt ntmfix}. See the routine {\tt tm3todm} for the tensor
1858 ! moment indexing convention.
1859 !
1860 ! \block{tmwrite}{
1861 ! {\tt tmwrite} & set to {\tt .true.} if the tensor moments and the
1862 ! corresponding decomposition of DFT+$U$ energy should be calculated
1863 ! at every loop of the self-consistent cycle & logical & {\tt .false.}}
1864 ! This variable is useful to check the convergence of the tensor moments in
1865 ! DFT+$U$ calculations. Alternatively, with {\tt task} equal to 400, one can
1866 ! calculate the tensor moments and corresponding DFT+$U$ energy contributions
1867 ! from a given density matrix and set of Slater parameters at the end of the
1868 ! self-consistent cycle.
1869 !
1870 ! \block{trdbfcr}{
1871 ! {\tt trdbfcr} & set to {\tt .true.} if the ultra long-range, real-space
1872 ! external magnetic field in Cartesian coordinates should be read in from
1873 ! {\tt BFCR.OUT} & logical & {\tt .false.}}
1874 !
1875 ! \block{trdvclr}{
1876 ! {\tt trdvclr} & set to {\tt .true.} if the ultra long-range, real-space
1877 ! external Coulomb potential should be read in from {\tt VCLR.OUT} &
1878 ! logical & {\tt .false.}}
1879 !
1880 ! \block{tsediag}{
1881 ! {\tt tsediag} & set to {\tt .true.} if the self-energy matrix should be
1882 ! treated as diagonal & logical & {\tt .false.}}
1883 ! When this variable is {\tt .true.}, the self-energy used in a $GW$
1884 ! calculation $\Sigma_{ij}({\bf k},\omega)$ is taken to be diagonal in the
1885 ! Kohn-Sham state indices $i$ and $j$. When {\tt tsediag} is {\tt .false.},
1886 ! the entire matrix is used.
1887 !
1888 ! \block{tshift}{
1889 ! {\tt tshift} & set to {\tt .true.} if the crystal can be shifted so that the
1890 ! atom closest to the origin is exactly at the origin &
1891 ! logical & {\tt .true.}}
1892 !
1893 ! \block{tstime}{
1894 ! {\tt tstime} & total simulation time of time evolution run & real &
1895 ! $1000.0$}
1896 ! See also {\tt dtimes}.
1897 !
1898 ! \block{vhmat}{
1899 ! {\tt vhmat(1)} & matrix row 1 & real(3) & $(1.0,0.0,0.0)$ \\
1900 ! \hline
1901 ! {\tt vhmat(2)} & matrix row 2 & real(3) & $(0.0,1.0,0.0)$ \\
1902 ! \hline
1903 ! {\tt vhmat(3)} & matrix row 3 & real(3) & $(0.0,0.0,1.0)$}
1904 ! This is the transformation matrix $M$ applied to every vector $\bf H$ in the
1905 ! structure factor output files {\tt SFACRHO.OUT} and {\tt SFACMAG.OUT}. It is
1906 ! stored in the usual row-column setting and applied directly as
1907 ! ${\bf H}'=M{\bf H}$ to every vector but {\em only} when writing the output
1908 ! files. See also {\tt hmaxvr} and {\tt reduceh}.
1909 !
1910 ! \block{vhighq}{
1911 ! {\tt vhighq} & {\tt .true.} if a very high-quality parameter set should be
1912 ! used & logical & {\tt .false.}}
1913 ! Setting this to {\tt .true.} results in some default parameters being
1914 ! changed to ensure excellent convergence in most situations. See also
1915 ! {\tt highq}.
1916 !
1917 ! \block{vklem}{
1918 ! {\tt vklem} & the $k$-point in lattice coordinates at which to compute the
1919 ! effective mass tensors & real(3) & $(0.0,0.0,0.0)$}
1920 ! See {\tt deltaem} and {\tt ndspem}.
1921 !
1922 ! \block{vkloff}{
1923 ! {\tt vkloff } & the $k$-point offset vector in lattice coordinates &
1924 ! real(3) & $(0.0,0.0,0.0)$}
1925 ! See {\tt ngridk}.
1926 !
1927 ! \block{vqlss}{
1928 ! {\tt vqlss} & the ${\bf q}$-vector of the spin-spiral state in lattice
1929 ! coordinates & real(3) & $(0.0,0.0,0.0)$}
1930 ! Spin-spirals arise from spinor states assumed to be of the form
1931 ! $$ \Psi^{\bf q}_{\bf k}({\bf r})=
1932 ! \left( \begin{array}{c}
1933 ! U^{{\bf q}\uparrow}_{\bf k}({\bf r})e^{i({\bf k+q/2})\cdot{\bf r}} \\
1934 ! U^{{\bf q}\downarrow}_{\bf k}({\bf r})e^{i({\bf k-q/2})\cdot{\bf r}} \\
1935 ! \end{array} \right). $$
1936 ! These are determined using a second-variational approach, and give rise to a
1937 ! magnetisation density of the form
1938 ! $$ {\bf m}^{\bf q}({\bf r})=(m_x({\bf r})\cos({\bf q \cdot r}),
1939 ! m_y({\bf r})\sin({\bf q \cdot r}),m_z({\bf r})), $$
1940 ! where $m_x$, $m_y$ and $m_z$ are lattice periodic. See also {\tt spinsprl}.
1941 !
1942 ! \block{wmaxgw}{
1943 ! {\tt wmaxgw} & maximum Matsubara frequency for $GW$ calculations & real &
1944 ! $-5.0$}
1945 ! This defines the cut-off of the Matsubara frequencies on the imaginary
1946 ! axis for calculating the $GW$ self-energy and solving the Dyson equation.
1947 ! If this number is negative then the cut-off is taken to be
1948 ! $|{\tt wmaxgw}|\times\Delta\epsilon$, where $\Delta\epsilon$ is the
1949 ! difference between the largest and smallest Kohn-Sham valence eigenvalues.
1950 !
1951 ! \block{wplot}{
1952 ! {\tt nwplot} & number of frequency/energy points in the DOS or optics plot &
1953 ! integer & $500$ \\
1954 ! {\tt ngrkf} & fine $k$-point grid size used for integrating functions in the
1955 ! Brillouin zone & integer & $100$ \\
1956 ! {\tt nswplot} & level of smoothing applied to DOS/optics output & integer &
1957 ! $1$ \\
1958 ! \hline
1959 ! {\tt wplot} & frequency/energy window for the DOS or optics plot & real(2) &
1960 ! $(-0.5,0.5)$}
1961 ! DOS and optics plots require integrals of the kind
1962 ! $$ g(\omega_i)=\frac{\Omega}{(2\pi)^3}\int_{\rm BZ} f({\bf k})
1963 ! \delta(\omega_i-e({\bf k}))d{\bf k}. $$
1964 ! These are calculated by first interpolating the functions $e({\bf k})$ and
1965 ! $f({\bf k})$ with the trilinear method on a much finer mesh whose size is
1966 ! determined by {\tt ngrkf}. Then the $\omega$-dependent histogram of the
1967 ! integrand is accumulated over the fine mesh. If the output function is noisy
1968 ! then either {\tt ngrkf} should be increased or {\tt nwplot} decreased.
1969 ! Alternatively, the output function can be artificially smoothed up to a
1970 ! level given by {\tt nswplot}. This is the number of successive 3-point
1971 ! averages to be applied to the function $g$.
1972 !
1973 ! \block{wsfac}{
1974 ! {\tt wsfac} & energy window to be used when calculating density or magnetic
1975 ! structure factors & real(2) & $(-10^6,10^6)$}
1976 ! Only those states with eigenvalues within this window will contribute to the
1977 ! density or magnetisation. See also {\tt hmaxvr} and {\tt vhmat}.
1978 !
1979 ! \block{xctype}{
1980 ! {\tt xctype} & integers defining the type of exchange-correlation functional
1981 ! to be used & integer(3) & $(3,0,0)$}
1982 ! Normally only the first value is used to define the functional type. The
1983 ! other value may be used for external libraries. Currently implemented are:
1984 ! \vskip 6pt
1985 ! \begin{tabularx}{\textwidth}[h]{lX}
1986 ! $-n$ & Exact-exchange optimised effective potential (EXX-OEP) method with
1987 ! correlation energy and potential given by functional number $n$ \\
1988 ! 1 & No exchange-correlation funtional ($E_{\rm xc}\equiv 0$) \\
1989 ! 2 & LDA, Perdew-Zunger/Ceperley-Alder, {\it Phys. Rev. B} {\bf 23}, 5048
1990 ! (1981) \\
1991 ! 3 & LSDA, Perdew-Wang/Ceperley-Alder, {\it Phys. Rev. B} {\bf 45}, 13244
1992 ! (1992) \\
1993 ! 4 & LDA, X-alpha approximation, J. C. Slater, {\it Phys. Rev.} {\bf 81}, 385
1994 ! (1951) \\
1995 ! 5 & LSDA, von Barth-Hedin, {\it J. Phys. C} {\bf 5}, 1629 (1972) \\
1996 ! 20 & GGA, Perdew-Burke-Ernzerhof, {\it Phys. Rev. Lett.} {\bf 77}, 3865
1997 ! (1996) \\
1998 ! 21 & GGA, Revised PBE, Zhang-Yang, {\it Phys. Rev. Lett.} {\bf 80}, 890
1999 ! (1998) \\
2000 ! 22 & GGA, PBEsol, Phys. Rev. Lett. 100, 136406 (2008) \\
2001 ! 26 & GGA, Wu-Cohen exchange (WC06) with PBE correlation, {\it Phys. Rev. B}
2002 ! {\bf 73}, 235116 (2006) \\
2003 ! 30 & GGA, Armiento-Mattsson (AM05) spin-unpolarised functional,
2004 ! {\it Phys. Rev. B} {\bf 72}, 085108 (2005) \\
2005 ! 100 & Libxc functionals; the second and third values of {\tt xctype} define
2006 ! the exchange and correlation functionals in the Libxc library,
2007 ! respectively
2008 ! \end{tabularx}
2009 !
2010 ! \section{Contributing to Elk}
2011 ! Please bear in mind when writing code for the Elk project that it should be
2012 ! an exercise in physics and not software engineering. All code should
2013 ! therefore be kept as simple and concise as possible, and above all it should
2014 ! be easy for anyone to locate and follow the Fortran representation of the
2015 ! original mathematics. We would also appreciate the following conventions
2016 ! being adhered to:
2017 ! \begin{itemize}
2018 ! \item Strict Fortran 2008 should be used. Features which are marked as
2019 ! obsolescent in Fortran 2008 should be avoided. These include assigned
2020 ! format specifiers, labeled do-loops, computed goto statements and statement
2021 ! functions.
2022 ! \item Modules should be used in place of common blocks for declaring
2023 ! global variables. Use the existing modules to declare new global variables.
2024 ! \item Any code should be written in lower-case free form style, starting
2025 ! from column one. Try and keep the length of each line to fewer than 80
2026 ! characters using the \& character for line continuation.
2027 ! \item Every function or subroutine, no matter how small, should be in its
2028 ! own file named {\tt routine.f90}, where {\tt routine} is the function or
2029 ! subroutine name. It is recommended that the routines are named so as to
2030 ! make their purpose apparent from the name alone.
2031 ! \item Use of {\tt implicit none} is mandatory. Remember also to define the
2032 ! {\tt intent} of any passed arguments.
2033 ! \item Local allocatable arrays must be deallocated on exit of the routine to
2034 ! prevent memory leakage. Use of automatic arrays should be limited to arrays
2035 ! of small size.
2036 ! \item Every function or subroutine must be documented with the Protex source
2037 ! code documentation system. This should include a short \LaTeX\ description
2038 ! of the algorithms and methods involved. Equations which need to be
2039 ! referenced should be labeled with {\tt routine\_1}, {\tt routine\_2}, etc.
2040 ! The authorship of each new piece of code or modification should be
2041 ! indicated in the {\tt REVISION HISTORY} part of the header. See the Protex
2042 ! documentation for details.
2043 ! \item Ensure as much as possible that a routine will terminate the program
2044 ! when given improper input instead of continuing with erroneous results.
2045 ! Specifically, functions should have a well-defined domain for which they
2046 ! return accurate results. Input outside that domain should result in an
2047 ! error message and termination.
2048 ! \item Report errors prior to termination with a short description, for
2049 ! example:
2050 ! \begin{verbatim}
2051 ! write(*,*)
2052 ! write(*,'("Error(readinput): natoms < 1 : ",I8)') natoms(is)
2053 ! write(*,'(" for species ",I4)') is
2054 ! write(*,*)
2055 ! stop
2056 ! \end{verbatim}
2057 ! \item Wherever possible, real numbers outputted as ASCII data should be
2058 ! formatted with the {\tt G18.10} specifier.
2059 ! \item Avoid redundant or repeated code: check to see if the routine you need
2060 ! already exists, before writing a new one.
2061 ! \item All reading in of ASCII data should be done in the subroutine
2062 ! {\tt readinput}. For binary data, separate routines for reading and writing
2063 ! should be used (for example {\tt writestate} and {\tt readstate}).
2064 ! \item Input filenames should be in lowercase and have the extension
2065 ! {\tt .in} . All output filenames should be in uppercase with the extension
2066 ! {\tt .OUT} .
2067 ! \item All internal units should be atomic. Input and output units should be
2068 ! atomic by default and clearly stated otherwise. Rydbergs should not be used
2069 ! under any circumstances.
2070 ! \end{itemize}
2071 ! \subsection{Licensing}
2072 ! Routines which constitute the main part of the code are released under the
2073 ! GNU General Public License (GPL). Library routines are released under the
2074 ! less restrictive GNU Lesser General Public License (LGPL). Both licenses
2075 ! are contained in the file {\tt COPYING}. Any contribution to the code must
2076 ! be licensed at the authors' discretion under either the GPL or LGPL.
2077 ! Author(s) of the code retain the copyrights. Copyright and (L)GPL
2078 ! information must be included at the beginning of every file, and no code
2079 ! will be accepted without this.
2080 !
2081 !EOI
2082 
subroutine writepmat
Definition: writepmat.f90:10
subroutine wxcplot
Definition: wxcplot.f90:7
subroutine gwsefm
Definition: gwsefm.f90:7
logical batch
Definition: modvars.f90:11
integer maxthd1
Definition: modomp.f90:13
subroutine gndstate
Definition: gndstate.f90:10
subroutine dielectric_tdrt
subroutine maguplot
Definition: maguplot.f90:7
subroutine gndstulr
Definition: gndstulr.f90:7
subroutine checkrst
Definition: checkrst.f90:7
integer task
Definition: modmain.f90:1299
logical mp_mpi
Definition: modmpi.f90:17
integer maxthd
Definition: modomp.f90:11
subroutine writeefieldw
Definition: writeefieldw.f90:7
subroutine writestrain
Definition: writestrain.f90:7
subroutine phononsc
Definition: phononsc.f90:7
subroutine writew90
Definition: writew90.f90:7
program elk
Definition: elk.f90:7
logical ramdisk
Definition: modramdisk.f90:9
subroutine ephdos
Definition: ephdos.f90:7
Definition: modomp.f90:6
subroutine ephcouple
Definition: ephcouple.f90:7
subroutine dbxcplot
Definition: dbxcplot.f90:7
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
subroutine moldyn
Definition: moldyn.f90:7
subroutine writeemd
Definition: writeemd.f90:7
subroutine tddft
Definition: tddft.f90:7
subroutine writeevbse
Definition: writeevbse.f90:7
integer maxthdmkl
Definition: modomp.f90:15
subroutine bandstr
Definition: bandstr.f90:10
subroutine phonon
Definition: phonon.f90:7
subroutine effmass
Definition: effmass.f90:7
subroutine phdos
Definition: phdos.f90:7
subroutine tddftlr
Definition: tddftlr.f90:7
subroutine spiralsc
Definition: spiralsc.f90:7
subroutine writephn
Definition: writephn.f90:7
integer ntasks
Definition: modmain.f90:1293
subroutine rdstatus
Definition: modramdisk.f90:297
subroutine batchdv
Definition: batchdv.f90:7
integer np_mpi
Definition: modmpi.f90:13
subroutine nonlinopt
Definition: nonlinopt.f90:10
subroutine hartfock
Definition: hartfock.f90:7
subroutine nesting
Definition: nesting.f90:7
subroutine omp_init
Definition: modomp.f90:32
subroutine tdhfc
Definition: tdhfc.f90:7
subroutine gwdmat
Definition: gwdmat.f90:7
integer maxlvl
Definition: modomp.f90:17
integer itask
Definition: modmain.f90:1295
subroutine dielectric
Definition: dielectric.f90:10
subroutine vecplot
Definition: vecplot.f90:10
subroutine writeevsp
Definition: writeevsp.f90:7
subroutine writehmlbse
Definition: writehmlbse.f90:7
subroutine writesf
Definition: writesf.f90:7
subroutine piezoelt
Definition: piezoelt.f90:7
subroutine wfcrplot
Definition: wfcrplot.f90:7
subroutine emdplot
Definition: emdplot.f90:7
subroutine writetm
Definition: writetm.f90:7
subroutine mossbauer
Definition: mossbauer.f90:10
logical wrtdsk
Definition: modramdisk.f90:15
subroutine geomopt
Definition: geomopt.f90:7
subroutine fermisurf
Definition: fermisurf.f90:7
subroutine elnes
Definition: elnes.f90:7
subroutine bandstrulr
Definition: bandstrulr.f90:7
subroutine potplot
Definition: potplot.f90:10
subroutine writeepsinv
Definition: writeepsinv.f90:7
subroutine gndsteph
Definition: gndsteph.f90:7
subroutine writewfpw
Definition: writewfpw.f90:7
subroutine mae
Definition: mae.f90:7
subroutine writelsj
Definition: writelsj.f90:7
subroutine elfplot
Definition: elfplot.f90:10
logical wrtvars
Definition: modvars.f90:9
subroutine delfile(fname)
Definition: moddelf.f90:15
subroutine writedosu
Definition: writedosu.f90:7
subroutine dielectric_bse
Definition: modmpi.f90:6
subroutine rhoplot
Definition: rhoplot.f90:10
subroutine wfplot
Definition: wfplot.f90:7
subroutine writeefg
Definition: writeefg.f90:10
subroutine writebox(fnum, str)
Definition: writebox.f90:7
logical trestart
Definition: modmain.f90:1057
subroutine rhosplot
Definition: rhosplot.f90:7
subroutine gwbandstr
Definition: gwbandstr.f90:7
subroutine moke
Definition: moke.f90:7
subroutine mkl_init
Definition: mkl_init.f90:7
subroutine omp_reset
Definition: modomp.f90:71
integer, dimension(3), parameter version
Definition: modmain.f90:1289
integer lp_mpi
Definition: modmpi.f90:15
subroutine rdmft
Definition: rdmft.f90:10
subroutine sfacmag
Definition: sfacmag.f90:10
subroutine phdisp
Definition: phdisp.f90:7
subroutine initrd
Definition: modramdisk.f90:37
subroutine readinput
Definition: readinput.f90:10
subroutine writeafpdt
Definition: writeafpdt.f90:7
subroutine bornechg
Definition: bornechg.f90:7
subroutine writegwefm
Definition: writegwefm.f90:7
subroutine alpha2f
Definition: alpha2f.f90:7
subroutine geomplot
Definition: geomplot.f90:7
subroutine eliashberg
Definition: eliashberg.f90:10
subroutine genafieldt
Definition: genafieldt.f90:10
subroutine bornecdyn
Definition: bornecdyn.f90:7
subroutine potuplot
Definition: potuplot.f90:7
subroutine writedos
Definition: writedos.f90:7
subroutine writeexpmat
Definition: writeexpmat.f90:7
subroutine aceplot
Definition: aceplot.f90:7
subroutine gwspecf
Definition: gwspecf.f90:7
integer, dimension(maxtasks) tasks
Definition: modmain.f90:1297
subroutine phlwidth
Definition: phlwidth.f90:7
subroutine tddftsplr
Definition: tddftsplr.f90:7
subroutine rhouplot
Definition: rhouplot.f90:7
subroutine writestress
Definition: writestress.f90:7
subroutine testcheck
Definition: testcheck.f90:7
subroutine sfacrho
Definition: sfacrho.f90:10
integer mpicom
Definition: modmpi.f90:11
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition: modvars.f90:16
subroutine fermisurfbxsf
subroutine magnetoelt
Definition: magnetoelt.f90:7
subroutine torque
Definition: torque.f90:7
integer ierror
Definition: modmpi.f90:19
subroutine jprplot
Definition: jprplot.f90:7