The Elk Code
readinput.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2008 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: readinput
8 ! !INTERFACE:
9 subroutine readinput
10 ! !USES:
11 use modmain
12 use moddftu
13 use modrdm
14 use modphonon
15 use modtest
16 use modrandom
17 use modpw
18 use modtddft
19 use modulr
20 use modvars
21 use modgw
22 use modbog
23 use modw90
24 use modtdhfc
25 use modmpi
26 use modomp
27 use modramdisk
28 ! !DESCRIPTION:
29 ! Reads in the input parameters from the file {\tt elk.in}. Also sets default
30 ! values for the input parameters.
31 !
32 ! !REVISION HISTORY:
33 ! Created September 2002 (JKD)
34 !EOP
35 !BOC
36 implicit none
37 ! local variables
38 logical lv
39 integer is,ia,ias,ios
40 integer i,j,k,l
41 real(8) sc,sc1,sc2,sc3
42 real(8) scx,scy,scz
43 real(8) scu,scu1,scu2,scu3
44 real(8) solscf,zn
45 real(8) axang(4),rot(3,3)
46 real(8) rndavec,v(3),t1
47 character(256) block,symb,str
48 
49 !------------------------!
50 ! default values !
51 !------------------------!
52 ntasks=0
53 avec(:,:)=0.d0
54 avec(1,1)=1.d0
55 avec(2,2)=1.d0
56 avec(3,3)=1.d0
57 davec(:,:)=0.d0
58 sc=1.d0
59 sc1=1.d0
60 sc2=1.d0
61 sc3=1.d0
62 scx=1.d0
63 scy=1.d0
64 scz=1.d0
65 epslat=1.d-6
66 primcell=.false.
67 tshift=.true.
68 ngridk(:)=1
69 dngridk(:)=0
70 vkloff(:)=0.d0
71 autokpt=.false.
72 radkpt=40.d0
73 reducek=1
74 ngridq(:)=-1
75 reduceq=1
76 rgkmax=7.d0
77 drgkmax=0.d0
78 gmaxvr=12.d0
79 dgmaxvr=0.d0
80 lmaxapw=8
81 dlmaxapw=0
82 lmaxo=6
83 dlmaxo=0
84 lmaxi=1
85 fracinr=0.01d0
86 trhonorm=.true.
87 xctype(1)=3
88 xctype(2:3)=0
89 xctsp(1)=3
90 xctsp(2:3)=0
91 ktype(1)=52
92 ktype(2:3)=0
93 stype=3
94 swidth=0.001d0
95 autoswidth=.false.
96 mstar=10.d0
97 epsocc=1.d-10
98 epschg=1.d-3
99 nempty0=4.d0
100 dnempty0=0.d0
101 maxscl=200
102 mixtype=3
103 mixsave=.false.
104 amixpm(1)=0.05d0
105 amixpm(2)=1.d0
106 ! Broyden parameters recommended by M. Meinert
107 mixsdb=5
108 broydpm(1)=0.4d0
109 broydpm(2)=0.15d0
110 mixrho=.false.
111 epspot=1.d-6
112 epsengy=1.d-4
113 epsforce=5.d-3
114 epsstress=2.d-3
115 molecule=.false.
116 nspecies=0
117 natoms(:)=0
118 atposl(:,:,:)=0.d0
119 datposl(:,:,:)=0.d0
120 atposc(:,:,:)=0.d0
121 bfcmt0(:,:,:)=0.d0
122 sppath=''
123 scrpath=''
124 nvp1d=2
125 if (allocated(vvlp1d)) deallocate(vvlp1d)
126 allocate(vvlp1d(3,nvp1d))
127 vvlp1d(:,1)=0.d0
128 vvlp1d(:,2)=1.d0
129 npp1d=200
130 ip01d=1
131 vclp2d(:,:)=0.d0
132 vclp2d(1,1)=1.d0
133 vclp2d(2,2)=1.d0
134 np2d(:)=40
135 vclp3d(:,:)=0.d0
136 vclp3d(1,1)=1.d0
137 vclp3d(2,2)=1.d0
138 vclp3d(3,3)=1.d0
139 np3d(:)=20
140 nwplot=500
141 ngrkf=100
142 nswplot=1
143 wplot(1)=-0.5d0
144 wplot(2)=0.5d0
145 dosocc=.false.
146 tpdos=.true.
147 dosmsum=.false.
148 dosssum=.false.
149 lmirep=.true.
150 spinpol=.false.
151 spinorb=.false.
152 socscf=1.d0
153 bforb=.false.
154 bfdmag=.false.
155 atpopt=1
156 maxatpstp=200
157 tau0atp=0.2d0
158 deltast=0.005d0
159 avecref(:,:)=0.d0
160 latvopt=0
161 maxlatvstp=30
162 tau0latv=0.2d0
163 lradstp=4
164 chgexs=0.d0
165 dchgexs=0.d0
166 scissor=0.d0
167 noptcomp=1
168 ! list of all optical tensor components
169 do k=1,3; do j=1,3; do i=1,3
170  l=(k-1)*9+(j-1)*3+i
171  optcomp(:,l)=[i,j,k]
172 end do; end do; end do
173 optcomp(:,1)=1
174 intraband=.false.
175 evaltol=-1.d0
176 epsband=1.d-12
177 demaxbnd=2.5d0
178 autolinengy=.false.
179 dlefe=-0.1d0
180 autodlefe=.true.
181 deapw=0.2d0
182 delorb=0.05d0
183 bfieldc0(:)=0.d0
184 dbfieldc0(:)=0.d0
185 efieldc(:)=0.d0
186 dmaxefc=1.d6
187 afieldc(:)=0.d0
188 dafieldc(:)=0.d0
189 afspc(:,:)=0.d0
190 dafspc(:,:)=0.d0
191 fsmtype=0
192 momfix(:)=0.d0
193 momfixm=0.d0
194 dmomfix(:)=0.d0
195 mommtfix(:,:,:)=1.d6
196 mommtfixm(:,:)=-1.d0
197 taufsm=0.01d0
198 rmtdelta=0.05d0
199 isgkmax=-1
200 symtype=1
201 deltaph=0.01d0
202 nphwrt=1
203 if (allocated(vqlwrt)) deallocate(vqlwrt)
204 allocate(vqlwrt(3,nphwrt))
205 vqlwrt(:,:)=0.d0
206 notelns=0
207 tforce=.false.
208 maxitoep=400
209 tau0oep=0.1d0
210 nkstlist=1
211 kstlist(:,1)=1
212 vklem(:)=0.d0
213 deltaem=0.025d0
214 ndspem=1
215 nosource=.false.
216 spinsprl=.false.
217 ssdph=.true.
218 vqlss(:)=0.d0
219 dvqlss(:)=0.d0
220 nwrite=0
221 dftu=0
222 inpdftu=1
223 ndftu=0
224 ujdu(:,:)=0.d0
225 fdu(:,:)=0.d0
226 edu(:,:)=0.d0
227 lamdu(:)=0.d0
228 udufix(:)=0.d0
229 dudufix(:)=0.d0
230 tmwrite=.false.
231 rdmxctype=2
232 rdmmaxscl=2
233 maxitn=200
234 maxitc=0
235 taurdmn=0.5d0
236 taurdmc=0.25d0
237 rdmalpha=0.656d0
238 rdmtemp=0.d0
239 reducebf=1.d0
240 ptnucl=.true.
241 tefvr=.true.
242 tefvit=.false.
243 nefvit=2
244 befvit=0.05d0
245 vecql(:)=0.d0
246 mustar=0.15d0
247 sqaxis(1:2)=0.d0
248 sqaxis(3)=1.d0
249 test=.false.
250 spincore=.false.
251 solscf=1.d0
252 emaxelnes=-1.2d0
253 wsfac(1)=-1.1d6; wsfac(2)=1.1d6
254 vhmat(:,:)=0.d0
255 vhmat(1,1)=1.d0
256 vhmat(2,2)=1.d0
257 vhmat(3,3)=1.d0
258 reduceh=.true.
259 hybrid0=.false.
260 hybridc=1.d0
261 ecvcut=-3.5d0
262 esccut=-0.4d0
263 gmaxrf=3.d0
264 mbwgrf=-1
265 emaxrf=1.d6
266 ntemp=40
267 nvbse0=2
268 ncbse0=3
269 nvxbse=0
270 ncxbse=0
271 bsefull=.false.
272 hxbse=.true.
273 hdbse=.true.
274 fxctype=-1
275 fxclrc(1)=0.d0
276 fxclrc(2)=0.d0
277 rndatposc=0.d0
278 rndbfcmt=0.d0
279 rndavec=0.d0
280 c_tb09=0.d0
281 tc_tb09=.false.
282 hmaxvr=20.d0
283 hkmax=12.d0
284 lorbcnd=.false.
285 lorbordc=3
286 nrmtscf=1.d0
287 dnrmtscf=0.d0
288 lmaxdb=3
289 epsdev=0.0025d0
290 npmae0=-1
291 wrtvars=.false.
292 ftmtype=0
293 ntmfix=0
294 tauftm=0.1d0
295 cmagz=.false.
296 axang(:)=0.d0
297 dncgga=1.d-8
298 tstime=1000.d0
299 dtimes=0.1d0
300 npulse=0
301 nramp=0
302 nstep=0
303 ntswrite(1)=500
304 ntswrite(2)=1
305 nxoapwlo=0
306 nxlo=0
307 tdrho1d=.false.
308 tdrho2d=.false.
309 tdrho3d=.false.
310 tdmag1d=.false.
311 tdmag2d=.false.
312 tdmag3d=.false.
313 tdjr1d=.false.
314 tdjr2d=.false.
315 tdjr3d=.false.
316 tddos=.false.
317 tdlsj=.false.
318 tdjtk=.false.
319 tdxrmk=.false.
320 rndevt0=0.d0
321 sxcscf=1.d0
322 dsxcscf=0.d0
323 avecu(:,:)=0.d0
324 avecu(1,1)=1.d0
325 avecu(2,2)=1.d0
326 avecu(3,3)=1.d0
327 scu=1.d0
328 scu1=1.d0
329 scu2=1.d0
330 scu3=1.d0
331 q0cut=0.d0
332 rndbfcu=0.d0
333 bfieldcu(:)=0.d0
334 efieldcu(:)=0.d0
335 tplotq0=.true.
336 trdvclr=.false.
337 trdbfcr=.false.
338 wmaxgw=-10.d0
339 tsediag=.false.
340 actype=10
341 npole=3
342 nspade=100
343 tfav0=.true.
344 rmtscf=1.d0
345 mrmtav=0
346 rmtall=-1.d0
347 maxthd=0
348 maxthd1=0
349 maxthdmkl=0
350 maxlvl=4
351 tdphi=0.d0
352 thetamld=45.d0*pi/180.d0
353 ntsbackup=0
354 ! Wannier90 variables
355 seedname='wannier'
356 num_wann=0
357 num_bands=0
358 num_iter=500
359 dis_num_iter=500
360 trial_step=1.d-3
361 nxlwin=0
362 wrtunk=.false.
363 tbdip=.false.
364 tjr=.false.
365 tauefm=0.01d0
366 epsefm=1.d-6
367 ehfb=1.d0
368 t0gclq0=.false.
369 tafindt=.false.
370 afindpm(:)=0.d0
371 afindpm(2)=1.d0
372 nkspolar=4
373 ntsforce=100
374 wphcut=1.d-6
375 ephscf(1)=8.d0
376 ephscf(2)=0.02d0
377 anomalous=.false.
378 tephde=.false.
379 bdiag=.false.
380 ecutb=0.001d0
381 ediag=.false.
382 pwxpsn=2
383 ramdisk=.true.
384 wrtdsk=.true.
385 epsdmat=1.d-8
386 tm3vdl=.false.
387 batch=.false.
388 tafspt=.false.
389 tbaspat=.false.
390 trdatdv=.false.
391 atdfc=0.d0
392 maxforce=-1.d0
393 msmgmt=0
394 ntsorth=1000
395 deltabf=0.5d0
396 jtconst0=.false.
397 trmt0=.true.
398 ksgwrho=.false.
399 npfftg=4
400 npfftgc=4
401 npfftq=4
402 npfftw=4
403 tphnat=.false.
404 ecutthc=0.01d0
405 tbdipu=.false.
406 bdipscf=1.d0
407 
408 !--------------------------!
409 ! read from elk.in !
410 !--------------------------!
411 open(50,file='elk.in',status='OLD',form='FORMATTED',iostat=ios)
412 if (ios /= 0) then
413  write(*,*)
414  write(*,'("Error(readinput): error opening elk.in")')
415  write(*,*)
416  stop
417 end if
418 10 continue
419 read(50,*,end=30) block
420 ! check for a comment
421 if ((block(1:1) == '!').or.(block(1:1) == '#')) goto 10
422 select case(trim(block))
423 case('tasks')
424  do i=1,maxtasks
425  read(50,'(A)',err=20) str
426  if (trim(str) == '') then
427  if (i == 1) then
428  write(*,*)
429  write(*,'("Error(readinput): no tasks to perform")')
430  write(*,*)
431  stop
432  end if
433  ntasks=i-1
434  goto 10
435  end if
436  read(str,*,iostat=ios) tasks(i)
437  if (ios /= 0) then
438  write(*,*)
439  write(*,'("Error(readinput): error reading tasks")')
440  write(*,'("(blank line required after tasks block)")')
441  write(*,*)
442  stop
443  end if
444  end do
445  write(*,*)
446  write(*,'("Error(readinput): too many tasks")')
447  write(*,'("Adjust maxtasks in modmain and recompile code")')
448  write(*,*)
449  stop
450 case('species')
451 ! generate a species file
452  call genspecies(50)
453 case('fspecies')
454 ! generate fractional species files
455  do is=1,maxspecies
456  read(50,'(A)',err=20) str
457  if (trim(str) == '') goto 10
458  read(str,*,iostat=ios) zn,symb
459  if (ios /= 0) then
460  write(*,*)
461  write(*,'("Error(readinput): error reading fractional species")')
462  write(*,'("(blank line required after fspecies block)")')
463  write(*,*)
464  stop
465  end if
466  if (zn > 0.d0) then
467  write(*,*)
468  write(*,'("Error(readinput): fractional nuclear Z > 0 : ",G18.10)') zn
469  write(*,*)
470  stop
471  end if
472  call genfspecies(zn,symb)
473  end do
474  write(*,*)
475  write(*,'("Error(readinput): too many fractional nucleus species")')
476  write(*,*)
477  stop
478 case('avec')
479  do i=1,3
480  read(50,'(A)',err=20) str
481  read(str,*,err=20) avec(:,i)
482  read(str,*,iostat=ios) avec(:,i),davec(:,i)
483  end do
484 case('scale')
485  read(50,*,err=20) sc
486 case('scale1')
487  read(50,*,err=20) sc1
488 case('scale2')
489  read(50,*,err=20) sc2
490 case('scale3')
491  read(50,*,err=20) sc3
492 case('scalex')
493  read(50,*,err=20) scx
494 case('scaley')
495  read(50,*,err=20) scy
496 case('scalez')
497  read(50,*,err=20) scz
498 case('epslat')
499  read(50,*,err=20) epslat
500  if (epslat <= 0.d0) then
501  write(*,*)
502  write(*,'("Error(readinput): epslat <= 0 : ",G18.10)') epslat
503  write(*,*)
504  stop
505  end if
506 case('primcell')
507  read(50,*,err=20) primcell
508 case('tshift')
509  read(50,*,err=20) tshift
510 case('autokpt')
511  read(50,*,err=20) autokpt
512 case('radkpt')
513  read(50,*,err=20) radkpt
514  if (radkpt <= 0.d0) then
515  write(*,*)
516  write(*,'("Error(readinput): radkpt <= 0 : ",G18.10)') radkpt
517  write(*,*)
518  stop
519  end if
520 case('ngridk')
521  read(50,'(A)',err=20) str
522  read(str,*,err=20) ngridk(:)
523  read(str,*,iostat=ios) ngridk(:),dngridk(:)
524  if (any(ngridk(:) < 1)) then
525  write(*,*)
526  write(*,'("Error(readinput): invalid ngridk : ",3I8)') ngridk
527  write(*,*)
528  stop
529  end if
530  autokpt=.false.
531 case('vkloff')
532  read(50,*,err=20) vkloff(:)
533  if (any(vkloff(:) < 0.d0).or.any(vkloff(:) >= 1.d0)) then
534  write(*,*)
535  write(*,'("Error(readinput): vkloff components should be in [0,1) : ",&
536  &3G18.10)') vkloff
537  write(*,*)
538  stop
539  end if
540 case('reducek')
541  read(50,*,err=20) reducek
542 case('ngridq')
543  read(50,*,err=20) ngridq(:)
544  if (any(ngridq(:) < 1)) then
545  write(*,*)
546  write(*,'("Error(readinput): invalid ngridq : ",3I8)') ngridq
547  write(*,*)
548  stop
549  end if
550 case('reduceq')
551  read(50,*,err=20) reduceq
552 case('rgkmax')
553  read(50,'(A)',err=20) str
554  read(str,*,err=20) rgkmax
555  read(str,*,iostat=ios) rgkmax,drgkmax
556  if (rgkmax <= 0.d0) then
557  write(*,*)
558  write(*,'("Error(readinput): rgkmax <= 0 : ",G18.10)') rgkmax
559  write(*,*)
560  stop
561  end if
562 case('gmaxvr')
563  read(50,'(A)',err=20) str
564  read(str,*,err=20) gmaxvr
565  read(str,*,iostat=ios) gmaxvr,dgmaxvr
566 case('lmaxapw')
567  read(50,'(A)',err=20) str
568  read(str,*,err=20) lmaxapw
569  read(str,*,iostat=ios) lmaxapw,dlmaxapw
570  if (lmaxapw < 0) then
571  write(*,*)
572  write(*,'("Error(readinput): lmaxapw < 0 : ",I8)') lmaxapw
573  write(*,*)
574  stop
575  end if
576  if (lmaxapw >= maxlapw) then
577  write(*,*)
578  write(*,'("Error(readinput): lmaxapw too large : ",I8)') lmaxapw
579  write(*,'("Adjust maxlapw in modmain and recompile code")')
580  write(*,*)
581  stop
582  end if
583 case('lmaxo','lmaxvr')
584  read(50,'(A)',err=20) str
585  read(str,*,err=20) lmaxo
586  read(str,*,iostat=ios) lmaxo,dlmaxo
587  if (lmaxo < 3) then
588  write(*,*)
589  write(*,'("Error(readinput): lmaxo < 3 : ",I8)') lmaxo
590  write(*,*)
591  stop
592  end if
593 case('lmaxi','lmaxinr')
594  read(50,*,err=20) lmaxi
595  if (lmaxi < 1) then
596  write(*,*)
597  write(*,'("Error(readinput): lmaxi < 1 : ",I8)') lmaxi
598  write(*,*)
599  stop
600  end if
601 case('lmaxmat')
602  read(50,*,err=20)
603  write(*,'("Info(readinput): variable ''lmaxmat'' is no longer used")')
604 case('fracinr')
605  read(50,*,err=20) fracinr
606 case('trhonorm')
607  read(50,*,err=20) trhonorm
608 case('spinpol')
609  read(50,*,err=20) spinpol
610 case('spinorb')
611  read(50,*,err=20) spinorb
612 case('socscf')
613  read(50,*,err=20) socscf
614  if (socscf < 0.d0) then
615  write(*,*)
616  write(*,'("Error(readinput): socscf < 0 : ",G18.10)') socscf
617  write(*,*)
618  stop
619  end if
620 case('bforb')
621  read(50,*,err=20) bforb
622 case('bfdmag')
623  read(50,*,err=20) bfdmag
624 case('xctype')
625  read(50,'(A)',err=20) str
626  str=trim(str)//' 0 0'
627  read(str,*,err=20) xctype(:)
628 case('xctsp')
629  read(50,'(A)',err=20) str
630  str=trim(str)//' 0 0'
631  read(str,*,err=20) xctsp(:)
632 case('ktype')
633  read(50,'(A)',err=20) str
634  str=trim(str)//' 0 0'
635  read(str,*,err=20) ktype(:)
636  if (ktype(3) /= 0) then
637  write(*,*)
638  write(*,'("Error(readinput): ktype(3) should be zero : ",I8)') ktype(3)
639  write(*,*)
640  stop
641  end if
642 case('stype')
643  read(50,*,err=20) stype
644 case('swidth')
645  read(50,*,err=20) swidth
646  if (swidth < 1.d-9) then
647  write(*,*)
648  write(*,'("Error(readinput): swidth too small or negative : ",G18.10)') &
649  swidth
650  write(*,*)
651  stop
652  end if
653 case('autoswidth')
654  read(50,*,err=20) autoswidth
655 case('mstar')
656  read(50,*,err=20) mstar
657  if (mstar <= 0.d0) then
658  write(*,*)
659  write(*,'("Error(readinput): mstar <= 0 : ",G18.10)') mstar
660  write(*,*)
661  stop
662  end if
663 case('epsocc')
664  read(50,*,err=20) epsocc
665  if (epsocc <= 0.d0) then
666  write(*,*)
667  write(*,'("Error(readinput): epsocc <= 0 : ",G18.10)') epsocc
668  write(*,*)
669  stop
670  end if
671 case('epschg')
672  read(50,*,err=20) epschg
673  if (epschg <= 0.d0) then
674  write(*,*)
675  write(*,'("Error(readinput): epschg <= 0 : ",G18.10)') epschg
676  write(*,*)
677  stop
678  end if
679 case('nempty','nempty0')
680  read(50,'(A)',err=20) str
681  read(str,*,err=20) nempty0
682  read(str,*,iostat=ios) nempty0,dnempty0
683  if (nempty0 <= 0.d0) then
684  write(*,*)
685  write(*,'("Error(readinput): nempty <= 0 : ",G18.10)') nempty0
686  write(*,*)
687  stop
688  end if
689 case('mixtype')
690  read(50,*,err=20) mixtype
691 case('mixsave')
692  read(50,*,err=20) mixsave
693 case('amixpm','beta0','betamax')
694  if (trim(block) == 'amixpm') then
695  read(50,*,err=20) amixpm(:)
696  else if (trim(block) == 'beta0') then
697  read(50,*,err=20) amixpm(1)
698  else
699  read(50,*,err=20) amixpm(2)
700  end if
701  if (amixpm(1) < 0.d0) then
702  write(*,*)
703  write(*,'("Error(readinput): beta0 [amixpm(1)] < 0 : ",G18.10)') amixpm(1)
704  write(*,*)
705  stop
706  end if
707  if ((amixpm(2) < 0.d0).or.(amixpm(2) > 1.d0)) then
708  write(*,*)
709  write(*,'("Error(readinput): betamax [amixpm(2)] not in [0,1] : ",G18.10)')&
710  amixpm(2)
711  write(*,*)
712  stop
713  end if
714 case('mixsdb')
715  read(50,*,err=20) mixsdb
716  if (mixsdb < 2) then
717  write(*,*)
718  write(*,'("Error(readinput): mixsdb < 2 : ",I8)') mixsdb
719  write(*,*)
720  stop
721  end if
722 case('broydpm')
723  read(50,*,err=20) broydpm(:)
724  if ((broydpm(1) < 0.d0).or.(broydpm(1) > 1.d0).or. &
725  (broydpm(2) < 0.d0).or.(broydpm(2) > 1.d0)) then
726  write(*,*)
727  write(*,'("Error(readinput): invalid Broyden mixing parameters : ",&
728  &2G18.10)') broydpm
729  write(*,*)
730  stop
731  end if
732 case('mixrho')
733  read(50,*,err=20) mixrho
734 case('maxscl')
735  read(50,*,err=20) maxscl
736  if (maxscl < 0) then
737  write(*,*)
738  write(*,'("Error(readinput): maxscl < 0 : ",I8)') maxscl
739  write(*,*)
740  stop
741  end if
742 case('epspot')
743  read(50,*,err=20) epspot
744 case('epsengy')
745  read(50,*,err=20) epsengy
746 case('epsforce')
747  read(50,*,err=20) epsforce
748 case('epsstress')
749  read(50,*,err=20) epsstress
750 case('sppath')
751  read(50,*,err=20) sppath
752  sppath=adjustl(sppath)
753 case('scrpath')
754  read(50,*,err=20) scrpath
755 case('molecule')
756  read(50,*,err=20) molecule
757 case('atoms')
758  read(50,*,err=20) nspecies
759  if (nspecies < 1) then
760  write(*,*)
761  write(*,'("Error(readinput): nspecies < 1 : ",I8)') nspecies
762  write(*,*)
763  stop
764  end if
765  if (nspecies > maxspecies) then
766  write(*,*)
767  write(*,'("Error(readinput): nspecies too large : ",I8)') nspecies
768  write(*,'("Adjust maxspecies in modmain and recompile code")')
769  write(*,*)
770  stop
771  end if
772  do is=1,nspecies
773  read(50,*,err=20) spfname(is)
774  spfname(is)=adjustl(spfname(is))
775  read(50,*,err=20) natoms(is)
776  if (natoms(is) < 1) then
777  write(*,*)
778  write(*,'("Error(readinput): natoms < 1 : ",I8)') natoms(is)
779  write(*,'(" for species ",I4)') is
780  write(*,*)
781  stop
782  end if
783  if (natoms(is) > maxatoms) then
784  write(*,*)
785  write(*,'("Error(readinput): natoms too large : ",I8)') natoms(is)
786  write(*,'(" for species ",I4)') is
787  write(*,'("Adjust maxatoms in modmain and recompile code")')
788  write(*,*)
789  stop
790  end if
791  do ia=1,natoms(is)
792  read(50,'(A)',err=20) str
793  read(str,*,err=20) atposl(:,ia,is)
794  read(str,*,iostat=ios) atposl(:,ia,is),bfcmt0(:,ia,is),datposl(:,ia,is)
795  end do
796  end do
797 case('plot1d')
798  read(50,*,err=20) nvp1d,npp1d
799  if (nvp1d < 1) then
800  write(*,*)
801  write(*,'("Error(readinput): nvp1d < 1 : ",I8)') nvp1d
802  write(*,*)
803  stop
804  end if
805  if (npp1d < nvp1d) then
806  write(*,*)
807  write(*,'("Error(readinput): npp1d < nvp1d : ",2I8)') npp1d,nvp1d
808  write(*,*)
809  stop
810  end if
811  if (allocated(vvlp1d)) deallocate(vvlp1d)
812  allocate(vvlp1d(3,nvp1d))
813  do i=1,nvp1d
814  read(50,*,err=20) vvlp1d(:,i)
815  end do
816 case('ip01d')
817  read(50,*,err=20) ip01d
818  if (ip01d < 1) then
819  write(*,*)
820  write(*,'("Error(readinput): ip01d < 1 : ",I8)') ip01d
821  write(*,*)
822  stop
823  end if
824 case('plot2d')
825  read(50,*,err=20) vclp2d(:,0)
826  read(50,*,err=20) vclp2d(:,1)
827  read(50,*,err=20) vclp2d(:,2)
828  read(50,*,err=20) np2d(:)
829  if ((np2d(1) < 1).or.(np2d(2) < 1)) then
830  write(*,*)
831  write(*,'("Error(readinput): np2d < 1 : ",2I8)') np2d
832  write(*,*)
833  stop
834  end if
835 case('plot3d')
836  read(50,*,err=20) vclp3d(:,0)
837  read(50,*,err=20) vclp3d(:,1)
838  read(50,*,err=20) vclp3d(:,2)
839  read(50,*,err=20) vclp3d(:,3)
840  read(50,*,err=20) np3d(:)
841  if ((np3d(1) < 1).or.(np3d(2) < 1).or.(np3d(3) < 1)) then
842  write(*,*)
843  write(*,'("Error(readinput): np3d < 1 : ",3I8)') np3d
844  write(*,*)
845  stop
846  end if
847 case('wplot','dos')
848  read(50,*,err=20) nwplot,ngrkf,nswplot
849  if (nwplot < 2) then
850  write(*,*)
851  write(*,'("Error(readinput): nwplot < 2 : ",I8)') nwplot
852  write(*,*)
853  stop
854  end if
855  if (ngrkf < 1) then
856  write(*,*)
857  write(*,'("Error(readinput): ngrkf < 1 : ",I8)') ngrkf
858  write(*,*)
859  stop
860  end if
861  if (nswplot < 0) then
862  write(*,*)
863  write(*,'("Error(readinput): nswplot < 0 : ",I8)') nswplot
864  write(*,*)
865  stop
866  end if
867  read(50,*,err=20) wplot(:)
868  if (wplot(1) > wplot(2)) then
869  write(*,*)
870  write(*,'("Error(readinput): wplot(1) > wplot(2) : ",2G18.10)') wplot
871  write(*,*)
872  stop
873  end if
874 case('dosocc')
875  read(50,*,err=20) dosocc
876 case('tpdos')
877  read(50,*,err=20) tpdos
878 case('dosmsum')
879  read(50,*,err=20) dosmsum
880 case('dosssum')
881  read(50,*,err=20) dosssum
882 case('lmirep')
883  read(50,*,err=20) lmirep
884 case('atpopt')
885  read(50,*,err=20) atpopt
886 case('maxatpstp','maxatmstp')
887  read(50,*,err=20) maxatpstp
888  if (maxatpstp < 1) then
889  write(*,*)
890  write(*,'("Error(readinput): maxatpstp < 1 : ",I8)') maxatpstp
891  write(*,*)
892  stop
893  end if
894 case('tau0atp','tau0atm')
895  read(50,*,err=20) tau0atp
896 case('deltast')
897  read(50,*,err=20) deltast
898  if (deltast <= 0.d0) then
899  write(*,*)
900  write(*,'("Error(readinput): deltast <= 0 : ",G18.10)') deltast
901  write(*,*)
902  stop
903  end if
904 case('avecref')
905  read(50,*,err=20) avecref(:,1)
906  read(50,*,err=20) avecref(:,2)
907  read(50,*,err=20) avecref(:,3)
908 case('latvopt')
909  read(50,*,err=20) latvopt
910 case('maxlatvstp')
911  read(50,*,err=20) maxlatvstp
912  if (maxlatvstp < 1) then
913  write(*,*)
914  write(*,'("Error(readinput): maxlatvstp < 1 : ",I8)') maxlatvstp
915  write(*,*)
916  stop
917  end if
918 case('tau0latv')
919  read(50,*,err=20) tau0latv
920 case('nstfsp')
921  read(50,*,err=20)
922  write(*,'("Info(readinput): variable ''nstfsp'' is no longer used")')
923 case('lradstp')
924  read(50,*,err=20) lradstp
925  if (lradstp < 1) then
926  write(*,*)
927  write(*,'("Error(readinput): lradstp < 1 : ",I8)') lradstp
928  write(*,*)
929  stop
930  end if
931 case('chgexs')
932  read(50,'(A)',err=20) str
933  read(str,*,err=20) chgexs
934  read(str,*,iostat=ios) chgexs,dchgexs
935 case('nprad')
936  read(50,*,err=20)
937  write(*,'("Info(readinput): variable ''nprad'' is no longer used")')
938 case('scissor')
939  read(50,*,err=20) scissor
940 case('noptcomp')
941  read(50,*,err=20) noptcomp
942  if ((noptcomp < 1).or.(noptcomp > 27)) then
943  write(*,*)
944  write(*,'("Error(readinput): noptcomp should be from 1 to 27 : ",I8)') &
945  noptcomp
946  write(*,*)
947  stop
948  end if
949 case('optcomp')
950  do i=1,27
951  read(50,'(A)',err=20) str
952  if (trim(str) == '') then
953  if (i == 1) then
954  write(*,*)
955  write(*,'("Error(readinput): empty optical component list")')
956  write(*,*)
957  stop
958  end if
959  noptcomp=i-1
960  goto 10
961  end if
962  str=trim(str)//' 1 1'
963  read(str,*,iostat=ios) optcomp(:,i)
964  if (ios /= 0) then
965  write(*,*)
966  write(*,'("Error(readinput): error reading optical component list")')
967  write(*,'("(blank line required after optcomp block)")')
968  write(*,*)
969  stop
970  end if
971  if (any(optcomp(:,i) < 1).or.any(optcomp(:,i) > 3)) then
972  write(*,*)
973  write(*,'("Error(readinput): invalid optcomp : ",3I8)') optcomp(:,i)
974  write(*,*)
975  stop
976  end if
977  end do
978  write(*,*)
979  write(*,'("Error(readinput): optical component list too long")')
980  write(*,*)
981  stop
982 case('intraband')
983  read(50,*,err=20) intraband
984 case('evaltol')
985  read(50,*,err=20) evaltol
986 case('deband')
987  read(50,*,err=20)
988  write(*,'("Info(readinput): variable ''deband'' is no longer used")')
989 case('epsband')
990  read(50,*,err=20) epsband
991  if (epsband <= 0.d0) then
992  write(*,*)
993  write(*,'("Error(readinput): epsband <= 0 : ",G18.10)') epsband
994  write(*,*)
995  stop
996  end if
997 case('demaxbnd')
998  read(50,*,err=20) demaxbnd
999  if (demaxbnd <= 0.d0) then
1000  write(*,*)
1001  write(*,'("Error(readinput): demaxbnd <= 0 : ",G18.10)') demaxbnd
1002  write(*,*)
1003  stop
1004  end if
1005 case('autolinengy')
1006  read(50,*,err=20) autolinengy
1007 case('dlefe')
1008  read(50,*,err=20) dlefe
1009 case('autodlefe')
1010  read(50,*,err=20) autodlefe
1011 case('deapw')
1012  read(50,*,err=20) deapw
1013  if (abs(deapw) < 1.d-8) then
1014  write(*,*)
1015  write(*,'("Error(readinput): invalid deapw : ",G18.10)') deapw
1016  write(*,*)
1017  stop
1018  end if
1019 case('delorb')
1020  read(50,*,err=20) delorb
1021  if (abs(delorb) < 1.d-8) then
1022  write(*,*)
1023  write(*,'("Error(readinput): invalid delorb : ",G18.10)') delorb
1024  write(*,*)
1025  stop
1026  end if
1027 case('bfieldc')
1028  read(50,'(A)',err=20) str
1029  read(str,*,err=20) bfieldc0(:)
1030  read(str,*,iostat=ios) bfieldc0(:),dbfieldc0(:)
1031 case('efieldc')
1032  read(50,*,err=20) efieldc(:)
1033 case('dmaxefc')
1034  read(50,*,err=20) dmaxefc
1035  if (dmaxefc < 0) then
1036  write(*,*)
1037  write(*,'("Error(readinput): dmaxefc < 0 : ",G18.10)') dmaxefc
1038  write(*,*)
1039  stop
1040  end if
1041 case('afieldc')
1042  read(50,'(A)',err=20) str
1043  read(str,*,err=20) afieldc(:)
1044  read(str,*,iostat=ios) afieldc(:),dafieldc(:)
1045 case('afspc')
1046  do i=1,3
1047  read(50,'(A)',err=20) str
1048  read(str,*,err=20) afspc(i,:)
1049  read(str,*,iostat=ios) afspc(i,:),dafspc(i,:)
1050  end do
1051 case('fsmtype','fixspin')
1052  read(50,*,err=20) fsmtype
1053 case('momfix')
1054  read(50,'(A)',err=20) str
1055  read(str,*,err=20) momfix(:)
1056  read(str,*,iostat=ios) momfix(:),dmomfix(:)
1057 case('momfixm')
1058  read(50,*,err=20) momfixm
1059  if (momfixm < 0.d0) then
1060  write(*,*)
1061  write(*,'("Error(readinput): momfixm < 0 : ",G18.10)') momfixm
1062  write(*,*)
1063  stop
1064  end if
1065 case('mommtfix')
1066  do ias=1,maxspecies*maxatoms
1067  read(50,'(A)',err=20) str
1068  if (trim(str) == '') goto 10
1069  read(str,*,iostat=ios) is,ia,mommtfix(:,ia,is)
1070  if (ios /= 0) then
1071  write(*,*)
1072  write(*,'("Error(readinput): error reading muffin-tin fixed spin &
1073  &moments")')
1074  write(*,'("(blank line required after mommtfix block)")')
1075  write(*,*)
1076  stop
1077  end if
1078  end do
1079 case('mommtfixm')
1080  do ias=1,maxspecies*maxatoms
1081  read(50,'(A)',err=20) str
1082  if (trim(str) == '') goto 10
1083  read(str,*,iostat=ios) is,ia,mommtfixm(ia,is)
1084  if (ios /= 0) then
1085  write(*,*)
1086  write(*,'("Error(readinput): error reading muffin-tin fixed spin &
1087  &moment magnitudes")')
1088  write(*,'("(blank line required after mommtfixm block)")')
1089  write(*,*)
1090  stop
1091  end if
1092  end do
1093 case('taufsm')
1094  read(50,*,err=20) taufsm
1095  if (taufsm < 0.d0) then
1096  write(*,*)
1097  write(*,'("Error(readinput): taufsm < 0 : ",G18.10)') taufsm
1098  write(*,*)
1099  stop
1100  end if
1101 case('autormt')
1102  read(50,*,err=20)
1103  write(*,'("Info(readinput): variable ''autormt'' is no longer used")')
1104 case('rmtdelta')
1105  read(50,*,err=20) rmtdelta
1106  if (rmtdelta < 0.d0) then
1107  write(*,*)
1108  write(*,'("Warning(readinput): rmtdelta < 0 : ",G18.10)') rmtdelta
1109  end if
1110 case('isgkmax')
1111  read(50,*,err=20) isgkmax
1112 case('nosym')
1113  read(50,*,err=20) lv
1114  if (lv) symtype=0
1115 case('symtype')
1116  read(50,*,err=20) symtype
1117  if ((symtype < 0).or.(symtype > 2)) then
1118  write(*,*)
1119  write(*,'("Error(readinput): symtype not defined : ",I8)') symtype
1120  write(*,*)
1121  stop
1122  end if
1123 case('deltaph')
1124  read(50,*,err=20) deltaph
1125  if (deltaph <= 0.d0) then
1126  write(*,*)
1127  write(*,'("Error(readinput): deltaph <= 0 : ",G18.10)') deltaph
1128  write(*,*)
1129  stop
1130  end if
1131 case('phwrite')
1132  read(50,*,err=20) nphwrt
1133  if (nphwrt < 1) then
1134  write(*,*)
1135  write(*,'("Error(readinput): nphwrt < 1 : ",I8)') nphwrt
1136  write(*,*)
1137  stop
1138  end if
1139  if (allocated(vqlwrt)) deallocate(vqlwrt)
1140  allocate(vqlwrt(3,nphwrt))
1141  do i=1,nphwrt
1142  read(50,*,err=20) vqlwrt(:,i)
1143  end do
1144 case('notes')
1145  if (allocated(notes)) deallocate(notes)
1146  allocate(notes(0))
1147  notelns=0
1148  do
1149  read(50,'(A)') str
1150  if (trim(str) == '') goto 10
1151  notes=[notes(1:notelns),str]
1152  notelns=notelns+1
1153  end do
1154 case('tforce')
1155  read(50,*,err=20) tforce
1156 case('tfibs')
1157  read(50,*,err=20)
1158  write(*,'("Info(readinput): variable ''tfibs'' is no longer used")')
1159 case('maxitoep')
1160  read(50,*,err=20) maxitoep
1161  if (maxitoep < 1) then
1162  write(*,*)
1163  write(*,'("Error(readinput): maxitoep < 1 : ",I8)') maxitoep
1164  write(*,*)
1165  stop
1166  end if
1167 case('tauoep')
1168  read(50,*,err=20)
1169  write(*,'("Info(readinput): variable ''tauoep'' is no longer used")')
1170 case('tau0oep')
1171  read(50,*,err=20) tau0oep
1172  if (tau0oep < 0.d0) then
1173  write(*,*)
1174  write(*,'("Error(readinput): tau0oep < 0 : ",G18.10)') tau0oep
1175  write(*,*)
1176  stop
1177  end if
1178 case('kstlist')
1179  do i=1,maxkst
1180  read(50,'(A)',err=20) str
1181  if (trim(str) == '') then
1182  if (i == 1) then
1183  write(*,*)
1184  write(*,'("Error(readinput): empty k-point and state list")')
1185  write(*,*)
1186  stop
1187  end if
1188  nkstlist=i-1
1189  goto 10
1190  end if
1191  str=trim(str)//' 1'
1192  read(str,*,iostat=ios) kstlist(:,i)
1193  if (ios /= 0) then
1194  write(*,*)
1195  write(*,'("Error(readinput): error reading k-point and state list")')
1196  write(*,'("(blank line required after kstlist block)")')
1197  write(*,*)
1198  stop
1199  end if
1200  end do
1201  write(*,*)
1202  write(*,'("Error(readinput): k-point and state list too long")')
1203  write(*,*)
1204  stop
1205 case('vklem')
1206  read(50,*,err=20) vklem
1207 case('deltaem')
1208  read(50,*,err=20) deltaem
1209  if (deltaem <= 0.d0) then
1210  write(*,*)
1211  write(*,'("Error(readinput): deltaem <= 0 : ",G18.10)') deltaem
1212  write(*,*)
1213  stop
1214  end if
1215 case('ndspem')
1216  read(50,*,err=20) ndspem
1217  if ((ndspem < 1).or.(ndspem > 4)) then
1218  write(*,*)
1219  write(*,'("Error(readinput): ndspem out of range : ",I8)') ndspem
1220  write(*,*)
1221  stop
1222  end if
1223 case('nosource')
1224  read(50,*,err=20) nosource
1225 case('spinsprl')
1226  read(50,*,err=20) spinsprl
1227 case('ssdph')
1228  read(50,*,err=20) ssdph
1229 case('vqlss')
1230  read(50,'(A)',err=20) str
1231  read(str,*,err=20) vqlss
1232  read(str,*,iostat=ios) vqlss,dvqlss
1233 case('nwrite')
1234  read(50,*,err=20) nwrite
1235 case('DFT+U','dft+u','lda+u')
1236  read(50,*,err=20) dftu,inpdftu
1237  do i=1,maxdftu
1238  read(50,'(A)',err=20) str
1239  if (trim(str) == '') then
1240  ndftu=i-1
1241  goto 10
1242  end if
1243  select case(inpdftu)
1244  case(1)
1245  read(str,*,iostat=ios) is,l,ujdu(1:2,i)
1246  case(2)
1247  read(str,*,iostat=ios) is,l,(fdu(k,i),k=0,2*l,2)
1248  case(3)
1249  read(str,*,iostat=ios) is,l,(edu(k,i),k=0,l)
1250  case(4)
1251  read(str,*,iostat=ios) is,l,lamdu(i)
1252  case(5)
1253  read(str,*,iostat=ios) is,l,udufix(i),dudufix(i)
1254  read(str,*,iostat=ios) is,l,udufix(i)
1255  case default
1256  write(*,*)
1257  write(*,'("Error(readinput): invalid inpdftu : ",I8)') inpdftu
1258  write(*,*)
1259  stop
1260  end select
1261  if (ios /= 0) then
1262  write(*,*)
1263  write(*,'("Error(readinput): error reading DFT+U parameters")')
1264  write(*,'("(blank line required after dft+u block)")')
1265  write(*,*)
1266  stop
1267  end if
1268  if ((is < 1).or.(is >= maxspecies)) then
1269  write(*,*)
1270  write(*,'("Error(readinput): invalid species number in dft+u block : ", &
1271  &I8)') is
1272  write(*,*)
1273  stop
1274  end if
1275  if (l < 0) then
1276  write(*,*)
1277  write(*,'("Error(readinput): l < 0 in dft+u block : ",I8)') l
1278  write(*,*)
1279  stop
1280  end if
1281  if (l > lmaxdm) then
1282  write(*,*)
1283  write(*,'("Error(readinput): l > lmaxdm in dft+u block : ",2I8)') l,lmaxdm
1284  write(*,*)
1285  stop
1286  end if
1287 ! check for repeated entries
1288  do j=1,i-1
1289  if ((is == isldu(1,j)).and.(l == isldu(2,j))) then
1290  write(*,*)
1291  write(*,'("Error(readinput): repeated entry in DFT+U block")')
1292  write(*,*)
1293  stop
1294  end if
1295  end do
1296  isldu(1,i)=is
1297  isldu(2,i)=l
1298  end do
1299  write(*,*)
1300  write(*,'("Error(readinput): too many DFT+U entries")')
1301  write(*,'("Adjust maxdftu in modmain and recompile code")')
1302  write(*,*)
1303  stop
1304 case('tmwrite','tmomlu')
1305  read(50,*,err=20) tmwrite
1306 case('readadu','readalu')
1307  read(50,*,err=20)
1308  write(*,'("Info(readinput): variable ''readadu'' is no longer used")')
1309 case('rdmxctype')
1310  read(50,*,err=20) rdmxctype
1311 case('rdmmaxscl')
1312  read(50,*,err=20) rdmmaxscl
1313  if (rdmmaxscl < 0) then
1314  write(*,*)
1315  write(*,'("Error(readinput): rdmmaxscl < 0 : ",I8)') rdmmaxscl
1316  write(*,*)
1317  end if
1318 case('maxitn')
1319  read(50,*,err=20) maxitn
1320 case('maxitc')
1321  read(50,*,err=20) maxitc
1322 case('taurdmn')
1323  read(50,*,err=20) taurdmn
1324  if (taurdmn < 0.d0) then
1325  write(*,*)
1326  write(*,'("Error(readinput): taurdmn < 0 : ",G18.10)') taurdmn
1327  write(*,*)
1328  stop
1329  end if
1330 case('taurdmc')
1331  read(50,*,err=20) taurdmc
1332  if (taurdmc < 0.d0) then
1333  write(*,*)
1334  write(*,'("Error(readinput): taurdmc < 0 : ",G18.10)') taurdmc
1335  write(*,*)
1336  stop
1337  end if
1338 case('rdmalpha')
1339  read(50,*,err=20) rdmalpha
1340  if ((rdmalpha <= 0.d0).or.(rdmalpha >= 1.d0)) then
1341  write(*,*)
1342  write(*,'("Error(readinput): rdmalpha not in (0,1) : ",G18.10)') rdmalpha
1343  write(*,*)
1344  stop
1345  end if
1346 case('rdmtemp')
1347  read(50,*,err=20) rdmtemp
1348  if (rdmtemp < 0.d0) then
1349  write(*,*)
1350  write(*,'("Error(readinput): rdmtemp < 0 : ",G18.10)') rdmtemp
1351  write(*,*)
1352  stop
1353  end if
1354 case('reducebf')
1355  read(50,*,err=20) reducebf
1356  if ((reducebf < 0.5d0).or.(reducebf > 1.d0)) then
1357  write(*,*)
1358  write(*,'("Error(readinput): reducebf not in [0.5,1] : ",G18.10)') reducebf
1359  write(*,*)
1360  stop
1361  end if
1362 case('ptnucl')
1363  read(50,*,err=20) ptnucl
1364 case('tefvr','tseqr')
1365  read(50,*,err=20) tefvr
1366 case('tefvit','tseqit')
1367  read(50,*,err=20) tefvit
1368 case('minitefv','minseqit')
1369  read(50,*,err=20)
1370  write(*,'("Info(readinput): variable ''minitefv'' is no longer used")')
1371 case('nefvit','maxitefv','maxseqit','nseqit')
1372  read(50,*,err=20) nefvit
1373  if (nefvit < 1) then
1374  write(*,*)
1375  write(*,'("Error(readinput): nefvit < 1 : ",I8)') nefvit
1376  write(*,*)
1377  stop
1378  end if
1379 case('befvit','bseqit')
1380  read(50,*,err=20) befvit
1381  if (befvit <= 0.d0) then
1382  write(*,*)
1383  write(*,'("Error(readinput): befvit <= 0 : ",G18.10)') befvit
1384  write(*,*)
1385  stop
1386  end if
1387 case('epsefvit','epsseqit')
1388  read(50,*,err=20)
1389  write(*,'("Info(readinput): variable ''epsefvit'' is no longer used")')
1390 case('tauseq')
1391  read(50,*,err=20)
1392  write(*,'("Info(readinput): variable ''tauseq'' is no longer used")')
1393 case('vecql')
1394  read(50,*,err=20) vecql(:)
1395 case('mustar')
1396  read(50,*,err=20) mustar
1397 case('sqaxis','sqados')
1398  read(50,*,err=20) sqaxis(:)
1399 case('test')
1400  read(50,*,err=20) test
1401 case('frozencr')
1402  read(50,*,err=20)
1403  write(*,'("Info(readinput): variable ''frozencr'' is no longer used")')
1404 case('spincore')
1405  read(50,*,err=20) spincore
1406 case('solscf')
1407  read(50,*,err=20) solscf
1408  if (solscf < 0.d0) then
1409  write(*,*)
1410  write(*,'("Error(readinput): solscf < 0 : ",G18.10)') solscf
1411  write(*,*)
1412  stop
1413  end if
1414 case('emaxelnes')
1415  read(50,*,err=20) emaxelnes
1416 case('wsfac')
1417  read(50,*,err=20) wsfac(:)
1418 case('vhmat')
1419  read(50,*,err=20) vhmat(1,:)
1420  read(50,*,err=20) vhmat(2,:)
1421  read(50,*,err=20) vhmat(3,:)
1422 case('reduceh')
1423  read(50,*,err=20) reduceh
1424 case('hybrid')
1425  read(50,*,err=20) hybrid0
1426 case('hybridc','hybmix')
1427  read(50,*,err=20) hybridc
1428  if ((hybridc < 0.d0).or.(hybridc > 1.d0)) then
1429  write(*,*)
1430  write(*,'("Error(readinput): invalid hybridc : ",G18.10)') hybridc
1431  write(*,*)
1432  stop
1433  end if
1434 case('ecvcut')
1435  read(50,*,err=20) ecvcut
1436 case('esccut')
1437  read(50,*,err=20) esccut
1438 case('nvbse')
1439  read(50,*,err=20) nvbse0
1440  if (nvbse0 < 0) then
1441  write(*,*)
1442  write(*,'("Error(readinput): nvbse < 0 : ",I8)') nvbse0
1443  write(*,*)
1444  stop
1445  end if
1446 case('ncbse')
1447  read(50,*,err=20) ncbse0
1448  if (ncbse0 < 0) then
1449  write(*,*)
1450  write(*,'("Error(readinput): ncbse < 0 : ",I8)') ncbse0
1451  write(*,*)
1452  stop
1453  end if
1454 case('istxbse')
1455  do i=1,maxxbse
1456  read(50,'(A)',err=20) str
1457  if (trim(str) == '') then
1458  if (i == 1) then
1459  write(*,*)
1460  write(*,'("Error(readinput): empty BSE extra valence state list")')
1461  write(*,*)
1462  stop
1463  end if
1464  nvxbse=i-1
1465  goto 10
1466  end if
1467  read(str,*,iostat=ios) istxbse(i)
1468  if (ios /= 0) then
1469  write(*,*)
1470  write(*,'("Error(readinput): error reading BSE valence state list")')
1471  write(*,'("(blank line required after istxbse block)")')
1472  write(*,*)
1473  stop
1474  end if
1475  end do
1476  write(*,*)
1477  write(*,'("Error(readinput): BSE extra valence state list too long")')
1478  write(*,*)
1479  stop
1480 case('jstxbse')
1481  do i=1,maxxbse
1482  read(50,'(A)',err=20) str
1483  if (trim(str) == '') then
1484  if (i == 1) then
1485  write(*,*)
1486  write(*,'("Error(readinput): empty BSE extra conduction state list")')
1487  write(*,*)
1488  stop
1489  end if
1490  ncxbse=i-1
1491  goto 10
1492  end if
1493  read(str,*,iostat=ios) jstxbse(i)
1494  if (ios /= 0) then
1495  write(*,*)
1496  write(*,'("Error(readinput): error reading BSE conduction state list")')
1497  write(*,'("(blank line required after jstxbse block)")')
1498  write(*,*)
1499  stop
1500  end if
1501  end do
1502  write(*,*)
1503  write(*,'("Error(readinput): BSE extra conduction state list too long")')
1504  write(*,*)
1505  stop
1506 case('bsefull')
1507  read(50,*,err=20) bsefull
1508 case('hxbse')
1509  read(50,*,err=20) hxbse
1510 case('hdbse')
1511  read(50,*,err=20) hdbse
1512 case('gmaxrf','gmaxrpa')
1513  read(50,*,err=20) gmaxrf
1514  if (gmaxrf < 0.d0) then
1515  write(*,*)
1516  write(*,'("Error(readinput): gmaxrf < 0 : ",G18.10)') gmaxrf
1517  write(*,*)
1518  stop
1519  end if
1520 case('mbwgrf')
1521  read(50,*,err=20) mbwgrf
1522 case('emaxrf')
1523  read(50,*,err=20) emaxrf
1524  if (emaxrf < 0.d0) then
1525  write(*,*)
1526  write(*,'("Error(readinput): emaxrf < 0 : ",G18.10)') emaxrf
1527  write(*,*)
1528  stop
1529  end if
1530 case('fxctype')
1531  read(50,'(A)',err=20) str
1532  str=trim(str)//' 0 0'
1533  read(str,*,err=20) fxctype
1534 case('fxclrc')
1535  read(50,'(A)',err=20) str
1536  str=trim(str)//' 0.0'
1537  read(str,*,err=20) fxclrc(:)
1538 case('ntemp')
1539  read(50,*,err=20) ntemp
1540  if (ntemp < 1) then
1541  write(*,*)
1542  write(*,'("Error(readinput): ntemp < 1 : ",I8)') ntemp
1543  write(*,*)
1544  stop
1545  end if
1546 case('trimvg')
1547  write(*,'("Info(readinput): variable ''trimvg'' is no longer used")')
1548  read(50,*,err=20)
1549 case('rndstate','rndseed')
1550  read(50,*,err=20) rndstate(0)
1551  rndstate(0)=abs(rndstate(0))
1552 case('rndatposc')
1553  read(50,*,err=20) rndatposc
1554 case('rndbfcmt')
1555  read(50,*,err=20) rndbfcmt
1556 case('rndavec')
1557  read(50,*,err=20) rndavec
1558 case('c_tb09')
1559  read(50,*,err=20) c_tb09
1560 ! set flag to indicate Tran-Blaha constant has been read in
1561  tc_tb09=.true.
1562 case('lowq','highq','vhighq','uhighq')
1563  read(50,*,err=20) lv
1564  if (lv) then
1565  if (trim(block) == 'lowq') then
1566  rgkmax=6.5d0
1567  gmaxvr=10.d0
1568  lmaxapw=7
1569  lmaxo=5
1570  nxlo=2
1571  lorbcnd=.true.
1572  radkpt=25.d0
1573  autokpt=.true.
1574  vkloff(:)=0.5d0
1575  nempty0=4.d0
1576  epspot=1.d-5
1577  epsengy=5.d-4
1578  epsforce=1.d-2
1579  epsstress=3.d-3
1580  autolinengy=.true.
1581  gmaxrf=2.5d0
1582  lradstp=6
1583  else if (trim(block) == 'highq') then
1584 ! parameter set for high-quality calculation
1585  rgkmax=max(rgkmax,8.d0)
1586  gmaxvr=max(gmaxvr,16.d0)
1587  lmaxapw=max(lmaxapw,9)
1588  lmaxo=max(lmaxo,7)
1589  nrmtscf=max(nrmtscf,1.5d0)
1590  nxlo=max(nxlo,2)
1591  lorbcnd=.true.
1592  radkpt=max(radkpt,50.d0)
1593  autokpt=.true.
1594  vkloff(:)=0.d0
1595  nempty0=max(nempty0,10.d0)
1596  epspot=min(epspot,1.d-7)
1597  epsengy=min(epsengy,1.d-5)
1598  epsforce=min(epsforce,5.d-4)
1599  epsstress=min(epsstress,1.d-3)
1600  autolinengy=.true.
1601  gmaxrf=max(gmaxrf,4.d0)
1602  else if (trim(block) == 'vhighq') then
1603 ! parameter set for very high-quality calculation
1604  rgkmax=max(rgkmax,9.d0)
1605  gmaxvr=max(gmaxvr,18.d0)
1606  lmaxapw=max(lmaxapw,11)
1607  lmaxo=max(lmaxo,9)
1608  nrmtscf=max(nrmtscf,2.d0)
1609  nxlo=max(nxlo,3)
1610  lorbcnd=.true.
1611  radkpt=max(radkpt,90.d0)
1612  autokpt=.true.
1613  vkloff(:)=0.d0
1614  nempty0=max(nempty0,20.d0)
1615  epspot=min(epspot,1.d-7)
1616  epsengy=min(epsengy,1.d-6)
1617  epsforce=min(epsforce,2.d-4)
1618  epsstress=min(epsstress,5.d-4)
1619  autolinengy=.true.
1620  gmaxrf=max(gmaxrf,5.d0)
1621  else
1622 ! parameter set for ultra high-quality calculation
1623  rgkmax=max(rgkmax,10.d0)
1624  gmaxvr=max(gmaxvr,20.d0)
1625  lmaxapw=max(lmaxapw,12)
1626  lmaxo=max(lmaxo,9)
1627  nrmtscf=max(nrmtscf,4.d0)
1628  nxlo=max(nxlo,3)
1629  lorbcnd=.true.
1630  radkpt=max(radkpt,120.d0)
1631  autokpt=.true.
1632  vkloff(:)=0.d0
1633  nempty0=max(nempty0,40.d0)
1634  epspot=min(epspot,1.d-7)
1635  epsengy=min(epsengy,1.d-6)
1636  epsforce=min(epsforce,1.d-4)
1637  epsstress=min(epsstress,2.d-4)
1638  autolinengy=.true.
1639  gmaxrf=max(gmaxrf,6.d0)
1640  end if
1641  if (mp_mpi) then
1642  write(*,*)
1643  write(*,'("Info(readinput): parameters set by ",A," option")') trim(block)
1644  write(*,'(" rgkmax : ",G18.10)') rgkmax
1645  write(*,'(" gmaxvr : ",G18.10)') gmaxvr
1646  write(*,'(" lmaxapw : ",I4)') lmaxapw
1647  write(*,'(" lmaxo : ",I4)') lmaxo
1648  write(*,'(" nrmtscf : ",G18.10)') nrmtscf
1649  write(*,'(" nxlo : ",I4)') nxlo
1650  write(*,'(" lorbcnd : ",L1)') lorbcnd
1651  write(*,'(" radkpt : ",G18.10)') radkpt
1652  write(*,'(" autokpt : ",L1)') autokpt
1653  write(*,'(" vkloff : ",3G18.10)') vkloff
1654  write(*,'(" nempty0 : ",G18.10)') nempty0
1655  write(*,'(" epspot : ",G18.10)') epspot
1656  write(*,'(" epsengy : ",G18.10)') epsengy
1657  write(*,'(" epsforce : ",G18.10)') epsforce
1658  write(*,'(" epsstress : ",G18.10)') epsstress
1659  write(*,'(" autolinengy : ",L1)') autolinengy
1660  write(*,'(" gmaxrf : ",G18.10)') gmaxrf
1661  if (trim(block) == 'lowq') then
1662  write(*,'(" lradstp : ",I4)') lradstp
1663  end if
1664  end if
1665  end if
1666 case('hmaxvr')
1667  read(50,*,err=20) hmaxvr
1668  if (hmaxvr < 0.d0) then
1669  write(*,*)
1670  write(*,'("Error(readinput): hmaxvr < 0 : ",G18.10)') hmaxvr
1671  write(*,*)
1672  stop
1673  end if
1674 case('hkmax')
1675  read(50,*,err=20) hkmax
1676  if (hkmax <= 0.d0) then
1677  write(*,*)
1678  write(*,'("Error(readinput): hkmax <= 0 : ",G18.10)') hkmax
1679  write(*,*)
1680  stop
1681  end if
1682 case('lorbcnd')
1683  read(50,*,err=20) lorbcnd
1684 case('lorbordc')
1685  read(50,*,err=20) lorbordc
1686  if (lorbordc < 2) then
1687  write(*,*)
1688  write(*,'("Error(readinput): lorbordc < 2 : ",I8)') lorbordc
1689  write(*,*)
1690  stop
1691  end if
1692  if (lorbordc > maxlorbord) then
1693  write(*,*)
1694  write(*,'("Error(readinput): lorbordc too large : ",I8)') lorbordc
1695  write(*,'("Adjust maxlorbord in modmain and recompile code")')
1696  write(*,*)
1697  stop
1698  end if
1699 case('nrmtscf')
1700  read(50,'(A)',err=20) str
1701  read(str,*,err=20) nrmtscf
1702  read(str,*,iostat=ios) nrmtscf,dnrmtscf
1703  if (nrmtscf < 0.5d0) then
1704  write(*,*)
1705  write(*,'("Error(readinput): nrmtscf < 0.5 : ",G18.10)') nrmtscf
1706  write(*,*)
1707  stop
1708  end if
1709 case('lmaxdb','lmaxdos')
1710  read(50,*,err=20) lmaxdb
1711  if (lmaxdb < 0) then
1712  write(*,*)
1713  write(*,'("Error(readinput): lmaxdb < 0 : ",I8)') lmaxdb
1714  write(*,*)
1715  stop
1716  end if
1717 case('epsdev')
1718  read(50,*,err=20) epsdev
1719  if (epsdev <= 0.d0) then
1720  write(*,*)
1721  write(*,'("Error(readinput): epsdev <= 0 : ",G18.10)') epsdev
1722  write(*,*)
1723  stop
1724  end if
1725 case('msmooth')
1726  read(50,*,err=20)
1727  write(*,'("Info(readinput): variable ''msmooth'' is no longer used")')
1728 case('npmae')
1729  read(50,*,err=20) npmae0
1730 case('wrtvars')
1731  read(50,*,err=20) wrtvars
1732 case('ftmtype')
1733  read(50,*,err=20) ftmtype
1734 case('tmomfix')
1735  write(*,*)
1736  write(*,'("Error(readinput): variable ''tmomfix'' is no longer used")')
1737  write(*,'(" use tm3fix instead")')
1738  write(*,*)
1739  stop
1740 case('tm3fix')
1741  read(50,*,err=20) ntmfix
1742  if (ntmfix < 1) then
1743  write(*,*)
1744  write(*,'("Error(readinput): ntmfix < 1 : ",I8)') ntmfix
1745  write(*,*)
1746  stop
1747  end if
1748  if (allocated(itmfix)) deallocate(itmfix)
1749  allocate(itmfix(7,ntmfix))
1750  if (allocated(wkprfix)) deallocate(wkprfix)
1751  allocate(wkprfix(ntmfix))
1752  do i=1,ntmfix
1753  read(50,*,err=20) is,ia,l
1754  if ((is < 1).or.(ia < 1).or.(l < 0)) then
1755  write(*,*)
1756  write(*,'("Error(readinput): invalid is, ia or l in tm3fix block : ",&
1757  &3I8)') is,ia,l
1758  write(*,*)
1759  stop
1760  end if
1761  itmfix(1,i)=is
1762  itmfix(2,i)=ia
1763  itmfix(3,i)=l
1764 ! read k, p, r, t for the 3-index tensor
1765  read(50,*,err=20) itmfix(4:7,i)
1766 ! read 3-index tensor component with conventional normalisation
1767  read(50,*,err=20) wkprfix(i)
1768  end do
1769 case('tauftm')
1770  read(50,*,err=20) tauftm
1771  if (tauftm < 0.d0) then
1772  write(*,*)
1773  write(*,'("Error(readinput): tauftm < 0 : ",G18.10)') tauftm
1774  write(*,*)
1775  stop
1776  end if
1777 case('ftmstep')
1778  read(50,*,err=20)
1779  write(*,'("Info(readinput): variable ''ftmstep'' is no longer used")')
1780 case('cmagz','forcecmag')
1781  read(50,*,err=20) cmagz
1782 case('rotavec')
1783  read(50,*,err=20) axang(:)
1784 case('tstime')
1785  read(50,*,err=20) tstime
1786  if (tstime <= 0.d0) then
1787  write(*,*)
1788  write(*,'("Error(readinput): tstime <= 0 : ",G18.10)') tstime
1789  write(*,*)
1790  stop
1791  end if
1792 case('dtimes')
1793  read(50,*,err=20) dtimes
1794  if (dtimes <= 0.d0) then
1795  write(*,*)
1796  write(*,'("Error(readinput): dtimes <= 0 : ",G18.10)') dtimes
1797  write(*,*)
1798  stop
1799  end if
1800 case('pulse')
1801  read(50,*,err=20) npulse
1802  if (npulse < 1) then
1803  write(*,*)
1804  write(*,'("Error(readinput): npulse < 1 : ",I8)') npulse
1805  write(*,*)
1806  stop
1807  end if
1808  if (allocated(pulse)) deallocate(pulse)
1809  allocate(pulse(12,npulse))
1810  do i=1,npulse
1811  read(50,'(A)',err=20) str
1812  str=trim(str)//' 1.0 0.0 0.0 0.0'
1813  read(str,*,err=20) pulse(:,i)
1814  end do
1815 case('ramp')
1816  read(50,*,err=20) nramp
1817  if (nramp < 1) then
1818  write(*,*)
1819  write(*,'("Error(readinput): nramp < 1 : ",I8)') nramp
1820  write(*,*)
1821  stop
1822  end if
1823  if (allocated(ramp)) deallocate(ramp)
1824  allocate(ramp(12,nramp))
1825  do i=1,nramp
1826  read(50,'(A)',err=20) str
1827  str=trim(str)//' 1.0 0.0 0.0 0.0'
1828  read(str,*,err=20) ramp(:,i)
1829  end do
1830 case('step')
1831  read(50,*,err=20) nstep
1832  if (nstep < 1) then
1833  write(*,*)
1834  write(*,'("Error(readinput): nstep < 1 : ",I8)') nstep
1835  write(*,*)
1836  stop
1837  end if
1838  if (allocated(step)) deallocate(step)
1839  allocate(step(9,nstep))
1840  do i=1,nstep
1841  read(50,'(A)',err=20) str
1842  str=trim(str)//' 1.0 0.0 0.0 0.0'
1843  read(str,*,err=20) step(:,i)
1844  end do
1845 case('ncgga')
1846  read(50,*,err=20)
1847  write(*,'("Info(readinput): variable ''ncgga'' is no longer used")')
1848 case('dncgga')
1849  read(50,*,err=20) dncgga
1850  if (dncgga < 0.d0) then
1851  write(*,*)
1852  write(*,'("Error(readinput): dncgga < 0 : ",G18.10)') dncgga
1853  write(*,*)
1854  stop
1855  end if
1856 case('ntswrite')
1857  read(50,'(A)',err=20) str
1858  str=trim(str)//' 1'
1859  read(str,*,err=20) ntswrite(:)
1860 case('nxoapwlo','nxapwlo')
1861  read(50,*,err=20) nxoapwlo
1862  if (nxoapwlo < 0) then
1863  write(*,*)
1864  write(*,'("Error(readinput): nxoapwlo < 0 : ",I8)') nxoapwlo
1865  write(*,*)
1866  stop
1867  end if
1868 case('nxlo')
1869  read(50,*,err=20) nxlo
1870  if (nxlo < 0) then
1871  write(*,*)
1872  write(*,'("Error(readinput): nxlo < 0 : ",I8)') nxlo
1873  write(*,*)
1874  stop
1875  end if
1876 case('tdrho1d')
1877  read(50,*,err=20) tdrho1d
1878 case('tdrho2d')
1879  read(50,*,err=20) tdrho2d
1880 case('tdrho3d')
1881  read(50,*,err=20) tdrho3d
1882 case('tdmag1d')
1883  read(50,*,err=20) tdmag1d
1884 case('tdmag2d')
1885  read(50,*,err=20) tdmag2d
1886 case('tdmag3d')
1887  read(50,*,err=20) tdmag3d
1888 case('tdjr1d','tdcd1d')
1889  read(50,*,err=20) tdjr1d
1890 case('tdjr2d','tdcd2d')
1891  read(50,*,err=20) tdjr2d
1892 case('tdjr3d','tdcd3d')
1893  read(50,*,err=20) tdjr3d
1894 case('tddos')
1895  read(50,*,err=20) tddos
1896 case('tdlsj')
1897  read(50,*,err=20) tdlsj
1898 case('tdjtk')
1899  read(50,*,err=20) tdjtk
1900 case('tdxrmk')
1901  read(50,*,err=20) tdxrmk
1902 case('epseph')
1903  read(50,*,err=20)
1904  write(*,'("Info(readinput): variable ''epseph'' is no longer used")')
1905 case('rndevt0')
1906  read(50,*,err=20) rndevt0
1907 case('sxcscf','ssxc','rstsf')
1908  read(50,'(A)',err=20) str
1909  read(str,*,err=20) sxcscf
1910  read(str,*,iostat=ios) sxcscf,dsxcscf
1911 case('tempk')
1912  read(50,*,err=20) tempk
1913  if (tempk <= 0.d0) then
1914  write(*,*)
1915  write(*,'("Error(readinput): tempk <= 0 : ",G18.10)') tempk
1916  write(*,*)
1917  stop
1918  end if
1919 ! set Fermi-Dirac smearing
1920  stype=3
1921 ! set the smearing width
1923 case('avecu')
1924  read(50,*,err=20) avecu(:,1)
1925  read(50,*,err=20) avecu(:,2)
1926  read(50,*,err=20) avecu(:,3)
1927 case('scaleu')
1928  read(50,*,err=20) scu
1929 case('scaleu1')
1930  read(50,*,err=20) scu1
1931 case('scaleu2')
1932  read(50,*,err=20) scu2
1933 case('scaleu3')
1934  read(50,*,err=20) scu3
1935 case('q0cut')
1936  read(50,*,err=20) q0cut
1937 case('rndbfcu')
1938  read(50,*,err=20) rndbfcu
1939 case('bfieldcu','bfielduc')
1940  read(50,*,err=20) bfieldcu
1941 case('efieldcu','efielduc')
1942  read(50,*,err=20) efieldcu
1943 case('tplotq0')
1944  read(50,*,err=20) tplotq0
1945 case('trdvclr')
1946  read(50,*,err=20) trdvclr
1947 case('trdbfcr')
1948  read(50,*,err=20) trdbfcr
1949 case('evtype')
1950  read(50,*,err=20)
1951  write(*,'("Info(readinput): variable ''evtype'' is no longer used")')
1952 case('wmaxgw')
1953  read(50,*,err=20) wmaxgw
1954 case('twdiag')
1955  read(50,*,err=20)
1956  write(*,'("Info(readinput): variable ''twdiag'' is no longer used")')
1957 case('tsediag')
1958  read(50,*,err=20) tsediag
1959 case('actype')
1960  read(50,*,err=20) actype
1961 case('npole')
1962  read(50,*,err=20) npole
1963  if (npole < 1) then
1964  write(*,*)
1965  write(*,'("Error(readinput): npole < 1 : ",I8)') npole
1966  write(*,*)
1967  stop
1968  end if
1969 case('nspade')
1970  read(50,*,err=20) nspade
1971  if (nspade < 1) then
1972  write(*,*)
1973  write(*,'("Error(readinput): nspade < 1 : ",I8)') nspade
1974  write(*,*)
1975  stop
1976  end if
1977 case('tfav0')
1978  read(50,*,err=20) tfav0
1979 case('rmtscf')
1980  read(50,*,err=20) rmtscf
1981  if (rmtscf <= 0.d0) then
1982  write(*,*)
1983  write(*,'("Error(readinput): rmtscf <= 0 : ",G18.10)') rmtscf
1984  write(*,*)
1985  stop
1986  end if
1987 case('mrmtav')
1988  read(50,*,err=20) mrmtav
1989 case('rmtall')
1990  read(50,*,err=20) rmtall
1991 case('maxthd','omp_num_threads','OMP_NUM_THREADS')
1992  read(50,*,err=20) maxthd
1993 case('maxthd1')
1994  read(50,*,err=20) maxthd1
1995 case('maxthdmkl')
1996  read(50,*,err=20) maxthdmkl
1997 case('maxlvl','omp_max_active_levels','OMP_MAX_ACTIVE_LEVELS')
1998  read(50,*,err=20) maxlvl
1999  if (maxlvl < 1) then
2000  write(*,*)
2001  write(*,'("Error(readinput): maxlvl < 1 : ",I8)') maxlvl
2002  write(*,*)
2003  stop
2004  end if
2005 case('stable')
2006  read(50,*,err=20) lv
2007  if (lv) then
2008  autolinengy=.true.
2009  mrmtav=1
2010  lmaxapw=max(lmaxapw,10)
2011  gmaxvr=max(gmaxvr,24.d0)
2012  msmgmt=max(msmgmt,4)
2013  if (mp_mpi) then
2014  write(*,*)
2015  write(*,'("Info(readinput): parameters set by stable option")')
2016  write(*,'(" autolinengy : ",L1)') autolinengy
2017  write(*,'(" mrmtav : ",I4)') mrmtav
2018  write(*,'(" lmaxapw : ",I4)') lmaxapw
2019  write(*,'(" gmaxvr : ",G18.10)') gmaxvr
2020  write(*,'(" msmgmt : ",I4)') msmgmt
2021  end if
2022  end if
2023 case('metagga')
2024  read(50,*,err=20) lv
2025  if (lv) then
2026  lmaxi=max(lmaxi,2)
2027  gmaxvr=max(gmaxvr,16.d0)
2028  nrmtscf=max(nrmtscf,3.d0)
2029  msmgmt=max(msmgmt,4)
2030  epspot=1.d6
2031  epsengy=min(epsengy,1.d-6)
2032  if (mp_mpi) then
2033  write(*,*)
2034  write(*,'("Info(readinput): parameters set by metagga option")')
2035  write(*,'(" lmaxi : ",I4)') lmaxi
2036  write(*,'(" gmaxvr : ",G18.10)') gmaxvr
2037  write(*,'(" nrmtscf : ",G18.10)') nrmtscf
2038  write(*,'(" msmgmt : ",I4)') msmgmt
2039  write(*,'(" epspot : ",G18.10)') epspot
2040  write(*,'(" epsengy : ",G18.10)') epsengy
2041  end if
2042  end if
2043 case('t0tdlr')
2044  read(50,*,err=20)
2045  write(*,'("Info(readinput): variable ''t0tdlr'' is no longer used")')
2046 case('tdphi')
2047  read(50,*,err=20) tdphi
2048 ! convert phase from degrees to radians
2049  tdphi=tdphi*pi/180.d0
2050 case('thetamld')
2051  read(50,*,err=20) thetamld
2052 ! convert MLD angle from degrees to radians
2053  thetamld=thetamld*pi/180.d0
2054 case('ntsbackup')
2055  read(50,*,err=20) ntsbackup
2056 case('seedname')
2057  read(50,*,err=20) seedname
2058  seedname=adjustl(seedname)
2059 case('num_wann')
2060  read(50,*,err=20) num_wann
2061 case('idxw90','wann_bands')
2062  read(50,'(A)',err=20) str
2063  num_bands=1024
2064  if (allocated(idxw90)) deallocate(idxw90)
2065  allocate(idxw90(num_bands))
2066  call numlist(str,num_bands,idxw90)
2067 case('num_iter')
2068  read(50,*,err=20) num_iter
2069 case('dis_num_iter')
2070  read(50,*,err=20) dis_num_iter
2071 case('trial_step')
2072  read(50,*,err=20) trial_step
2073 case('xlwin','wannierExtra')
2074  if (allocated(xlwin)) deallocate(xlwin)
2075  allocate(xlwin(0))
2076  nxlwin=0
2077  do
2078  read(50,'(A)',err=20) str
2079  if (trim(str) == '') goto 10
2080  xlwin=[xlwin(1:nxlwin),str]
2081  nxlwin=nxlwin+1
2082  end do
2083 case('wrtunk')
2084  read(50,*,err=20) wrtunk
2085 case('tbdip')
2086  read(50,*,err=20) tbdip
2087 case('tjr','tcden')
2088  read(50,*,err=20) tjr
2089 case('tauefm')
2090  read(50,*,err=20) tauefm
2091 case('epsefm')
2092  read(50,*,err=20) epsefm
2093 case('ehfb')
2094  read(50,*,err=20) ehfb
2095 case('t0gclq0')
2096  read(50,*,err=20) t0gclq0
2097 case('tafindt')
2098  read(50,*,err=20) tafindt
2099 case('afindscf')
2100  read(50,*,err=20)
2101  write(*,'("Info(readinput): variable ''afindscf'' is no longer used")')
2102 case('afindpm')
2103  read(50,*,err=20) afindpm(:)
2104  if (afindpm(2) == 0.d0) then
2105  write(*,*)
2106  write(*,'("Error(readinput): afindpm(2) = 0")')
2107  write(*,*)
2108  stop
2109  end if
2110 case('nkspolar')
2111  read(50,*,err=20) nkspolar
2112  if (nkspolar < 1) then
2113  write(*,*)
2114  write(*,'("Error(readinput): nkspolar < 1 : ",I8)') nkspolar
2115  write(*,*)
2116  stop
2117  end if
2118 case('ntsforce')
2119  read(50,*,err=20) ntsforce
2120  if (ntsforce < 1) then
2121  write(*,*)
2122  write(*,'("Error(readinput): ntsforce < 1 : ",I8)') ntsforce
2123  write(*,*)
2124  stop
2125  end if
2126 case('wphcut')
2127  read(50,*,err=20) wphcut
2128  if (wphcut <= 0.d0) then
2129  write(*,*)
2130  write(*,'("Error(readinput): wphcut <= 0 : ",G18.10)') wphcut
2131  write(*,*)
2132  stop
2133  end if
2134 case('ephscf')
2135  read(50,*,err=20) ephscf(:)
2136 case('anomalous')
2137  read(50,*,err=20) anomalous
2138 case('tephde')
2139  read(50,*,err=20) tephde
2140 case('bdiag')
2141  read(50,*,err=20) bdiag
2142 case('ecutb')
2143  read(50,*,err=20) ecutb
2144  if (ecutb <= 0.d0) then
2145  write(*,*)
2146  write(*,'("Error(readinput): ecutb <= 0 : ",G18.10)') ecutb
2147  write(*,*)
2148  stop
2149  end if
2150 case('ediag')
2151  read(50,*,err=20) ediag
2152 case('pwxpsn')
2153  read(50,*,err=20) pwxpsn
2154  if (pwxpsn < 1) then
2155  write(*,*)
2156  write(*,'("Error(readinput): pwxpsn < 1 : ",I8)') pwxpsn
2157  write(*,*)
2158  stop
2159  end if
2160 case('ramdisk')
2161  read(50,*,err=20) ramdisk
2162 case('wrtdsk')
2163  read(50,*,err=20) wrtdsk
2164 case('epsdmat')
2165  read(50,*,err=20) epsdmat
2166 case('tm3vdl','tm3old')
2167  read(50,*,err=20) tm3vdl
2168 case('batch')
2169  read(50,*,err=20) batch
2170 case('tafspt')
2171  read(50,*,err=20) tafspt
2172 case('tbaspat')
2173  read(50,*,err=20) tbaspat
2174 case('trdatdv')
2175  read(50,*,err=20) trdatdv
2176 case('atdfc')
2177  read(50,*,err=20) atdfc
2178  if (atdfc < 0.d0) then
2179  write(*,*)
2180  write(*,'("Error(readinput): atdfc < 0 : ",G18.10)') atdfc
2181  write(*,*)
2182  stop
2183  end if
2184 case('maxforce')
2185  read(50,*,err=20) maxforce
2186 case('msmgmt','msmg2mt')
2187  read(50,*,err=20) msmgmt
2188 case('ntsorth')
2189  read(50,*,err=20) ntsorth
2190 case('deltabf')
2191  read(50,*,err=20) deltabf
2192  if (deltabf <= 0.d0) then
2193  write(*,*)
2194  write(*,'("Error(readinput): deltabf <= 0 : ",G18.10)') deltabf
2195  write(*,*)
2196  stop
2197  end if
2198 case('jtconst0')
2199  read(50,*,err=20) jtconst0
2200 case('trmt0')
2201  read(50,*,err=20) trmt0
2202 case('ksgwrho')
2203  read(50,*,err=20) ksgwrho
2204 case('npfftg')
2205  read(50,*,err=20) npfftg
2206 case('npfftgc')
2207  read(50,*,err=20) npfftgc
2208 case('npfftq')
2209  read(50,*,err=20) npfftq
2210 case('npfftw')
2211  read(50,*,err=20) npfftw
2212 case('tphnat')
2213  read(50,*,err=20) tphnat
2214 case('ecutthc')
2215  read(50,*,err=20) ecutthc
2216  if (ecutthc <= 0.d0) then
2217  write(*,*)
2218  write(*,'("Error(readinput): ecutthc <= 0 : ",G18.10)') ecutthc
2219  write(*,*)
2220  stop
2221  end if
2222 case('tbdipu')
2223  read(50,*,err=20) tbdipu
2224 case('bdipscf')
2225  read(50,*,err=20) bdipscf
2226 case('')
2227  goto 10
2228 case default
2229  write(*,*)
2230  write(*,'("Error(readinput): invalid block name : ",A)') trim(block)
2231  write(*,*)
2232  stop
2233 end select
2234 goto 10
2235 20 continue
2236 write(*,*)
2237 write(*,'("Error(readinput): error reading from elk.in")')
2238 write(*,'("Problem occurred in ''",A,"'' block")') trim(block)
2239 write(*,'("Check input convention in manual")')
2240 write(*,*)
2241 stop
2242 30 continue
2243 close(50)
2244 ! scale the speed of light
2245 solsc=sol*solscf
2246 ! scale and rotate the lattice vectors (not referenced again in code)
2247 avec(:,:)=sc*avec(:,:)
2248 avec(:,1)=sc1*avec(:,1)
2249 avec(:,2)=sc2*avec(:,2)
2250 avec(:,3)=sc3*avec(:,3)
2251 avec(1,:)=scx*avec(1,:)
2252 avec(2,:)=scy*avec(2,:)
2253 avec(3,:)=scz*avec(3,:)
2254 t1=axang(4)
2255 if (t1 /= 0.d0) then
2256  t1=t1*pi/180.d0
2257  call axangrot(axang(:),t1,rot)
2258  do i=1,3
2259  v(:)=avec(:,i)
2260  call r3mv(rot,v,avec(:,i))
2261  end do
2262 end if
2263 ! randomise lattice vectors if required
2264 if (rndavec > 0.d0) then
2265  do i=1,3
2266  do j=1,3
2267  t1=rndavec*(randomu()-0.5d0)
2268  avec(i,j)=avec(i,j)+t1
2269  end do
2270  end do
2271 end if
2272 ! check if reference lattice vectors should be used
2273 tavref=(any(abs(avecref(:,:)) > epslat))
2274 ! case of isolated molecule
2275 if (molecule) then
2276 ! convert atomic positions from Cartesian to lattice coordinates
2277  call r3minv(avec,ainv)
2278  do is=1,nspecies
2279  do ia=1,natoms(is)
2280  call r3mv(ainv,atposl(:,ia,is),v)
2281  atposl(:,ia,is)=v(:)
2282  end do
2283  end do
2284 end if
2285 ! randomise atomic positions if required
2286 if (rndatposc > 0.d0) then
2287  call r3minv(avec,ainv)
2288  do is=1,nspecies
2289  do ia=1,natoms(is)
2290  call r3mv(avec,atposl(:,ia,is),v)
2291  do i=1,3
2292  t1=rndatposc*(randomu()-0.5d0)
2293  v(i)=v(i)+t1
2294  end do
2295  call r3mv(ainv,v,atposl(:,ia,is))
2296  end do
2297  end do
2298 end if
2299 ! randomise the muffin-tin magnetic fields if required
2300 if (rndbfcmt > 0.d0) then
2301  do is=1,nspecies
2302  do ia=1,natoms(is)
2303  do i=1,3
2304  t1=rndbfcmt*(randomu()-0.5d0)
2305  bfcmt0(i,ia,is)=bfcmt0(i,ia,is)+t1
2306  end do
2307  end do
2308  end do
2309 end if
2310 ! set fxctype to fxctype if required
2311 if (fxctype(1) == -1) fxctype(:)=xctype(:)
2312 ! find primitive cell if required
2313 if (primcell) call findprimcell
2314 ! scale the ultracell vectors if required
2315 avecu(:,1)=scu1*avecu(:,1)
2316 avecu(:,2)=scu2*avecu(:,2)
2317 avecu(:,3)=scu3*avecu(:,3)
2318 avecu(:,:)=scu*avecu(:,:)
2319 ! read in atomic species data
2320 call readspecies
2321 return
2322 
2323 end subroutine
2324 !EOC
2325 
real(8) socscf
Definition: modmain.f90:232
integer ngrkf
Definition: modmain.f90:1075
real(8) ecvcut
Definition: modmain.f90:117
integer maxlatvstp
Definition: modmain.f90:1039
real(8), dimension(3, 3) dafspc
Definition: modmain.f90:331
real(8) epsband
Definition: modmain.f90:820
subroutine genfspecies(zn, symb)
Definition: genfspecies.f90:7
character(256) scrpath
Definition: modmain.f90:1303
real(8) dchgexs
Definition: modmain.f90:724
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition: modmain.f90:275
integer, parameter maxspecies
Definition: modmain.f90:30
integer ncxbse
Definition: modmain.f90:1185
real(8) scissor
Definition: modmain.f90:911
real(8), dimension(3, 3) afspc
Definition: modmain.f90:331
subroutine findprimcell
real(8) dmaxefc
Definition: modmain.f90:318
real(8) epsforce
Definition: modmain.f90:1065
logical tjr
Definition: modmain.f90:620
integer, dimension(maxxbse) jstxbse
Definition: modmain.f90:1187
integer lorbordc
Definition: modmain.f90:837
integer maxthd1
Definition: modomp.f90:13
integer mixtype
Definition: modmain.f90:695
integer, dimension(3) ktype
Definition: modmain.f90:606
real(8), dimension(3, 3) davec
Definition: modmain.f90:12
real(8), dimension(3, maxatoms, maxspecies) datposl
Definition: modmain.f90:52
integer maxthd
Definition: modomp.f90:11
integer, dimension(3) xctype
Definition: modmain.f90:588
logical spinpol
Definition: modmain.f90:228
integer, parameter maxlapw
Definition: modmain.f90:195
integer npfftg
Definition: modmain.f90:404
integer mbwgrf
Definition: modmain.f90:1166
integer, dimension(3) xctsp
Definition: modmain.f90:142
logical tshift
Definition: modmain.f90:352
character(256) sppath
Definition: modmain.f90:72
real(8) rmtall
Definition: modmain.f90:158
real(8) taufsm
Definition: modmain.f90:265
logical tefvit
Definition: modmain.f90:873
real(8) reducebf
Definition: modmain.f90:279
logical autokpt
Definition: modmain.f90:444
integer msmgmt
Definition: modmain.f90:222
logical tfav0
Definition: modmain.f90:1003
real(8) dnrmtscf
Definition: modmain.f90:148
real(8) sxcscf
Definition: modmain.f90:668
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
real(8), dimension(2) amixpm
Definition: modmain.f90:702
logical lorbcnd
Definition: modmain.f90:835
real(8) rmtscf
Definition: modmain.f90:154
subroutine r3minv(a, b)
Definition: r3minv.f90:10
logical spinsprl
Definition: modmain.f90:283
Definition: modomp.f90:6
real(8), parameter kboltz
Definition: modmain.f90:1263
real(8) swidth
Definition: modmain.f90:895
logical hybrid0
Definition: modmain.f90:1152
real(8) epsocc
Definition: modmain.f90:903
real(8) dncgga
Definition: modmain.f90:604
character(256), dimension(maxspecies) spfname
Definition: modmain.f90:74
real(8) dnempty0
Definition: modmain.f90:883
real(8), dimension(3, 0:3) vclp3d
Definition: modmain.f90:1132
integer, dimension(maxxbse) istxbse
Definition: modmain.f90:1187
integer nvxbse
Definition: modmain.f90:1185
integer notelns
Definition: modmain.f90:1305
integer, parameter maxxbse
Definition: modmain.f90:1183
subroutine numlist(str, n, list)
Definition: numlist.f90:10
logical mixrho
Definition: modmain.f90:687
logical nosource
Definition: modmain.f90:664
real(8), dimension(3) vkloff
Definition: modmain.f90:450
integer maxthdmkl
Definition: modomp.f90:15
integer lmaxo
Definition: modmain.f90:201
real(8), dimension(3) vqlss
Definition: modmain.f90:293
logical trhonorm
Definition: modmain.f90:618
integer, parameter maxkst
Definition: modmain.f90:925
real(8) c_tb09
Definition: modmain.f90:678
real(8), parameter pi
Definition: modmain.f90:1232
real(8), dimension(3) vecql
Definition: modmain.f90:1104
real(8) nempty0
Definition: modmain.f90:883
integer nkstlist
Definition: modmain.f90:927
logical dosmsum
Definition: modmain.f90:1087
integer nxlo
Definition: modmain.f90:841
logical tc_tb09
Definition: modmain.f90:680
integer dlmaxo
Definition: modmain.f90:201
integer ntasks
Definition: modmain.f90:1293
integer lmaxapw
Definition: modmain.f90:197
real(8) nrmtscf
Definition: modmain.f90:148
integer symtype
Definition: modmain.f90:340
real(8) epschg
Definition: modmain.f90:712
real(8), dimension(2) broydpm
Definition: modmain.f90:706
logical tefvr
Definition: modmain.f90:870
real(8), dimension(3) momfix
Definition: modmain.f90:253
integer mixsdb
Definition: modmain.f90:704
logical bforb
Definition: modmain.f90:234
logical tforce
Definition: modmain.f90:989
real(8) rmtdelta
Definition: modmain.f90:160
real(8), dimension(3, 3) avecref
Definition: modmain.f90:1030
Definition: modrdm.f90:6
integer ncbse0
Definition: modmain.f90:1181
integer maxlvl
Definition: modomp.f90:17
real(8), dimension(2) wplot
Definition: modmain.f90:1079
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
real(8) emaxrf
Definition: modmain.f90:1162
real(8) gmaxrf
Definition: modmain.f90:1160
real(8) delorb
Definition: modmain.f90:802
logical lmirep
Definition: modmain.f90:1098
real(8) esccut
Definition: modmain.f90:119
Definition: modw90.f90:6
real(8), dimension(3) sqaxis
Definition: modmain.f90:1101
logical intraband
Definition: modmain.f90:1095
real(8) rndatposc
Definition: modmain.f90:56
real(8) radkpt
Definition: modmain.f90:446
integer, parameter maxlorbord
Definition: modmain.f90:782
integer nvbse0
Definition: modmain.f90:1181
Definition: modbog.f90:6
real(8) deltaem
Definition: modmain.f90:481
integer lmaxdb
Definition: modmain.f90:1081
real(8), dimension(3) dbfieldc0
Definition: modmain.f90:271
pure subroutine axangrot(v, th, rot)
Definition: axangrot.f90:10
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
real(8) rndbfcmt
Definition: modmain.f90:277
real(8), dimension(2) wsfac
Definition: modmain.f90:1108
real(8) tempk
Definition: modmain.f90:684
integer, dimension(3, 27) optcomp
Definition: modmain.f90:1093
integer, dimension(2, maxkst) kstlist
Definition: modmain.f90:929
real(8) rgkmax
Definition: modmain.f90:493
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8) maxforce
Definition: modmain.f90:1000
logical cmagz
Definition: modmain.f90:242
integer lradstp
Definition: modmain.f90:171
integer, parameter maxatoms
Definition: modmain.f90:32
real(8), dimension(3) dafieldc
Definition: modmain.f90:325
integer, parameter maxtasks
Definition: modmain.f90:1291
real(8), parameter sol
Definition: modmain.f90:1251
logical tpdos
Definition: modmain.f90:1085
real(8), dimension(3) dvqlss
Definition: modmain.f90:293
real(8), dimension(3) afieldc
Definition: modmain.f90:325
integer latvopt
Definition: modmain.f90:1037
real(8), dimension(3, 0:2) vclp2d
Definition: modmain.f90:1128
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition: modmain.f90:259
logical mixsave
Definition: modmain.f90:700
integer ip01d
Definition: modmain.f90:1118
logical tavref
Definition: modmain.f90:1032
integer npfftq
Definition: modmain.f90:535
real(8) solsc
Definition: modmain.f90:1253
real(8) tau0atp
Definition: modmain.f90:1011
real(8) evaltol
Definition: modmain.f90:919
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8) momfixm
Definition: modmain.f90:255
real(8) drgkmax
Definition: modmain.f90:493
logical t0gclq0
Definition: modmain.f90:555
Definition: modgw.f90:6
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nwrite
Definition: modmain.f90:1059
logical ssdph
Definition: modmain.f90:285
integer noptcomp
Definition: modmain.f90:1091
real(8) epslat
Definition: modmain.f90:24
integer maxitoep
Definition: modmain.f90:1142
real(8), dimension(maxatoms, maxspecies) mommtfixm
Definition: modmain.f90:261
integer ndspem
Definition: modmain.f90:483
integer stype
Definition: modmain.f90:891
logical hxbse
Definition: modmain.f90:1209
Definition: modmpi.f90:6
subroutine readspecies
Definition: readspecies.f90:7
logical dosssum
Definition: modmain.f90:1089
real(8) mstar
Definition: modmain.f90:899
real(8) epspot
Definition: modmain.f90:1061
real(8) dsxcscf
Definition: modmain.f90:668
integer nefvit
Definition: modmain.f90:875
logical ptnucl
Definition: modmain.f90:83
Definition: modpw.f90:6
real(8), dimension(:,:), allocatable vvlp1d
Definition: modmain.f90:1120
real(8) dlefe
Definition: modmain.f90:830
logical autolinengy
Definition: modmain.f90:828
integer nswplot
Definition: modmain.f90:1077
logical tbdip
Definition: modmain.f90:643
real(8) deapw
Definition: modmain.f90:764
logical spinorb
Definition: modmain.f90:230
integer maxatpstp
Definition: modmain.f90:1009
real(8), dimension(3) bfieldc0
Definition: modmain.f90:271
real(8) tau0oep
Definition: modmain.f90:1144
integer, dimension(3) dngridk
Definition: modmain.f90:448
subroutine readinput
Definition: readinput.f90:10
integer nspecies
Definition: modmain.f90:34
integer, dimension(2) np2d
Definition: modmain.f90:1130
logical bsefull
Definition: modmain.f90:1206
real(8), dimension(3) dmomfix
Definition: modmain.f90:253
logical autodlefe
Definition: modmain.f90:833
integer, dimension(3) np3d
Definition: modmain.f90:1134
integer reducek
Definition: modmain.f90:455
real(8) emaxelnes
Definition: modmain.f90:1106
real(8) tau0latv
Definition: modmain.f90:1041
real(8) bdipscf
Definition: modmain.f90:645
logical autoswidth
Definition: modmain.f90:897
integer nkspolar
Definition: modmain.f90:485
integer reduceq
Definition: modmain.f90:519
character(256), dimension(:), allocatable notes
Definition: modmain.f90:1307
real(8) befvit
Definition: modmain.f90:877
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
Definition: modulr.f90:6
integer nwplot
Definition: modmain.f90:1073
real(8) deltabf
Definition: modmain.f90:281
integer mrmtav
Definition: modmain.f90:156
real(8) epsstress
Definition: modmain.f90:1067
real(8) chgexs
Definition: modmain.f90:724
real(8) demaxbnd
Definition: modmain.f90:823
integer nvp1d
Definition: modmain.f90:1114
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
subroutine genspecies(fnum)
Definition: genspecies.f90:7
logical hdbse
Definition: modmain.f90:1209
integer, dimension(maxtasks) tasks
Definition: modmain.f90:1297
integer npp1d
Definition: modmain.f90:1116
logical molecule
Definition: modmain.f90:47
real(8) epsengy
Definition: modmain.f90:1063
real(8), dimension(3) efieldc
Definition: modmain.f90:312
integer dlmaxapw
Definition: modmain.f90:197
integer npfftgc
Definition: modmain.f90:408
integer nxoapwlo
Definition: modmain.f90:839
real(8), dimension(3) vklem
Definition: modmain.f90:479
real(8) hybridc
Definition: modmain.f90:1154
logical dosocc
Definition: modmain.f90:1083
logical trmt0
Definition: modmain.f90:165
integer maxscl
Definition: modmain.f90:1049
real(8) dgmaxvr
Definition: modmain.f90:384
logical spincore
Definition: modmain.f90:943
integer npmae0
Definition: modmain.f90:302
integer isgkmax
Definition: modmain.f90:491
real(8) atdfc
Definition: modmain.f90:66
logical primcell
Definition: modmain.f90:49
integer atpopt
Definition: modmain.f90:1007
real(8) fracinr
Definition: modmain.f90:209
real(8) deltast
Definition: modmain.f90:1022
integer fsmtype
Definition: modmain.f90:251
logical trdatdv
Definition: modmain.f90:62
real(8) gmaxvr
Definition: modmain.f90:384
integer lmaxi
Definition: modmain.f90:205
logical bfdmag
Definition: modmain.f90:236