The Elk Code
writeinfo.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2009 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 !BOP
7 ! !ROUTINE: writeinfo
8 ! !INTERFACE:
9 subroutine writeinfo(fnum)
10 ! !USES:
11 use modmain
12 use moddftu
13 use modrdm
14 use modgw
15 use modxcifc
16 use modmpi
17 ! !INPUT/OUTPUT PARAMETERS:
18 ! fnum : unit specifier for INFO.OUT file (in,integer)
19 ! !DESCRIPTION:
20 ! Outputs basic information about the run to the file {\tt INFO.OUT}. Does not
21 ! close the file afterwards.
22 !
23 ! !REVISION HISTORY:
24 ! Created January 2003 (JKD)
25 ! Updated with DFT+U quantities July 2009 (FC)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer fnum
31 ! local variables
32 integer is,ia,k,l,i
33 real(8) t1
34 character(10) dat,tim
35 character(64) str
36 write(str,'("Elk code version ",I0,".",I0,".",I0)') version
37 call writebox(fnum,trim(str))
38 call date_and_time(date=dat,time=tim)
39 write(fnum,*)
40 write(fnum,'("Date (YYYY-MM-DD) : ",A4,"-",A2,"-",A2)') dat(1:4),dat(5:6), &
41  dat(7:8)
42 write(fnum,'("Time (hh:mm:ss) : ",A2,":",A2,":",A2)') tim(1:2),tim(3:4), &
43  tim(5:6)
44 if (np_mpi > 1) then
45  write(fnum,*)
46  write(fnum,'("Using MPI, number of processes : ",I0)') np_mpi
47 end if
48 if (notelns > 0) then
49  write(fnum,*)
50  write(fnum,'("Notes :")')
51  do i=1,notelns
52  write(fnum,'(A)') trim(notes(i))
53  end do
54 end if
55 write(fnum,*)
56 write(fnum,'("All units are atomic (Hartree, Bohr, etc.)")')
57 select case(task)
58 case(2,3)
59  if (trdstate) then
60  call writebox(fnum,"Geometry optimisation run resuming from STATE.OUT")
61  else
62  call writebox(fnum,"Geometry optimisation starting from atomic densities")
63  end if
64 case(5)
65  call writebox(fnum,"Ground-state Hartree-Fock run")
66 case(300)
67  call writebox(fnum,"Reduced density matrix functional theory run")
68 case default
69  if (trdstate) then
70  call writebox(fnum,"Ground-state run resuming from STATE.OUT")
71  else
72  call writebox(fnum,"Ground-state run starting from atomic densities")
73  end if
74 end select
75 write(fnum,*)
76 write(fnum,'("Lattice vectors :")')
77 write(fnum,'(3G18.10)') avec(1,1),avec(2,1),avec(3,1)
78 write(fnum,'(3G18.10)') avec(1,2),avec(2,2),avec(3,2)
79 write(fnum,'(3G18.10)') avec(1,3),avec(2,3),avec(3,3)
80 write(fnum,*)
81 write(fnum,'("Reciprocal lattice vectors :")')
82 write(fnum,'(3G18.10)') bvec(1,1),bvec(2,1),bvec(3,1)
83 write(fnum,'(3G18.10)') bvec(1,2),bvec(2,2),bvec(3,2)
84 write(fnum,'(3G18.10)') bvec(1,3),bvec(2,3),bvec(3,3)
85 write(fnum,*)
86 write(fnum,'("Unit cell volume : ",G18.10)') omega
87 write(fnum,'("Brillouin zone volume : ",G18.10)') omegabz
88 write(fnum,*)
89 write(fnum,'("Muffin-tin inner radius fraction : ",G18.10)') fracinr
90 write(fnum,*)
91 if (ptnucl) then
92  write(fnum,'("Nuclei treated as point charges")')
93 else
94  write(fnum,'("Nuclei treated as charged spheres")')
95 end if
96 do is=1,nspecies
97  write(fnum,*)
98  write(fnum,'("Species : ",I0," (",A,")")') is,trim(spsymb(is))
99  write(fnum,'(" parameters loaded from : ",A)') trim(spfname(is))
100  write(fnum,'(" name : ",A)') trim(spname(is))
101  write(fnum,'(" nuclear charge : ",G18.10)') spzn(is)
102  write(fnum,'(" electronic charge : ",G18.10)') spze(is)
103  write(fnum,'(" atomic mass : ",G18.10)') spmass(is)
104  write(fnum,'(" muffin-tin radius : ",G18.10)') rmt(is)
105  write(fnum,'(" number of radial points in muffin-tin : ",I6)') nrmt(is)
106  write(fnum,'(" number on inner part of muffin-tin : ",I6)') nrmti(is)
107  write(fnum,'(" approximate nuclear radius : ",G18.10)') rnucl(is)
108  write(fnum,'(" number of mesh points to nuclear radius : ",I6)') nrnucl(is)
109  write(fnum,'(" atomic positions (lattice), magnetic fields (Cartesian) :")')
110  do ia=1,natoms(is)
111  write(fnum,'(" ",I0," : ",3F12.8," ",3F12.8)') ia,atposl(:,ia,is), &
112  bfcmt(:,ia,is)
113  end do
114 end do
115 write(fnum,*)
116 write(fnum,'("Total number of atoms per unit cell : ",I4)') natmtot
117 write(fnum,'("Total muffin-tin volume : ",G18.10)') omegamt
118 write(fnum,'(" ratio of muffin-tin to unit cell volume : ",G18.10)') &
119  omegamt/omega
120 write(fnum,*)
121 write(fnum,'("Spin treatment :")')
122 if (spinpol) then
123  write(fnum,'(" spin-polarised")')
124 else
125  write(fnum,'(" spin-unpolarised")')
126 end if
127 if (spinorb) then
128  write(fnum,'(" spin-orbit coupling")')
129 end if
130 if (spincore) then
131  write(fnum,'(" spin-polarised core states")')
132 end if
133 if (spinpol) then
134  write(fnum,'(" global magnetic field (Cartesian) : ",3G18.10)') bfieldc
135  if (ncmag) then
136  write(fnum,'(" non-collinear magnetisation")')
137  else
138  write(fnum,'(" collinear magnetisation in z-direction")')
139  end if
140 end if
141 if (tbdip) then
142  write(fnum,'(" magnetic dipole field included")')
143  if (tjr) then
144  write(fnum,'(" spin and current contribution")')
145  else
146  write(fnum,'(" only spin contribution")')
147  end if
148 end if
149 if (spinsprl) then
150  write(fnum,'(" spin-spiral state assumed")')
151  write(fnum,'(" q-vector (lattice) : ",3G18.10)') vqlss
152  write(fnum,'(" q-vector (Cartesian) : ",3G18.10)') vqcss
153  write(fnum,'(" q-vector length : ",G18.10)') sqrt(vqcss(1)**2 &
154  +vqcss(2)**2+vqcss(3)**2)
155 end if
156 if (fsmtype /= 0) then
157  write(fnum,'(" fixed spin moment (FSM) calculation, type : ",I4)') fsmtype
158  if (fsmtype < 0) then
159  write(fnum,'(" only moment direction is fixed")')
160  end if
161 end if
162 if ((abs(fsmtype) == 1).or.(abs(fsmtype) == 3)) then
163  write(fnum,'(" fixing total moment to (Cartesian) :")')
164  write(fnum,'(" ",3G18.10)') momfix
165 end if
166 if ((abs(fsmtype) == 2).or.(abs(fsmtype) == 3)) then
167  write(fnum,'(" fixing local muffin-tin moments to (Cartesian) :")')
168  do is=1,nspecies
169  write(fnum,'(" species : ",I0," (",A,")")') is,trim(spsymb(is))
170  do ia=1,natoms(is)
171  write(fnum,'(" ",I0,3G18.10)') ia,mommtfix(:,ia,is)
172  end do
173  end do
174 end if
175 if (tssxc) then
176  write(fnum,'(" scaled spin exchange-correlation enabled")')
177  write(fnum,'(" scaling factor : ",G18.10)') sxcscf
178 end if
179 if (ftmtype /= 0) then
180  write(fnum,*)
181  write(fnum,'(" fixed tensor moment (FTM) calculation, type : ",I4)') ftmtype
182 end if
183 if (tefield) then
184  write(fnum,*)
185  write(fnum,'("Constant electric field applied across unit cell")')
186  write(fnum,'(" field strength : ",3G18.10)') efieldc
187  t1=norm2(efieldc(1:3))
188  write(fnum,'(" magnitude : ",G18.10)') t1
189  write(fnum,'(" volts/nanometer : ",G18.10)') t1*ef_si/1.d9
190  write(fnum,'(" maximum distance from center over which E-field is &
191  &applied : ",G18.10)') dmaxefc
192  write(fnum,'(" potential at maximum distance : ",G18.10)') vmaxefc
193 end if
194 if (tafield) then
195  write(fnum,*)
196  write(fnum,'("Constant A-field applied across unit cell")')
197  write(fnum,'(" field strength : ",3G18.10)') afieldc
198 end if
199 write(fnum,*)
200 write(fnum,'("Number of Bravais lattice symmetries : ",I4)') nsymlat
201 write(fnum,'("Number of crystal symmetries : ",I4)') nsymcrys
202 if (tsyminv) then
203  write(fnum,'("Crystal has inversion symmetry")')
204 else
205  write(fnum,'("Crystal has no inversion symmetry")')
206 end if
207 if (tefvr) then
208  write(fnum,'("Real symmetric eigensolver will be used")')
209 else
210  write(fnum,'("Complex Hermitian eigensolver will be used")')
211 end if
212 write(fnum,*)
213 if (autokpt) then
214  write(fnum,'("Radius of sphere used to determine k-point grid density : ",&
215  &G18.10)') radkpt
216 end if
217 write(fnum,'("k-point grid : ",3I6)') ngridk
218 write(fnum,'("k-point offset : ",3G18.10)') vkloff
219 if (reducek == 0) then
220  write(fnum,'("k-point set is not reduced")')
221 else if (reducek == 1) then
222  write(fnum,'("k-point set is reduced with full crystal symmetry group")')
223 else if (reducek == 2) then
224  write(fnum,'("k-point set is reduced with symmorphic symmetries only")')
225 else
226  write(*,*)
227  write(*,'("Error(writeinfo): undefined k-point reduction type : ",I0)') &
228  reducek
229  write(*,*)
230  stop
231 end if
232 write(fnum,'("Total number of k-points : ",I8)') nkpt
233 write(fnum,*)
234 write(fnum,'("Muffin-tin radius times maximum |G+k| : ",G18.10)') rgkmax
235 select case(isgkmax)
236 case(:-4)
237  write(fnum,'(" using largest radius")')
238 case(-3)
239  write(fnum,'(" using smallest radius")')
240 case(-2)
241  write(fnum,'(" using gkmax = rgkmax / 2")')
242 case(-1)
243  write(fnum,'(" using average radius")')
244 case(1:)
245  if (isgkmax <= nspecies) then
246  write(fnum,'(" using radius of species ",I0," (",A,")")') isgkmax, &
247  trim(spsymb(isgkmax))
248  end if
249 end select
250 write(fnum,'("Maximum |G+k| for APW functions : ",G18.10)') gkmax
251 write(fnum,'("Maximum (1/2)|G+k|² : ",G18.10)') 0.5d0*gkmax**2
252 write(fnum,'("Maximum |G| for potential and density : ",G18.10)') gmaxvr
253 write(fnum,'("Constant for pseudocharge density : ",I4)') npsd
254 write(fnum,'("Radial integration step length : ",I4)') lradstp
255 write(fnum,*)
256 write(fnum,'("G-vector grid sizes : ",3I6)') ngridg(:)
257 write(fnum,'("Number of G-vectors : ",I8)') ngvec
258 write(fnum,*)
259 write(fnum,'("Maximum angular momentum used for")')
260 write(fnum,'(" APW functions : ",I4)') lmaxapw
261 write(fnum,'(" outer part of muffin-tin : ",I4)') lmaxo
262 write(fnum,'(" inner part of muffin-tin : ",I4)') lmaxi
263 write(fnum,*)
264 write(fnum,'("Total nuclear charge : ",G18.10)') chgzn
265 write(fnum,'("Total core charge : ",G18.10)') chgcrtot
266 write(fnum,'("Total valence charge : ",G18.10)') chgval
267 write(fnum,'("Total excess charge : ",G18.10)') chgexs
268 write(fnum,'("Total electronic charge : ",G18.10)') chgtot
269 write(fnum,*)
270 write(fnum,'("Effective Wigner radius, rₛ : ",G18.10)') rwigner
271 write(fnum,*)
272 write(fnum,'("Number of empty states : ",I6)') nempty
273 write(fnum,'("Total number of valence states : ",I6)') nstsv
274 write(fnum,'("Total number of core states : ",I6)') nstcr
275 write(fnum,*)
276 if (lorbcnd) then
277  write(fnum,'("Conduction state local-orbitals added automatically")')
278 end if
279 write(fnum,'("Total number of local-orbitals : ",I6)') nlotot
280 if (tefvs) then
281  write(fnum,*)
282  write(fnum,'("First-variational eigenvalue equation solved in subspace")')
283  write(fnum,'(" subspace size parameter, mefvs : ",I6)') mefvs
284 end if
285 write(fnum,*)
286 if (task == 5) then
287  write(fnum,'("Hartree-Fock calculation using Kohn-Sham states")')
288  if (hybrid) then
289  write(fnum,'(" hybrid functional, coefficient : ",G18.10)') hybridc
290  end if
291 end if
292 if (xctype(1) == 100) then
293  write(fnum,'("Using Libxc version ",I0,".",I0,".",I0)') libxcv(:)
294 end if
295 if (xctype(1) < 0) then
296  write(fnum,'("Optimised effective potential (OEP) and exact exchange (EXX)")')
297  write(fnum,'(" Phys. Rev. B 53, 7024 (1996)")')
298  write(fnum,'("Correlation functional : ",3I6)') abs(xctype(1)),xctype(2:3)
299  write(fnum,'(" ",A)') trim(xcdescr)
300 else
301  write(fnum,'("Exchange-correlation functional : ",3I6)') xctype(:)
302  write(fnum,'(" ",A)') trim(xcdescr)
303  write(fnum,'(" gradient requirement : ",I4)') xcgrad
304 end if
305 if (xcgrad == 0) then
306  write(fnum,'(" local density approximation (LDA)")')
307 else if ((xcgrad == 1).or.(xcgrad == 2)) then
308  write(fnum,'(" generalised gradient approximation (GGA)")')
309 else if (any(xcgrad == [3,4,5,6])) then
310  write(fnum,'(" meta-GGA")')
311  if (xcgrad == 3) then
312  write(fnum,'(" fully deorbitalised functional")')
313  else if (xcgrad == 6) then
314  write(fnum,'(" potential-only functional")')
315  else
316  write(fnum,'(" partially deorbitalised functional")')
317  write(fnum,'(" using Kohn-Sham orbital kinetic energy density τ(r)")')
318  write(fnum,'(" kinetic energy functional used for δτ(r'')/δρ(r) : ",3I6)') &
319  ktype(:)
320  write(fnum,'(" ",A)') trim(kdescr)
321  write(fnum,'(" gradient requirement : ",I4)') kgrad
322  end if
323 end if
324 if (ksgwrho) then
325  write(fnum,*)
326  write(fnum,'("Kohn-Sham density determined via the GW Green''s function")')
327 end if
328 if (dftu /= 0) then
329  write(fnum,*)
330  write(fnum,'("DFT+U calculation")')
331  if (dftu == 1) then
332  write(fnum,'(" fully localised limit (FLL)")')
333  write(fnum,'(" see Phys. Rev. B 52, R5467 (1995)")')
334  else if (dftu == 2) then
335  write(fnum,'(" around mean field (AMF)")')
336  write(fnum,'(" see Phys. Rev. B 49, 14211 (1994)")')
337  else
338  write(*,*)
339  write(*,'("Error(writeinfo): dftu not defined : ",I0)') dftu
340  write(*,*)
341  stop
342  end if
343  do i=1,ndftu
344  is=isldu(1,i)
345  l=isldu(2,i)
346  if (inpdftu == 1) then
347  write(fnum,'(" species : ",I0," (",A,")",", l = ",I0,", U = ",F12.8, &
348  &", J = ",F12.8)') is,trim(spsymb(is)),l,ujdu(1,i),ujdu(2,i)
349  else if (inpdftu == 2) then
350  write(fnum,'(" species : ",I0," (",A,")",", l = ",I0)') is, &
351  trim(spsymb(is)),l
352  write(fnum,'(" Slater integrals are provided as input")')
353  do k=0,2*l,2
354  write(fnum,'(" F^(",I1,") = ",F12.8)') k,fdu(k,i)
355  end do
356  else if (inpdftu == 3) then
357  write(fnum,'(" species : ",I0," (",A,")",", l = ",I0)') is, &
358  trim(spsymb(is)),l
359  write(fnum,'(" Racah parameters are provided as input")')
360  do k=0,l
361  write(fnum,'(" E^(",I1,") = ",F12.8)') k,edu(k,i)
362  end do
363  else if (inpdftu == 4) then
364  write(fnum,'(" species : ",I0," (",A,")",", l = ",I0)') is, &
365  trim(spsymb(is)),l
366  write(fnum,'(" Slater integrals are calculated by means of Yukawa &
367  &potential")')
368  write(fnum,'(" Yukawa potential screening length (a.u.⁻¹) : ",F12.8)') &
369  lamdu(i)
370  else if(inpdftu == 5) then
371  write(fnum,'(" species : ",I0," (",A,")",", l = ",I0)') is, &
372  trim(spsymb(is)),l
373  write(fnum,'(" Slater integrals are calculated by means of Yukawa &
374  &potential")')
375  write(fnum,'(" Yukawa potential screening length corresponds to U = ",&
376  &F12.8)') udufix(i)
377  end if
378  end do
379 end if
380 if (task == 300) then
381  write(fnum,*)
382  write(fnum,'("RDMFT calculation")')
383  write(fnum,'(" see Phys. Rev. B 78, 201103 (2008)")')
384  write(fnum,'(" RDMFT exchange-correlation type : ",I4)') rdmxctype
385  if (rdmxctype == 1) then
386  write(fnum,'(" Hartree-Fock functional")')
387  else if (rdmxctype == 2) then
388  write(fnum,'(" Power functional, exponent : ",G18.10)') rdmalpha
389  end if
390 end if
391 write(fnum,*)
392 write(fnum,'("Smearing type : ",I4)') stype
393 write(fnum,'(" ",A)') trim(sdescr)
394 if (autoswidth) then
395  write(fnum,'("Automatic determination of smearing width")')
396 else
397  write(fnum,'("Smearing width : ",G18.10)') swidth
398  write(fnum,'("Effective electronic temperature (K) : ",G18.10)') tempk
399 end if
400 write(fnum,*)
401 write(fnum,'("Mixing type : ",I4)') mixtype
402 write(fnum,'(" ",A)') trim(mixdescr)
403 flush(fnum)
404 end subroutine
405 !EOC
406 
real(8), dimension(maxspecies) rnucl
Definition: modmain.f90:85
real(8) dmaxefc
Definition: modmain.f90:318
logical tjr
Definition: modmain.f90:618
integer mixtype
Definition: modmain.f90:693
integer, dimension(3) ktype
Definition: modmain.f90:604
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer task
Definition: modmain.f90:1291
integer, dimension(3) xctype
Definition: modmain.f90:586
real(8) rwigner
Definition: modmain.f90:734
logical spinpol
Definition: modmain.f90:228
logical hybrid
Definition: modmain.f90:1144
integer nkpt
Definition: modmain.f90:461
logical autokpt
Definition: modmain.f90:444
real(8), dimension(0:2 *lmaxdm, maxdftu) fdu
Definition: moddftu.f90:57
integer nlotot
Definition: modmain.f90:788
real(8) sxcscf
Definition: modmain.f90:666
subroutine writeinfo(fnum)
Definition: writeinfo.f90:10
real(8) omega
Definition: modmain.f90:20
logical lorbcnd
Definition: modmain.f90:833
integer nsymcrys
Definition: modmain.f90:358
integer nstcr
Definition: modmain.f90:129
logical spinsprl
Definition: modmain.f90:283
real(8) swidth
Definition: modmain.f90:889
character(256), dimension(maxspecies) spfname
Definition: modmain.f90:74
logical tssxc
Definition: modmain.f90:664
character(264) kdescr
Definition: modmain.f90:606
integer notelns
Definition: modmain.f90:1297
character(264) xcdescr
Definition: modmain.f90:588
real(8), dimension(3) vkloff
Definition: modmain.f90:450
integer lmaxo
Definition: modmain.f90:201
real(8), dimension(3) vqlss
Definition: modmain.f90:293
real(8) vmaxefc
Definition: modmain.f90:320
integer lmaxapw
Definition: modmain.f90:197
integer kgrad
Definition: modmain.f90:608
integer np_mpi
Definition: modmpi.f90:13
logical tefvr
Definition: modmain.f90:866
real(8), dimension(3) momfix
Definition: modmain.f90:253
real(8) chgcrtot
Definition: modmain.f90:716
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:49
Definition: modrdm.f90:6
logical tefield
Definition: modmain.f90:310
integer nstsv
Definition: modmain.f90:883
logical tsyminv
Definition: modmain.f90:354
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
integer ndftu
Definition: moddftu.f90:47
real(8), dimension(maxdftu) lamdu
Definition: moddftu.f90:61
real(8) radkpt
Definition: modmain.f90:446
logical tefvs
Definition: modmain.f90:869
real(8) omegamt
Definition: modmain.f90:169
real(8), dimension(0:lmaxdm, maxdftu) edu
Definition: moddftu.f90:59
integer ftmtype
Definition: moddftu.f90:79
integer nsymlat
Definition: modmain.f90:342
integer inpdftu
Definition: moddftu.f90:43
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
real(8) tempk
Definition: modmain.f90:682
real(8), dimension(maxspecies) spmass
Definition: modmain.f90:101
real(8) rgkmax
Definition: modmain.f90:493
integer, dimension(3) ngridk
Definition: modmain.f90:448
integer ngvec
Definition: modmain.f90:396
real(8) chgzn
Definition: modmain.f90:712
integer lradstp
Definition: modmain.f90:171
real(8), dimension(3) afieldc
Definition: modmain.f90:325
integer mefvs
Definition: modmain.f90:871
real(8) rdmalpha
Definition: modrdm.f90:29
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition: modmain.f90:259
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
real(8), dimension(maxspecies) spze
Definition: modmain.f90:99
real(8), dimension(3, maxatoms, maxspecies) bfcmt
Definition: modmain.f90:273
character(64) mixdescr
Definition: modmain.f90:695
real(8), dimension(3) vqcss
Definition: modmain.f90:295
Definition: modgw.f90:6
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer npsd
Definition: modmain.f90:624
real(8), dimension(maxdftu) udufix
Definition: moddftu.f90:69
integer stype
Definition: modmain.f90:885
Definition: modmpi.f90:6
character(64), dimension(maxspecies) spname
Definition: modmain.f90:76
character(64) sdescr
Definition: modmain.f90:887
integer dftu
Definition: moddftu.f90:36
real(8) chgval
Definition: modmain.f90:720
integer, dimension(maxspecies) nrnucl
Definition: modmain.f90:89
logical ptnucl
Definition: modmain.f90:83
subroutine writebox(fnum, str)
Definition: writebox.f90:7
integer rdmxctype
Definition: modrdm.f90:21
logical tbdip
Definition: modmain.f90:641
logical spinorb
Definition: modmain.f90:230
integer, dimension(3), parameter version
Definition: modmain.f90:1281
real(8), parameter ef_si
Definition: modmain.f90:1265
integer nspecies
Definition: modmain.f90:34
real(8) chgtot
Definition: modmain.f90:724
real(8) gkmax
Definition: modmain.f90:495
logical trdstate
Definition: modmain.f90:680
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer reducek
Definition: modmain.f90:455
logical autoswidth
Definition: modmain.f90:891
integer natmtot
Definition: modmain.f90:40
real(8), dimension(2, maxdftu) ujdu
Definition: moddftu.f90:51
real(8) omegabz
Definition: modmain.f90:22
logical ksgwrho
Definition: modgw.f90:38
character(256), dimension(:), allocatable notes
Definition: modmain.f90:1299
logical ncmag
Definition: modmain.f90:240
real(8) chgexs
Definition: modmain.f90:722
logical tafield
Definition: modmain.f90:322
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), dimension(3) efieldc
Definition: modmain.f90:312
real(8) hybridc
Definition: modmain.f90:1146
integer xcgrad
Definition: modmain.f90:600
logical spincore
Definition: modmain.f90:935
character(64), dimension(maxspecies) spsymb
Definition: modmain.f90:78
integer isgkmax
Definition: modmain.f90:491
real(8) fracinr
Definition: modmain.f90:209
real(8), dimension(3) bfieldc
Definition: modmain.f90:269
integer fsmtype
Definition: modmain.f90:251
integer nempty
Definition: modmain.f90:879
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
real(8) gmaxvr
Definition: modmain.f90:384
integer lmaxi
Definition: modmain.f90:205