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