The Elk Code
init0.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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: init0
8 ! !INTERFACE:
9 subroutine init0
10 ! !USES:
11 use modmain
12 use modxcifc
13 use moddftu
14 use modtddft
15 use modphonon
16 use modulr
17 use modgw
18 use modtest
19 use modvars
20 use modmpi
21 use modomp
22 ! !DESCRIPTION:
23 ! Performs basic consistency checks as well as allocating and initialising
24 ! global variables not dependent on the $k$-point set.
25 !
26 ! !REVISION HISTORY:
27 ! Created January 2004 (JKD)
28 !EOP
29 !BOC
30 implicit none
31 ! local variables
32 logical hybrid_
33 integer is,ia,ias,npc
34 integer idu,ist,nr,l,n,i
35 integer xcspin_
36 real(8) hybridc_,t1
37 real(8) ts0,ts1
38 
39 !-------------------------------!
40 ! zero timing variables !
41 !-------------------------------!
42 timeinit=0.d0
43 timemat=0.d0
44 timefv=0.d0
45 timesv=0.d0
46 timerho=0.d0
47 timepot=0.d0
48 timefor=0.d0
49 call timesec(ts0)
50 
51 !------------------------------------!
52 ! angular momentum variables !
53 !------------------------------------!
54 if (lmaxo > lmaxapw) then
55  write(*,*)
56  write(*,'("Error(init0): lmaxo > lmaxapw :",2(X,I0))') lmaxo,lmaxapw
57  write(*,*)
58  stop
59 end if
60 lmaxi=min(lmaxi,lmaxo)
61 lmmaxapw=(lmaxapw+1)**2
62 lmmaxi=(lmaxi+1)**2
63 lmmaxo=(lmaxo+1)**2
64 ! check DOS and band structure angular momentum maximum is within range
65 lmaxdb=min(lmaxdb,lmaxo)
66 lmmaxdb=(lmaxdb+1)**2
67 ! write to VARIABLES.OUT
68 if (wrtvars) then
69  call writevars('lmaxapw',iv=lmaxapw)
70  call writevars('lmaxi',iv=lmaxi)
71  call writevars('lmaxo',iv=lmaxo)
72  call writevars('lmaxdb',iv=lmaxdb)
73 end if
74 
75 !------------------------------------!
76 ! index to atoms and species !
77 !------------------------------------!
78 natmmax=0
79 ias=0
80 do is=1,nspecies
81  do ia=1,natoms(is)
82  ias=ias+1
83  idxas(ia,is)=ias
84  idxis(ias)=is
85  idxia(ias)=ia
86  end do
87 ! maximum number of atoms over all species
88  natmmax=max(natmmax,natoms(is))
89 end do
90 ! total number of atoms
91 natmtot=ias
92 ! number of phonon branches
93 nbph=3*natmtot
94 ! write to VARIABLES.OUT
95 if (wrtvars) then
96  call writevars('nspecies',iv=nspecies)
97  call writevars('natoms',nv=nspecies,iva=natoms)
98  call writevars('spsymb',nv=nspecies,sva=spsymb)
99  call writevars('spname',nv=nspecies,sva=spname)
100  call writevars('spzn',nv=nspecies,rva=spzn)
101 end if
102 
103 !------------------------!
104 ! spin variables !
105 !------------------------!
106 if (spinsprl) then
107  spinpol=.true.
108  spinorb=.false.
109  if (any(task == [5,51,52,53,61,62,63,700,701])) then
110  write(*,*)
111  write(*,'("Error(init0): spin-spirals do not work with task ",I0)') task
112  write(*,*)
113  stop
114  end if
115  if (xctype(1) < 0) then
116  write(*,*)
117  write(*,'("Error(init0): spin-spirals do not work with the OEP method")')
118  write(*,*)
119  stop
120  end if
121  if (tefvit) then
122  write(*,*)
123  write(*,'("Error(init0): spin-spirals do not work with iterative &
124  &diagonalisation")')
125  write(*,*)
126  stop
127  end if
128 end if
129 ! de-phasing required only for spin-spirals
130 if (.not.spinsprl) ssdph=.false.
131 ! spin-orbit coupling, B-field-orbit coupling, fixed spin moment, spin spirals
132 ! or spin-polarised cores requires a spin-polarised calculation
133 if (spinorb.or.bforb.or.(fsmtype /= 0).or.spinsprl.or.spincore) spinpol=.true.
134 ! number of spinor components and maximum allowed occupancy
135 if (spinpol) then
136  nspinor=2
137  occmax=1.d0
138 else
139  nspinor=1
140  occmax=2.d0
141 end if
142 ! number of spin-dependent first-variational functions per state and map from
143 ! second- to first-variational spin index
144 if (spinsprl) then
145  nspnfv=2
146  jspnfv(1)=1
147  jspnfv(2)=2
148 else
149  nspnfv=1
150  jspnfv(1)=1
151  jspnfv(2)=1
152 end if
153 ! no calculation of second-variational eigenvectors by default
154 tevecsv=.false.
155 ! spin-polarised calculations require second-variational eigenvectors
156 if (spinpol) tevecsv=.true.
157 ! Hartree-Fock/RDMFT/TDDFT/GW/ULR requires second-variational eigenvectors
158 if (any(task == [5,10,170,300,460,461,462,463,600,601,620,700,701,720,725]).or.&
159  ksgwrho) then
160  tevecsv=.true.
161 end if
162 ! get exchange-correlation functional data
164 if ((spinpol).and.(xcspin == 0)) then
165  write(*,*)
166  write(*,'("Error(init0): requested spin-polarised run with &
167  &spin-unpolarised")')
168  write(*,'(" exchange-correlation functional")')
169  write(*,*)
170  stop
171 end if
172 ! set flag for hybrid functional
173 if (task == 5) then
175 else
176  hybrid=.false.
177 end if
178 ! check for collinearity in the z-direction and set the dimension of the
179 ! magnetisation and exchange-correlation vector fields
180 if (spinpol) then
181  ndmag=1
182  if ((abs(bfieldc0(1)) > epslat).or.(abs(bfieldc0(2)) > epslat)) ndmag=3
183  do is=1,nspecies
184  do ia=1,natoms(is)
185  if ((abs(bfcmt0(1,ia,is)) > epslat).or. &
186  (abs(bfcmt0(2,ia,is)) > epslat)) ndmag=3
187  end do
188  end do
189 ! spin-orbit coupling is non-collinear in general
190  if (spinorb) ndmag=3
191 ! source-free fields and spin-spirals must be non-collinear
192  if (nosource.or.spinsprl) then
193  ndmag=3
194  cmagz=.false.
195  end if
196 ! force collinear magnetism along the z-axis if required
197  if (cmagz) ndmag=1
198 else
199  ndmag=0
200 end if
201 ! set the non-collinear flag
202 ncmag=(ndmag == 3)
203 ! check for meta-GGA with non-collinearity
204 if (any(xcgrad == [3,4,5,6]).and.ncmag) then
205  write(*,*)
206  write(*,'("Error(init0): meta-GGA is not valid for non-collinear magnetism")')
207  write(*,*)
208  stop
209 end if
210 if (tbdip.and.(.not.ncmag)) then
211  write(*,*)
212  write(*,'("Error(init0): non-collinear magnetism required for inclusion of &
213  &the dipole field")')
214  write(*,*)
215  stop
216 end if
217 ! spin-polarised cores
218 if (.not.spinpol) spincore=.false.
219 if (fsmtype /= 0) then
220 ! set fixed spin moment effective field to zero
221  bfsmc(:)=0.d0
222 ! set muffin-tin FSM fields to zero
223  if (allocated(bfsmcmt)) deallocate(bfsmcmt)
224  allocate(bfsmcmt(3,natmtot))
225  bfsmcmt(:,:)=0.d0
226  if (mixtype > 1) then
227  mixtype=1
228  if (mp_mpi) then
229  write(*,'("Info(init0): mixtype changed to 1 for FSM calculation")')
230  end if
231  end if
232 end if
233 ! number of independent spin components of the f_xc spin tensor
234 if (spinpol) then
235  if (ncmag) then
236  nscfxc=10
237  else
238  nscfxc=3
239  end if
240 else
241  nscfxc=1
242 end if
243 ! set the magnetic fields to the initial values
244 bfieldc(:)=bfieldc0(:)
245 bfcmt(:,:,:)=bfcmt0(:,:,:)
246 if (tmwrite.or.(ftmtype /= 0).or.(task == 400)) then
247  if (.not.spinorb) then
248  write(*,*)
249  write(*,'("Error(init0): tensor moments require spin-orbit coupling &
250  &enabled")')
251  write(*,'(" set spinorb=.true.")')
252  write(*,*)
253  stop
254  end if
255 end if
256 ! generate the fixed tensor moment density matrices if required
257 if (ftmtype /= 0) call gendmftm
258 ! if reducebf < 1 then reduce the external magnetic fields and tensor moments
259 ! immediately for non-self-consistent calculations
260 if (reducebf < 1.d0-1.d-4) then
261  if (all(task /= [0,1,2,3,28,200,201,208,350,351,380,420,421,440])) then
262  bfieldc(:)=0.d0
263  bfcmt(:,:,:)=0.d0
264  if (ftmtype < 0) dmftm(:,:,:,:,:)=0.d0
265  end if
266 end if
267 ! write to VARIABLES.OUT
268 if (wrtvars) then
269  call writevars('nspinor',iv=nspinor)
270  call writevars('ndmag',iv=ndmag)
271  call writevars('sxcscf',rv=sxcscf)
272  call writevars('bfieldc0',nv=3,rva=bfieldc0)
273 end if
274 
275 !----------------------------------!
276 ! crystal structure set up !
277 !----------------------------------!
278 ! use reference lattice vectors if required
279 if (tavref) then
280  avec0(:,:)=avec(:,:)
281  avec(:,:)=avecref(:,:)
282 end if
283 ! generate the reciprocal lattice vectors, unit cell and Brillouin zone volumes
285 ! inverse of the lattice vector matrix
286 call r3minv(avec,ainv)
287 ! inverse of the reciprocal vector matrix
288 call r3minv(bvec,binv)
289 ! Cartesian coordinates of the spin-spiral vector
290 call r3mv(bvec,vqlss,vqcss)
291 do is=1,nspecies
292  do ia=1,natoms(is)
293 ! map atomic lattice coordinates to [0,1)
294  call r3frac(epslat,atposl(:,ia,is))
295 ! determine atomic Cartesian coordinates
296  call r3mv(avec,atposl(:,ia,is),atposc(:,ia,is))
297  end do
298 end do
299 ! check for overlapping muffin-tins and adjust radii if required
300 call checkmt
301 ! compute the total muffin-tin volume (M. Meinert)
302 omegamt=0.d0
303 do is=1,nspecies
304  omegamt=omegamt+dble(natoms(is))*(fourpi/3.d0)*rmt(is)**3
305 end do
306 ! input q-vector in Cartesian coordinates
307 call r3mv(bvec,vecql,vecqc)
308 
309 !-------------------------------!
310 ! vector fields E and A !
311 !-------------------------------!
312 ! static electric field
313 tefield=.false.
314 if (any(abs(efieldc(:)) > epslat)) then
315 ! no shift of the atomic positions
316  tshift=.false.
317 ! electric field vector in lattice coordinates
318  call r3mv(ainv,efieldc,efieldl)
319 ! potential at maximum distance
320  vmaxefc=dmaxefc*norm2(efieldc(1:3))
321 ! allocate array for average electric field in each muffin-tin
322  if (allocated(efcmt)) deallocate(efcmt)
323  allocate(efcmt(3,natmtot))
324 ! set the E-field flag
325  tefield=.true.
326 end if
327 ! static vector potential
328 tafield=.false.
329 if (any(abs(afieldc(:)) > epslat)) then
330  tafield=.true.
331 ! A-field in lattice coordinates
332  call r3mv(ainv,afieldc,afieldl)
333 ! vector potential added in second-variational step
334  tevecsv=.true.
335 end if
336 ! static spin-dependent vector potential
337 tafsp=.false.
338 if (spinpol.and.(any(abs(afspc(:,:)) > epslat))) then
339  tafsp=.true.
340  tevecsv=.true.
341 end if
342 ! time-dependent vector potential
343 tafieldt=.false.
344 if (any(task == [460,461,462,463,480,481,485])) then
345 ! read time-dependent A-field from file
346  call readafieldt
347  tafieldt=.true.
348 ! zero the induced A-field and its time derivative
349  afindt(:,:)=0.d0
350 end if
351 ! write to VARIABLES.OUT
352 if (wrtvars) then
353  call writevars('efieldc',nv=3,rva=efieldc)
354  call writevars('afieldc',nv=3,rva=afieldc)
355  call writevars('afspc',nv=9,rva=afspc)
356 end if
357 
358 !---------------------------------!
359 ! crystal symmetry set up !
360 !---------------------------------!
361 call symmetry
362 
363 !-----------------------!
364 ! radial meshes !
365 !-----------------------!
366 nrmtmax=1
367 nrcmtmax=1
368 do is=1,nspecies
369 ! make the muffin-tin mesh commensurate with lradstp
370  nrmt(is)=nrmt(is)-mod(nrmt(is)-1,lradstp)
371  nrmtmax=max(nrmtmax,nrmt(is))
372 ! number of coarse radial mesh points
373  nrcmt(is)=(nrmt(is)-1)/lradstp+1
374  nrcmtmax=max(nrcmtmax,nrcmt(is))
375 end do
376 ! set up atomic and muffin-tin radial meshes
377 call genrmesh
378 ! number of points in packed muffin-tins
379 npmtmax=1
380 npcmtmax=1
381 do is=1,nspecies
382  npmti(is)=lmmaxi*nrmti(is)
383  npmt(is)=npmti(is)+lmmaxo*(nrmt(is)-nrmti(is))
384  npmtmax=max(npmtmax,npmt(is))
385  npcmti(is)=lmmaxi*nrcmti(is)
386  npcmt(is)=npcmti(is)+lmmaxo*(nrcmt(is)-nrcmti(is))
387  npcmtmax=max(npcmtmax,npcmt(is))
388 end do
389 ! index to first muffin-tin point in packed array over all atoms
390 if (allocated(ipcmt)) deallocate(ipcmt)
391 allocate(ipcmt(natmtot))
392 ipcmt(1)=1
393 npc=npcmt(1)
394 npcmttot=npc
395 do ias=2,natmtot
396  is=idxis(ias)
397  ipcmt(ias)=ipcmt(ias)+npc
398  npc=npcmt(is)
399  npcmttot=npcmttot+npc
400 end do
401 
402 !--------------------------------------!
403 ! charges and number of states !
404 !--------------------------------------!
405 chgzn=0.d0
406 chgcrtot=0.d0
407 chgval=0.d0
408 nstspmax=0
409 nstcr=0
410 do is=1,nspecies
411 ! nuclear charge
412  chgzn=chgzn+spzn(is)*natoms(is)
413 ! find the maximum number of atomic states
414  nstspmax=max(nstspmax,nstsp(is))
415 ! compute the electronic charge for each species, as well as the total core and
416 ! valence charge
417  spze(is)=0.d0
418  chgcr(is)=0.d0
419  do ist=1,nstsp(is)
420  spze(is)=spze(is)+occsp(ist,is)
421  if (spcore(ist,is)) then
422  chgcr(is)=chgcr(is)+occsp(ist,is)
423  nstcr=nstcr+2*ksp(ist,is)*natoms(is)
424  else
425  chgval=chgval+occsp(ist,is)*natoms(is)
426  end if
427  end do
428  chgcrtot=chgcrtot+chgcr(is)*natoms(is)
429 end do
430 ! add excess charge
432 ! total charge
434 if (chgtot < 1.d-8) then
435  write(*,*)
436  write(*,'("Error(init0): zero total charge")')
437  write(*,*)
438  stop
439 end if
440 ! subtract small charge to ensure consistent Fermi energy for insulators
442 ! effective Wigner radius
443 rwigner=(3.d0/(fourpi*(chgtot/omega)))**(1.d0/3.d0)
444 ! write to VARIABLES.OUT
445 if (wrtvars) then
446  call writevars('spze',nv=nspecies,rva=spze)
447  call writevars('chgcr',nv=nspecies,rva=chgcr)
448  call writevars('chgexs',rv=chgexs)
449  call writevars('chgval',rv=chgtot)
450 end if
451 
452 !-------------------------!
453 ! G-vector arrays !
454 !-------------------------!
455 ! determine gkmax from rgkmax
456 if (nspecies == 0) isgkmax=-2
457 select case(isgkmax)
458 case(:-4)
459 ! use largest muffin-tin radius
460  gkmax=rgkmax/maxval(rmt(1:nspecies))
461 case(-3)
462 ! use smallest muffin-tin radius
463  gkmax=rgkmax/minval(rmt(1:nspecies))
464 case(-2)
465 ! use the fixed value of 2.0
466  gkmax=rgkmax/2.d0
467 case(-1,0)
468 ! use average muffin-tin radius
469  t1=sum(natoms(1:nspecies)*rmt(1:nspecies))/dble(natmtot)
470  gkmax=rgkmax/t1
471 case(1:)
472 ! use user-specified muffin-tin radius
473  if (isgkmax <= nspecies) then
475  else
476  write(*,*)
477  write(*,'("Error(init0): isgkmax > nspecies :",2(X,I0))') isgkmax,nspecies
478  write(*,*)
479  stop
480  end if
481 end select
482 ! generate the G-vectors
483 call gengvec
484 ! apply strain to A-, B- and G-vectors if required
485 call strainabg
486 ! write number of G-vectors to test file
487 call writetest(900,'number of G-vectors',iv=ngvec)
488 ! Poisson solver pseudocharge density constant
489 if (nspecies > 0) then
490  t1=0.25d0*gmaxvr*maxval(rmt(1:nspecies))
491 else
492  t1=0.25d0*gmaxvr*2.d0
493 end if
494 npsd=max(nint(t1),1)
495 lnpsd=lmaxo+npsd+1
496 ! generate the Coulomb Green's function in G-space = 4π/G²
497 call gengclg
498 ! compute the spherical Bessel functions j_l(|G|R_mt)
499 if (allocated(jlgrmt)) deallocate(jlgrmt)
500 allocate(jlgrmt(0:lnpsd,ngvec,nspecies))
502 ! generate the spherical harmonics of the G-vectors
503 call genylmg
504 ! allocate structure factor array for G-vectors
505 if (allocated(sfacg)) deallocate(sfacg)
506 allocate(sfacg(ngvec,natmtot))
507 ! generate structure factors for G-vectors
509 ! generate the smooth step function form factors
510 if (allocated(ffacg)) deallocate(ffacg)
511 allocate(ffacg(ngtot,nspecies))
513 ! generate the smooth characteristic function
514 call gencfun
515 ! G-vector variables for coarse grid with |G| < 2 gkmax
516 call gengvc
517 ! generate the characteristic function on the coarse grid
518 call gencfrc
519 ! write to VARIABLES.OUT
520 if (wrtvars) then
521  call writevars('avec',nv=9,rva=avec)
522  call writevars('bvec',nv=9,rva=bvec)
523  call writevars('omega',rv=omega)
524  do is=1,nspecies
525  call writevars('atposl',n1=is,nv=3*natoms(is),rva=atposl(:,:,is))
526  end do
527  do is=1,nspecies
528  call writevars('atposc',n1=is,nv=3*natoms(is),rva=atposc(:,:,is))
529  end do
530  call writevars('vqlss',nv=3,rva=vqlss)
531  call writevars('vqcss',nv=3,rva=vqcss)
532  call writevars('gmaxvr',rv=gmaxvr)
533  call writevars('ngridg',nv=3,iva=ngridg)
534  call writevars('intgv',nv=6,iva=intgv)
535  call writevars('ngvec',iv=ngvec)
536  call writevars('ivg',nv=3*ngtot,iva=ivg)
537  call writevars('igfft',nv=ngtot,iva=igfft)
538 end if
539 
540 !-------------------------!
541 ! atoms and cores !
542 !-------------------------!
543 ! determine the nuclear Coulomb potential
544 if (allocated(vcln)) deallocate(vcln)
545 allocate(vcln(nrspmax,nspecies))
546 do is=1,nspecies
547  nr=nrsp(is)
548  call potnucl(ptnucl,nr,rsp(:,is),spzn(is),vcln(:,is))
549  vcln(1:nr,is)=vcln(1:nr,is)*y00i
550 end do
551 ! solve the Kohn-Sham-Dirac equations for all atoms
552 call allatoms
553 ! allocate core state occupancy and eigenvalue arrays and set to default
554 if (allocated(occcr)) deallocate(occcr)
555 allocate(occcr(nstspmax,natmtot))
556 if (allocated(evalcr)) deallocate(evalcr)
557 allocate(evalcr(nstspmax,natmtot))
558 do ias=1,natmtot
559  is=idxis(ias)
560  do ist=1,nstsp(is)
561  occcr(ist,ias)=occsp(ist,is)
562  evalcr(ist,ias)=evalsp(ist,is)
563  end do
564 end do
565 ! allocate core state radial wavefunction array
566 if (allocated(rwfcr)) deallocate(rwfcr)
567 allocate(rwfcr(nrspmax,2,nstspmax,natmtot))
568 ! number of core spin channels
569 if (spincore) then
570  nspncr=2
571 else
572  nspncr=1
573 end if
574 ! allocate core state charge density array
575 if (allocated(rhocr)) deallocate(rhocr)
576 allocate(rhocr(nrmtmax,natmtot,nspncr))
577 
578 !-------------------------------------------------------------!
579 ! charge density, potentials and exchange-correlation !
580 !-------------------------------------------------------------!
581 ! combined target array for density and magnetisation
582 if (allocated(rhmg)) deallocate(rhmg)
584 if (spinpol) n=n*(1+ndmag)
585 allocate(rhmg(n))
586 ! associate pointer arrays with target
587 rhomt(1:npmtmax,1:natmtot) => rhmg(1:)
588 i=size(rhomt)+1
589 rhoir(1:ngtot) => rhmg(i:)
590 if (spinpol) then
591  i=i+size(rhoir)
592  magmt(1:npmtmax,1:natmtot,1:ndmag) => rhmg(i:)
593  i=i+size(magmt)
594  magir(1:ngtot,1:ndmag) => rhmg(i:)
595 end if
596 if (any(task == [371,372,373]).or.tafield.or.tdjr1d.or.tdjr2d.or.tdjr3d) then
597  tjr=.true.
598 end if
599 ! allocate current density arrays
600 if (allocated(jrmt)) deallocate(jrmt)
601 if (allocated(jrir)) deallocate(jrir)
602 if (tjr) then
603  allocate(jrmt(npmtmax,natmtot,3),jrir(ngtot,3))
604 end if
605 ! Coulomb potential
606 if (allocated(vclmt)) deallocate(vclmt)
607 allocate(vclmt(npmtmax,natmtot))
608 if (allocated(vclir)) deallocate(vclir)
609 allocate(vclir(ngtot))
610 ! exchange energy density
611 if (allocated(exmt)) deallocate(exmt)
612 allocate(exmt(npmtmax,natmtot))
613 if (allocated(exir)) deallocate(exir)
614 allocate(exir(ngtot))
615 ! correlation energy density
616 if (allocated(ecmt)) deallocate(ecmt)
617 allocate(ecmt(npmtmax,natmtot))
618 if (allocated(ecir)) deallocate(ecir)
619 allocate(ecir(ngtot))
620 ! exchange-correlation potential
621 if (allocated(vxcmt)) deallocate(vxcmt)
622 allocate(vxcmt(npmtmax,natmtot))
623 if (allocated(vxcir)) deallocate(vxcir)
624 allocate(vxcir(ngtot))
625 ! exchange-correlation and dipole magnetic fields
626 if (allocated(bxcmt)) deallocate(bxcmt)
627 if (allocated(bxcir)) deallocate(bxcir)
628 if (allocated(bdmt)) deallocate(bdmt)
629 if (allocated(bdir)) deallocate(bdir)
630 if (allocated(bdmta)) deallocate(bdmta)
631 if (spinpol) then
633  if (tbdip) then
634  allocate(bdmt(npmtmax,natmtot,ndmag),bdir(ngtot,ndmag))
635  allocate(bdmta(ndmag,natmtot))
636  bdmta(1:ndmag,1:natmtot)=0.d0
637  end if
638 end if
639 ! combined target array for Kohn-Sham potential and magnetic field
640 if (allocated(vsbs)) deallocate(vsbs)
642 if (spinpol) n=n+(npcmtmax*natmtot+ngtc)*ndmag
643 allocate(vsbs(n))
644 ! associate pointer arrays with target
645 vsmt(1:npmtmax,1:natmtot) => vsbs(1:)
646 i=size(vsmt)+1
647 vsirc(1:ngtc) => vsbs(i:)
648 if (spinpol) then
649  i=i+size(vsirc)
650  bsmt(1:npcmtmax,1:natmtot,1:ndmag) => vsbs(i:)
651  i=i+size(bsmt)
652  bsirc(1:ngtc,1:ndmag) => vsbs(i:)
653 ! allocate the Kohn-Sham magnetic field on the intersitial grid
654  if (allocated(bsir)) deallocate(bsir)
655  allocate(bsir(ngtot,ndmag))
656 end if
657 ! interstitial Kohn-Sham potential
658 if (allocated(vsir)) deallocate(vsir)
659 allocate(vsir(ngtot))
660 ! interstitial Kohn-Sham potential in G-space
661 if (allocated(vsig)) deallocate(vsig)
662 allocate(vsig(ngvc))
663 ! kinetic energy density and meta-GGA exchange-correlation potential
664 if (allocated(taumt)) deallocate(taumt)
665 if (allocated(tauir)) deallocate(tauir)
666 if (allocated(taucr)) deallocate(taucr)
667 if (allocated(wxcmt)) deallocate(wxcmt)
668 if (allocated(wxcir)) deallocate(wxcir)
669 if (any(xcgrad == [3,4,5,6])) then
671  allocate(taucr(npmtmax,natmtot,nspinor))
672  allocate(wxcmt(npmtmax,natmtot),wxcir(ngtot))
673 ! approximate kinetic energy density functional used to compute the functional
674 ! derivative δτ(r')/δρ(r) for meta-GGA
675  call getxcdata(ktype,kdescr,xcspin_,kgrad,hybrid_,hybridc_)
676 end if
677 ! spin-orbit coupling radial function
678 if (allocated(socfr)) deallocate(socfr)
679 if (spinorb) then
680  allocate(socfr(nrcmtmax,natmtot))
681 end if
682 ! allocate muffin-tin charge and moment arrays
683 if (allocated(chgcrlk)) deallocate(chgcrlk)
684 allocate(chgcrlk(natmtot))
685 if (allocated(chgmt)) deallocate(chgmt)
686 allocate(chgmt(natmtot))
687 if (allocated(mommt)) deallocate(mommt)
688 allocate(mommt(3,natmtot))
689 ! check if scaled spin exchange-correlation should be used
690 tssxc=(abs(sxcscf-1.d0) > 1.d-6)
691 ! spin-spiral phase factors
692 if (ssdph) then
693  if (allocated(zqss)) deallocate(zqss)
694  allocate(zqss(natmtot))
695  do ias=1,natmtot
696  is=idxis(ias)
697  ia=idxia(ias)
698  t1=-0.5d0*dot_product(vqcss(1:3),atposc(1:3,ia,is))
699  zqss(ias)=cmplx(cos(t1),sin(t1),8)
700  end do
701 end if
702 ! mixing vector: either density/magnetisation or potential/magnetic field
703 if (mixrho) then
704  vmixer => rhmg
705 else
706  vmixer => vsbs
707 end if
708 ! zero the mixing vector
709 vmixer(:)=0.d0
710 
711 !-------------------------!
712 ! force variables !
713 !-------------------------!
714 if (tforce) then
715  if (allocated(forcehf)) deallocate(forcehf)
716  allocate(forcehf(3,natmtot))
717  if (allocated(forcetot)) deallocate(forcetot)
718  allocate(forcetot(3,natmtot))
719 end if
720 
721 !-------------------------------------------------!
722 ! DFT+U and fixed tensor moment variables !
723 !-------------------------------------------------!
724 if ((dftu /= 0).or.(ftmtype /= 0)) then
725 ! density matrix elements in each muffin-tin
726  if (allocated(dmatmt)) deallocate(dmatmt)
728 ! potential matrix elements in each muffin-tin
729  if (allocated(vmatmt)) deallocate(vmatmt)
731 ! zero the potential matrix
732  vmatmt(:,:,:,:,:)=0.d0
733 ! matrix elements in spherical coordinates for TDDFT+U
734  if (any(task == [460,461,462,463,478])) then
735  if (allocated(vmatmti)) deallocate(vmatmti)
737  if (allocated(vmatmto)) deallocate(vmatmto)
739  end if
740 ! require the potential matrix elements be calculated
741  tvmatmt=.true.
742 ! flags for non-zero muffin-tin potential matrices
743  if (allocated(tvmmt)) deallocate(tvmmt)
744  allocate(tvmmt(0:lmaxdm,natmtot))
745  tvmmt(:,:)=.false.
746 ! require second-variational eigenvectors
747  tevecsv=.true.
748 end if
749 if (dftu /= 0) then
750  if (any(task == [5,300,600,601,610,620,630,640])) then
751  write(*,*)
752  write(*,'("Error(init0): DFT+U does not work with task ",I0)') task
753  write(*,*)
754  stop
755  end if
756 ! DFT+U energy for each atom
757  if (allocated(engyadu)) deallocate(engyadu)
758  allocate(engyadu(natmmax,ndftu))
759 ! flag the muffin-tin potential matrices which are non-zero
760  do idu=1,ndftu
761  is=isldu(1,idu)
762  if (is > nspecies) then
763  write(*,*)
764  write(*,'("Error(init0): invalid species number : ",I0)') is
765  write(*,*)
766  stop
767  end if
768  l=isldu(2,idu)
769  do ia=1,natoms(is)
770  ias=idxas(ia,is)
771  tvmmt(l,ias)=.true.
772  end do
773  end do
774 ! zero the initial values of screening length
775  lamdu0(:)=0.d0
776 ! write to VARIABLES.OUT
777  if (wrtvars) then
778  call writevars('udufix',nv=ndftu,rva=udufix)
779  end if
780 end if
781 if (ftmtype /= 0) then
782 ! allocate and zero the fixed tensor moment potential array
783  if (allocated(vmftm)) deallocate(vmftm)
785  vmftm(:,:,:,:,:)=0.d0
786 ! flag the muffin-tin potential matrices which are non-zero
787  do i=1,ntmfix
788  is=itmfix(1,i)
789  ia=itmfix(2,i)
790  ias=idxas(ia,is)
791  l=itmfix(3,i)
792  tvmmt(l,ias)=.true.
793  end do
794 end if
795 
796 !-----------------------!
797 ! miscellaneous !
798 !-----------------------!
799 ! determine nuclear radii and volumes
800 call nuclei
801 ! determine the nuclear-nuclear energy
802 call energynn
803 ! get smearing function description
804 call getsdata(stype,sdescr)
805 ! get mixing type description
807 ! generate the spherical harmonic transform (SHT) matrices
808 call genshtmat
809 ! find the maximum size of the spherical Bessel function array over all species
810 call findnjcmax
811 ! allocate 1D plotting arrays
812 if (allocated(dvp1d)) deallocate(dvp1d)
813 allocate(dvp1d(nvp1d))
814 if (allocated(vplp1d)) deallocate(vplp1d)
815 allocate(vplp1d(3,npp1d))
816 if (allocated(dpp1d)) deallocate(dpp1d)
817 allocate(dpp1d(npp1d))
818 ! initial self-consistent loop number
819 iscl=1
820 tlast=.false.
821 ! set the Fermi energy to zero
822 efermi=0.d0
823 ! set the temperature from the smearing width
825 
826 call timesec(ts1)
827 timeinit=timeinit+ts1-ts0
828 
829 end subroutine
830 !EOC
831 
integer, dimension(:,:), allocatable itmfix
Definition: moddftu.f90:76
subroutine getsdata(stype, sdescr)
Definition: sdelta.f90:68
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition: modmain.f90:275
real(8), dimension(:), allocatable wxcir
Definition: modmain.f90:676
real(8) efermi
Definition: modmain.f90:905
subroutine potnucl(ptnucl, nr, r, zn, vn)
Definition: potnucl.f90:10
real(8), dimension(maxstsp, maxspecies) occsp
Definition: modmain.f90:133
integer npcmttot
Definition: modmain.f90:218
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
real(8), dimension(3, 3) afspc
Definition: modmain.f90:331
integer, dimension(maxstsp, maxspecies) ksp
Definition: modmain.f90:125
real(8) dmaxefc
Definition: modmain.f90:318
real(8), dimension(:,:), allocatable mommt
Definition: modmain.f90:744
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
logical tjr
Definition: modmain.f90:620
real(8), dimension(3) efieldl
Definition: modmain.f90:314
integer mixtype
Definition: modmain.f90:695
integer, dimension(3) ktype
Definition: modmain.f90:606
integer natmmax
Definition: modmain.f90:38
integer npcmtmax
Definition: modmain.f90:216
integer, dimension(3) ngridg
Definition: modmain.f90:386
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
subroutine gengvc
Definition: gengvc.f90:7
real(8), dimension(:), allocatable, target rhmg
Definition: modmain.f90:612
integer task
Definition: modmain.f90:1297
logical mp_mpi
Definition: modmpi.f90:17
subroutine reciplat(avec, bvec, omega, omegabz)
Definition: reciplat.f90:10
integer ngtot
Definition: modmain.f90:390
real(8), dimension(:), pointer, contiguous vmixer
Definition: modmain.f90:689
integer lmmaxo
Definition: modmain.f90:203
complex(4), dimension(:,:,:,:,:), allocatable vmatmti
Definition: moddftu.f90:22
real(8), dimension(:,:), allocatable vxcmt
Definition: modmain.f90:634
integer, dimension(3) xctype
Definition: modmain.f90:588
integer ngtc
Definition: modmain.f90:392
real(8) rwigner
Definition: modmain.f90:736
real(8), dimension(:,:), allocatable occcr
Definition: modmain.f90:933
real(8), dimension(:), allocatable ecir
Definition: modmain.f90:632
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer lmmaxapw
Definition: modmain.f90:199
logical hybrid
Definition: modmain.f90:1150
integer, parameter lmmaxdm
Definition: moddftu.f90:14
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
complex(8), dimension(:,:,:,:,:), allocatable dmftm
Definition: moddftu.f90:80
logical tshift
Definition: modmain.f90:352
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
logical tefvit
Definition: modmain.f90:873
real(8) reducebf
Definition: modmain.f90:279
integer ndmag
Definition: modmain.f90:238
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8) sxcscf
Definition: modmain.f90:668
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
real(8) omega
Definition: modmain.f90:20
logical tevecsv
Definition: modmain.f90:921
real(8), dimension(:), allocatable dpp1d
Definition: modmain.f90:1124
integer iscl
Definition: modmain.f90:1049
integer, dimension(:), allocatable ipcmt
Definition: modmain.f90:220
integer nstcr
Definition: modmain.f90:129
subroutine r3minv(a, b)
Definition: r3minv.f90:10
real(8), dimension(:,:,:), allocatable bdmt
Definition: modmain.f90:638
logical spinsprl
Definition: modmain.f90:283
Definition: modomp.f90:6
real(8), parameter kboltz
Definition: modmain.f90:1261
real(8) swidth
Definition: modmain.f90:893
real(8), dimension(:,:), allocatable evalcr
Definition: modmain.f90:935
subroutine gendmftm
Definition: gendmftm.f90:7
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
logical hybrid0
Definition: modmain.f90:1150
real(8) epsocc
Definition: modmain.f90:901
real(8), dimension(:,:), allocatable bdmta
Definition: modmain.f90:640
subroutine getmixdata(mtype, mixdescr)
Definition: mixerifc.f90:47
real(8), dimension(:,:,:), allocatable rhocr
Definition: modmain.f90:939
subroutine gencfun
Definition: gencfun.f90:10
logical tssxc
Definition: modmain.f90:666
real(8), dimension(:), allocatable chgmt
Definition: modmain.f90:732
character(264) kdescr
Definition: modmain.f90:608
real(8) timemat
Definition: modmain.f90:1215
logical mixrho
Definition: modmain.f90:687
character(264) xcdescr
Definition: modmain.f90:590
real(8), dimension(:,:), allocatable ecmt
Definition: modmain.f90:632
logical nosource
Definition: modmain.f90:664
real(8), dimension(:,:), allocatable forcetot
Definition: modmain.f90:991
integer lmaxo
Definition: modmain.f90:201
real(8), dimension(3) vqlss
Definition: modmain.f90:293
logical tmwrite
Definition: moddftu.f90:66
logical tafieldt
Definition: modtddft.f90:56
real(8), dimension(3) vecql
Definition: modmain.f90:1102
integer nrspmax
Definition: modmain.f90:109
real(8), dimension(3) bfsmc
Definition: modmain.f90:257
integer ngvc
Definition: modmain.f90:398
real(8) timefor
Definition: modmain.f90:1225
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition: moddftu.f90:16
real(8) vmaxefc
Definition: modmain.f90:320
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
integer xcspin
Definition: modmain.f90:592
integer lmaxapw
Definition: modmain.f90:197
real(8), dimension(:,:), allocatable exmt
Definition: modmain.f90:630
integer kgrad
Definition: modmain.f90:610
real(8) timerho
Definition: modmain.f90:1221
real(8), dimension(:), allocatable vxcir
Definition: modmain.f90:634
logical bforb
Definition: modmain.f90:234
real(8) chgcrtot
Definition: modmain.f90:718
integer lmmaxdb
Definition: modmain.f90:1079
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
complex(4), dimension(:,:,:,:,:), allocatable vmatmto
Definition: moddftu.f90:22
logical tforce
Definition: modmain.f90:987
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:40
real(8), dimension(3, 3) avecref
Definition: modmain.f90:1028
logical, dimension(maxstsp, maxspecies) spcore
Definition: modmain.f90:127
real(8) timeinit
Definition: modmain.f90:1213
integer nrcmtmax
Definition: modmain.f90:175
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
real(8) timepot
Definition: modmain.f90:1223
logical tefield
Definition: modmain.f90:310
real(8), dimension(:,:,:), allocatable taucr
Definition: modmain.f90:674
complex(8), dimension(:), allocatable zqss
Definition: modmain.f90:287
subroutine symmetry
Definition: symmetry.f90:7
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
logical tlast
Definition: modmain.f90:1051
real(8), dimension(:), allocatable dvp1d
Definition: modmain.f90:1120
integer nspncr
Definition: modmain.f90:943
subroutine gencfrc
Definition: gencfrc.f90:7
integer ndftu
Definition: moddftu.f90:38
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
real(8), dimension(:,:,:), allocatable bxcmt
Definition: modmain.f90:636
real(8), dimension(:,:), allocatable bsir
Definition: modmain.f90:658
subroutine getxcdata(xctype, xcdescr, xcspin, xcgrad, hybrid, hybridc)
Definition: modxcifc.f90:411
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
real(8) omegamt
Definition: modmain.f90:169
real(8) occmax
Definition: modmain.f90:899
subroutine readafieldt
Definition: readafieldt.f90:7
real(8), dimension(maxspecies) chgcr
Definition: modmain.f90:716
integer ftmtype
Definition: moddftu.f90:70
integer nbph
Definition: modphonon.f90:13
subroutine energynn
Definition: energynn.f90:7
integer lmaxdb
Definition: modmain.f90:1079
real(8), dimension(maxstsp, maxspecies) evalsp
Definition: modmain.f90:131
real(8), dimension(:,:), allocatable engyadu
Definition: moddftu.f90:44
real(8), dimension(:), allocatable chgcrlk
Definition: modmain.f90:720
subroutine genrmesh
Definition: genrmesh.f90:10
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
subroutine genshtmat
Definition: genshtmat.f90:10
real(8), dimension(3, 0:1) afindt
Definition: modtddft.f90:67
logical tdjr2d
Definition: modtddft.f90:93
real(8) tempk
Definition: modmain.f90:684
real(8), dimension(:,:), allocatable bdir
Definition: modmain.f90:638
real(8) rgkmax
Definition: modmain.f90:493
real(8), dimension(:,:), allocatable ffacg
Definition: modmain.f90:432
logical cmagz
Definition: modmain.f90:242
integer ngvec
Definition: modmain.f90:396
integer nscfxc
Definition: modtddft.f90:16
real(8) chgzn
Definition: modmain.f90:714
integer lradstp
Definition: modmain.f90:171
logical tvmatmt
Definition: moddftu.f90:24
real(8), dimension(:,:), allocatable forcehf
Definition: modmain.f90:989
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
integer nspinor
Definition: modmain.f90:267
subroutine gengvec
Definition: gengvec.f90:10
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:20
real(8), dimension(3) afieldc
Definition: modmain.f90:325
logical, dimension(:,:), allocatable tvmmt
Definition: moddftu.f90:26
pure subroutine r3frac(eps, v)
Definition: r3frac.f90:10
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
logical tavref
Definition: modmain.f90:1030
subroutine allatoms
Definition: allatoms.f90:10
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
real(8), dimension(3, 3) avec0
Definition: modmain.f90:12
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
complex(8), dimension(:,:,:,:,:), allocatable vmftm
Definition: moddftu.f90:82
subroutine checkmt
Definition: checkmt.f90:10
real(8), dimension(maxspecies) spze
Definition: modmain.f90:99
real(8), dimension(:), allocatable, target vsbs
Definition: modmain.f90:647
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
real(8), dimension(3, maxatoms, maxspecies) bfcmt
Definition: modmain.f90:273
real(8), dimension(:), pointer, contiguous vsirc
Definition: modmain.f90:653
character(64) mixdescr
Definition: modmain.f90:697
real(8), dimension(3) vqcss
Definition: modmain.f90:295
subroutine strainabg
Definition: strainabg.f90:7
Definition: modgw.f90:6
integer, dimension(maxspecies) nrsp
Definition: modmain.f90:107
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer, dimension(2, 3) intgv
Definition: modmain.f90:394
logical wrtvars
Definition: modvars.f90:9
logical ssdph
Definition: modmain.f90:285
integer, dimension(maxspecies) npcmti
Definition: modmain.f90:214
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(3, 3) binv
Definition: modmain.f90:18
real(8) epslat
Definition: modmain.f90:24
integer npsd
Definition: modmain.f90:626
real(8), parameter y00i
Definition: modmain.f90:1235
real(8), dimension(maxdftu) udufix
Definition: moddftu.f90:60
integer stype
Definition: modmain.f90:889
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
character(64), dimension(maxspecies) spname
Definition: modmain.f90:76
real(8), dimension(:,:), allocatable bfsmcmt
Definition: modmain.f90:263
subroutine timesec(ts)
Definition: timesec.f90:10
character(64) sdescr
Definition: modmain.f90:891
integer dftu
Definition: moddftu.f90:32
real(8) chgval
Definition: modmain.f90:722
logical ptnucl
Definition: modmain.f90:83
real(8) timefv
Definition: modmain.f90:1217
subroutine nuclei
Definition: nuclei.f90:7
pure subroutine genffacgp(ngp, gpc, ld, ffacgp)
Definition: genffacgp.f90:10
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
logical tbdip
Definition: modmain.f90:643
logical spinorb
Definition: modmain.f90:230
real(8), dimension(:,:), allocatable tauir
Definition: modmain.f90:672
real(8), dimension(3) bfieldc0
Definition: modmain.f90:271
real(8), dimension(maxdftu) lamdu0
Definition: moddftu.f90:54
integer lnpsd
Definition: modmain.f90:628
integer nspecies
Definition: modmain.f90:34
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
real(8), dimension(:,:), allocatable efcmt
Definition: modmain.f90:316
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
Definition: genjlgprmt.f90:10
real(8) chgtot
Definition: modmain.f90:726
logical tdjr3d
Definition: modtddft.f90:93
real(8) timesv
Definition: modmain.f90:1219
integer, dimension(maxspecies) nstsp
Definition: modmain.f90:113
real(8) gkmax
Definition: modmain.f90:495
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition: modmain.f90:937
logical tafsp
Definition: modmain.f90:329
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer npmtmax
Definition: modmain.f90:216
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
real(8), dimension(:,:,:), allocatable jrmt
Definition: modmain.f90:622
subroutine init0
Definition: init0.f90:10
real(8), parameter fourpi
Definition: modmain.f90:1232
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer natmtot
Definition: modmain.f90:40
real(8) omegabz
Definition: modmain.f90:22
logical ksgwrho
Definition: modgw.f90:38
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, dimension(maxspecies) npmti
Definition: modmain.f90:213
integer, parameter lmaxdm
Definition: moddftu.f90:13
logical tdjr1d
Definition: modtddft.f90:93
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
Definition: modulr.f90:6
integer nrmtmax
Definition: modmain.f90:152
subroutine genylmg
Definition: genylmg.f90:10
logical ncmag
Definition: modmain.f90:240
subroutine gengclg
Definition: gengclg.f90:7
real(8) chgexs
Definition: modmain.f90:724
complex(8), dimension(:), allocatable vsig
Definition: modmain.f90:662
logical tafield
Definition: modmain.f90:322
real(8), dimension(:,:), pointer, contiguous bsirc
Definition: modmain.f90:660
integer nvp1d
Definition: modmain.f90:1112
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
integer nstspmax
Definition: modmain.f90:115
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer npp1d
Definition: modmain.f90:1114
real(8), dimension(3) efieldc
Definition: modmain.f90:312
real(8), dimension(3) afieldl
Definition: modmain.f90:327
real(8), dimension(:,:,:), allocatable taumt
Definition: modmain.f90:672
subroutine findnjcmax
Definition: findnjcmax.f90:7
real(8), dimension(:,:), allocatable jrir
Definition: modmain.f90:622
real(8) hybridc
Definition: modmain.f90:1152
integer xcgrad
Definition: modmain.f90:602
real(8), dimension(:,:), allocatable vplp1d
Definition: modmain.f90:1122
integer, dimension(2) jspnfv
Definition: modmain.f90:291
logical spincore
Definition: modmain.f90:941
character(64), dimension(maxspecies) spsymb
Definition: modmain.f90:78
real(8), dimension(:), allocatable exir
Definition: modmain.f90:630
real(8), dimension(3) vecqc
Definition: modmain.f90:1102
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition: modvars.f90:16
integer isgkmax
Definition: modmain.f90:491
integer nspnfv
Definition: modmain.f90:289
real(8), dimension(3) bfieldc
Definition: modmain.f90:269
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition: modmain.f90:656
integer fsmtype
Definition: modmain.f90:251
real(8), dimension(:,:), allocatable socfr
Definition: modmain.f90:670
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
real(8) gmaxvr
Definition: modmain.f90:384
integer ntmfix
Definition: moddftu.f90:72
integer lmaxi
Definition: modmain.f90:205