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