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 : ",2I8)') 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 ",I4)') 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 (mp_mpi.and.(mixtype /= 1)) then
227  write(*,'("Info(init0): mixtype changed to 1 for FSM calculation")')
228  end if
229  mixtype=1
230 end if
231 ! number of independent spin components of the f_xc spin tensor
232 if (spinpol) then
233  if (ncmag) then
234  nscfxc=10
235  else
236  nscfxc=3
237  end if
238 else
239  nscfxc=1
240 end if
241 ! set the magnetic fields to the initial values
242 bfieldc(:)=bfieldc0(:)
243 bfcmt(:,:,:)=bfcmt0(:,:,:)
244 if (tmwrite.or.(ftmtype /= 0).or.(task == 400)) then
245  if (.not.spinorb) then
246  write(*,*)
247  write(*,'("Error(init0): tensor moments require spin-orbit coupling &
248  &enabled")')
249  write(*,'(" set spinorb=.true.")')
250  write(*,*)
251  stop
252  end if
253 end if
254 ! generate the fixed tensor moment density matrices if required
255 if (ftmtype /= 0) call gendmftm
256 ! if reducebf < 1 then reduce the external magnetic fields and tensor moments
257 ! immediately for non-self-consistent calculations
258 if (reducebf < 1.d0-1.d-4) then
259  if (all(task /= [0,1,2,3,28,200,201,208,350,351,380,420,421,440])) then
260  bfieldc(:)=0.d0
261  bfcmt(:,:,:)=0.d0
262  if (ftmtype < 0) dmftm(:,:,:,:,:)=0.d0
263  end if
264 end if
265 ! write to VARIABLES.OUT
266 if (wrtvars) then
267  call writevars('nspinor',iv=nspinor)
268  call writevars('ndmag',iv=ndmag)
269  call writevars('sxcscf',rv=sxcscf)
270  call writevars('bfieldc0',nv=3,rva=bfieldc0)
271 end if
272 
273 !----------------------------------!
274 ! crystal structure set up !
275 !----------------------------------!
276 ! use reference lattice vectors if required
277 if (tavref) then
278  avec0(:,:)=avec(:,:)
279  avec(:,:)=avecref(:,:)
280 end if
281 ! generate the reciprocal lattice vectors, unit cell and Brillouin zone volumes
283 ! inverse of the lattice vector matrix
284 call r3minv(avec,ainv)
285 ! inverse of the reciprocal vector matrix
286 call r3minv(bvec,binv)
287 ! Cartesian coordinates of the spin-spiral vector
288 call r3mv(bvec,vqlss,vqcss)
289 do is=1,nspecies
290  do ia=1,natoms(is)
291 ! map atomic lattice coordinates to [0,1)
292  call r3frac(epslat,atposl(:,ia,is))
293 ! determine atomic Cartesian coordinates
294  call r3mv(avec,atposl(:,ia,is),atposc(:,ia,is))
295  end do
296 end do
297 ! check for overlapping muffin-tins and adjust radii if required
298 call checkmt
299 ! compute the total muffin-tin volume (M. Meinert)
300 omegamt=0.d0
301 do is=1,nspecies
302  omegamt=omegamt+dble(natoms(is))*(fourpi/3.d0)*rmt(is)**3
303 end do
304 ! input q-vector in Cartesian coordinates
305 call r3mv(bvec,vecql,vecqc)
306 
307 !-------------------------------!
308 ! vector fields E and A !
309 !-------------------------------!
310 ! static electric field
311 tefield=.false.
312 if (any(abs(efieldc(:)) > epslat)) then
313 ! no shift of the atomic positions
314  tshift=.false.
315 ! electric field vector in lattice coordinates
316  call r3mv(ainv,efieldc,efieldl)
317 ! potential at maximum distance
318  vmaxefc=dmaxefc*norm2(efieldc(1:3))
319 ! allocate array for average electric field in each muffin-tin
320  if (allocated(efcmt)) deallocate(efcmt)
321  allocate(efcmt(3,natmtot))
322 ! set the E-field flag
323  tefield=.true.
324 end if
325 ! static vector potential
326 tafield=.false.
327 if (any(abs(afieldc(:)) > epslat)) then
328  tafield=.true.
329 ! A-field in lattice coordinates
330  call r3mv(ainv,afieldc,afieldl)
331 ! vector potential added in second-variational step
332  tevecsv=.true.
333 end if
334 ! static spin-dependent vector potential
335 tafsp=.false.
336 if (spinpol.and.(any(abs(afspc(:,:)) > epslat))) then
337  tafsp=.true.
338  tevecsv=.true.
339 end if
340 ! time-dependent vector potential
341 tafieldt=.false.
342 if (any(task == [460,461,462,463,480,481,485])) then
343 ! read time-dependent A-field from file
344  call readafieldt
345  tafieldt=.true.
346 ! zero the induced A-field and its time derivative
347  afindt(:,:)=0.d0
348 end if
349 ! write to VARIABLES.OUT
350 if (wrtvars) then
351  call writevars('efieldc',nv=3,rva=efieldc)
352  call writevars('afieldc',nv=3,rva=afieldc)
353  call writevars('afspc',nv=9,rva=afspc)
354 end if
355 
356 !---------------------------------!
357 ! crystal symmetry set up !
358 !---------------------------------!
359 call symmetry
360 
361 !-----------------------!
362 ! radial meshes !
363 !-----------------------!
364 nrmtmax=1
365 nrcmtmax=1
366 do is=1,nspecies
367 ! make the muffin-tin mesh commensurate with lradstp
368  nrmt(is)=nrmt(is)-mod(nrmt(is)-1,lradstp)
369  nrmtmax=max(nrmtmax,nrmt(is))
370 ! number of coarse radial mesh points
371  nrcmt(is)=(nrmt(is)-1)/lradstp+1
372  nrcmtmax=max(nrcmtmax,nrcmt(is))
373 end do
374 ! set up atomic and muffin-tin radial meshes
375 call genrmesh
376 ! number of points in packed muffin-tins
377 npmtmax=1
378 npcmtmax=1
379 do is=1,nspecies
380  npmti(is)=lmmaxi*nrmti(is)
381  npmt(is)=npmti(is)+lmmaxo*(nrmt(is)-nrmti(is))
382  npmtmax=max(npmtmax,npmt(is))
383  npcmti(is)=lmmaxi*nrcmti(is)
384  npcmt(is)=npcmti(is)+lmmaxo*(nrcmt(is)-nrcmti(is))
385  npcmtmax=max(npcmtmax,npcmt(is))
386 end do
387 ! index to first muffin-tin point in packed array over all atoms
388 if (allocated(ipcmt)) deallocate(ipcmt)
389 allocate(ipcmt(natmtot))
390 ipcmt(1)=1
391 npc=npcmt(1)
392 npcmttot=npc
393 do ias=2,natmtot
394  is=idxis(ias)
395  ipcmt(ias)=ipcmt(ias)+npc
396  npc=npcmt(is)
397  npcmttot=npcmttot+npc
398 end do
399 
400 !--------------------------------------!
401 ! charges and number of states !
402 !--------------------------------------!
403 chgzn=0.d0
404 chgcrtot=0.d0
405 chgval=0.d0
406 nstspmax=0
407 nstcr=0
408 do is=1,nspecies
409 ! nuclear charge
410  chgzn=chgzn+spzn(is)*natoms(is)
411 ! find the maximum number of atomic states
412  nstspmax=max(nstspmax,nstsp(is))
413 ! compute the electronic charge for each species, as well as the total core and
414 ! valence charge
415  spze(is)=0.d0
416  chgcr(is)=0.d0
417  do ist=1,nstsp(is)
418  spze(is)=spze(is)+occsp(ist,is)
419  if (spcore(ist,is)) then
420  chgcr(is)=chgcr(is)+occsp(ist,is)
421  nstcr=nstcr+2*ksp(ist,is)*natoms(is)
422  else
423  chgval=chgval+occsp(ist,is)*natoms(is)
424  end if
425  end do
426  chgcrtot=chgcrtot+chgcr(is)*natoms(is)
427 end do
428 ! add excess charge
430 ! total charge
432 if (chgtot < 1.d-8) then
433  write(*,*)
434  write(*,'("Error(init0): zero total charge")')
435  write(*,*)
436  stop
437 end if
438 ! subtract small charge to ensure consistent Fermi energy for insulators
440 ! effective Wigner radius
441 rwigner=(3.d0/(fourpi*(chgtot/omega)))**(1.d0/3.d0)
442 ! write to VARIABLES.OUT
443 if (wrtvars) then
444  call writevars('spze',nv=nspecies,rva=spze)
445  call writevars('chgcr',nv=nspecies,rva=chgcr)
446  call writevars('chgexs',rv=chgexs)
447  call writevars('chgval',rv=chgtot)
448 end if
449 
450 !-------------------------!
451 ! G-vector arrays !
452 !-------------------------!
453 ! determine gkmax from rgkmax
454 if (nspecies == 0) isgkmax=-2
455 select case(isgkmax)
456 case(:-4)
457 ! use largest muffin-tin radius
458  gkmax=rgkmax/maxval(rmt(1:nspecies))
459 case(-3)
460 ! use smallest muffin-tin radius
461  gkmax=rgkmax/minval(rmt(1:nspecies))
462 case(-2)
463 ! use the fixed value of 2.0
464  gkmax=rgkmax/2.d0
465 case(-1,0)
466 ! use average muffin-tin radius
467  t1=sum(natoms(1:nspecies)*rmt(1:nspecies))/dble(natmtot)
468  gkmax=rgkmax/t1
469 case(1:)
470 ! use user-specified muffin-tin radius
471  if (isgkmax <= nspecies) then
473  else
474  write(*,*)
475  write(*,'("Error(init0): isgkmax > nspecies : ",2I8)') isgkmax,nspecies
476  write(*,*)
477  stop
478  end if
479 end select
480 ! generate the G-vectors
481 call gengvec
482 ! apply strain to A-, B- and G-vectors if required
483 call strainabg
484 ! write number of G-vectors to test file
485 call writetest(900,'number of G-vectors',iv=ngvec)
486 ! Poisson solver pseudocharge density constant
487 if (nspecies > 0) then
488  t1=0.25d0*gmaxvr*maxval(rmt(1:nspecies))
489 else
490  t1=0.25d0*gmaxvr*2.d0
491 end if
492 npsd=max(nint(t1),1)
493 lnpsd=lmaxo+npsd+1
494 ! generate the Coulomb Green's function in G-space = 4π/G²
495 call gengclg
496 ! compute the spherical Bessel functions j_l(|G|R_mt)
497 if (allocated(jlgrmt)) deallocate(jlgrmt)
498 allocate(jlgrmt(0:lnpsd,ngvec,nspecies))
500 ! generate the spherical harmonics of the G-vectors
501 call genylmg
502 ! allocate structure factor array for G-vectors
503 if (allocated(sfacg)) deallocate(sfacg)
504 allocate(sfacg(ngvec,natmtot))
505 ! generate structure factors for G-vectors
507 ! generate the smooth step function form factors
508 if (allocated(ffacg)) deallocate(ffacg)
509 allocate(ffacg(ngtot,nspecies))
511 ! generate the smooth characteristic function
512 call gencfun
513 ! G-vector variables for coarse grid with |G| < 2 gkmax
514 call gengvc
515 ! generate the characteristic function on the coarse grid
516 call gencfrc
517 ! write to VARIABLES.OUT
518 if (wrtvars) then
519  call writevars('avec',nv=9,rva=avec)
520  call writevars('bvec',nv=9,rva=bvec)
521  call writevars('omega',rv=omega)
522  do is=1,nspecies
523  call writevars('atposl',n1=is,nv=3*natoms(is),rva=atposl(:,:,is))
524  end do
525  do is=1,nspecies
526  call writevars('atposc',n1=is,nv=3*natoms(is),rva=atposc(:,:,is))
527  end do
528  call writevars('vqlss',nv=3,rva=vqlss)
529  call writevars('vqcss',nv=3,rva=vqcss)
530  call writevars('gmaxvr',rv=gmaxvr)
531  call writevars('ngridg',nv=3,iva=ngridg)
532  call writevars('intgv',nv=6,iva=intgv)
533  call writevars('ngvec',iv=ngvec)
534  call writevars('ivg',nv=3*ngtot,iva=ivg)
535  call writevars('igfft',nv=ngtot,iva=igfft)
536 end if
537 
538 !-------------------------!
539 ! atoms and cores !
540 !-------------------------!
541 ! determine the nuclear Coulomb potential
542 if (allocated(vcln)) deallocate(vcln)
543 allocate(vcln(nrspmax,nspecies))
544 do is=1,nspecies
545  nr=nrsp(is)
546  call potnucl(ptnucl,nr,rsp(:,is),spzn(is),vcln(:,is))
547  vcln(1:nr,is)=vcln(1:nr,is)*y00i
548 end do
549 ! solve the Kohn-Sham-Dirac equations for all atoms
550 call allatoms
551 ! allocate core state occupancy and eigenvalue arrays and set to default
552 if (allocated(occcr)) deallocate(occcr)
553 allocate(occcr(nstspmax,natmtot))
554 if (allocated(evalcr)) deallocate(evalcr)
555 allocate(evalcr(nstspmax,natmtot))
556 do ias=1,natmtot
557  is=idxis(ias)
558  do ist=1,nstsp(is)
559  occcr(ist,ias)=occsp(ist,is)
560  evalcr(ist,ias)=evalsp(ist,is)
561  end do
562 end do
563 ! allocate core state radial wavefunction array
564 if (allocated(rwfcr)) deallocate(rwfcr)
565 allocate(rwfcr(nrspmax,2,nstspmax,natmtot))
566 ! number of core spin channels
567 if (spincore) then
568  nspncr=2
569 else
570  nspncr=1
571 end if
572 ! allocate core state charge density array
573 if (allocated(rhocr)) deallocate(rhocr)
574 allocate(rhocr(nrmtmax,natmtot,nspncr))
575 
576 !-------------------------------------------------------------!
577 ! charge density, potentials and exchange-correlation !
578 !-------------------------------------------------------------!
579 ! combined target array for density and magnetisation
580 if (allocated(rhmg)) deallocate(rhmg)
582 if (spinpol) n=n*(1+ndmag)
583 allocate(rhmg(n))
584 ! associate pointer arrays with target
585 rhomt(1:npmtmax,1:natmtot) => rhmg(1:)
586 i=size(rhomt)+1
587 rhoir(1:ngtot) => rhmg(i:)
588 if (spinpol) then
589  i=i+size(rhoir)
590  magmt(1:npmtmax,1:natmtot,1:ndmag) => rhmg(i:)
591  i=i+size(magmt)
592  magir(1:ngtot,1:ndmag) => rhmg(i:)
593 end if
594 if (any(task == [371,372,373]).or.tafield.or.tdjr1d.or.tdjr2d.or.tdjr3d) then
595  tjr=.true.
596 end if
597 ! allocate current density arrays
598 if (allocated(jrmt)) deallocate(jrmt)
599 if (allocated(jrir)) deallocate(jrir)
600 if (tjr) then
601  allocate(jrmt(npmtmax,natmtot,3),jrir(ngtot,3))
602 end if
603 ! Coulomb potential
604 if (allocated(vclmt)) deallocate(vclmt)
605 allocate(vclmt(npmtmax,natmtot))
606 if (allocated(vclir)) deallocate(vclir)
607 allocate(vclir(ngtot))
608 ! exchange energy density
609 if (allocated(exmt)) deallocate(exmt)
610 allocate(exmt(npmtmax,natmtot))
611 if (allocated(exir)) deallocate(exir)
612 allocate(exir(ngtot))
613 ! correlation energy density
614 if (allocated(ecmt)) deallocate(ecmt)
615 allocate(ecmt(npmtmax,natmtot))
616 if (allocated(ecir)) deallocate(ecir)
617 allocate(ecir(ngtot))
618 ! exchange-correlation potential
619 if (allocated(vxcmt)) deallocate(vxcmt)
620 allocate(vxcmt(npmtmax,natmtot))
621 if (allocated(vxcir)) deallocate(vxcir)
622 allocate(vxcir(ngtot))
623 ! exchange-correlation and dipole magnetic fields
624 if (allocated(bxcmt)) deallocate(bxcmt)
625 if (allocated(bxcir)) deallocate(bxcir)
626 if (allocated(bdmt)) deallocate(bdmt)
627 if (allocated(bdir)) deallocate(bdir)
628 if (allocated(bdmta)) deallocate(bdmta)
629 if (spinpol) then
631  if (tbdip) then
632  allocate(bdmt(npmtmax,natmtot,ndmag),bdir(ngtot,ndmag))
633  allocate(bdmta(ndmag,natmtot))
634  bdmta(1:ndmag,1:natmtot)=0.d0
635  end if
636 end if
637 ! combined target array for Kohn-Sham potential and magnetic field
638 if (allocated(vsbs)) deallocate(vsbs)
640 if (spinpol) n=n+(npcmtmax*natmtot+ngtc)*ndmag
641 allocate(vsbs(n))
642 ! associate pointer arrays with target
643 vsmt(1:npmtmax,1:natmtot) => vsbs(1:)
644 i=size(vsmt)+1
645 vsirc(1:ngtc) => vsbs(i:)
646 if (spinpol) then
647  i=i+size(vsirc)
648  bsmt(1:npcmtmax,1:natmtot,1:ndmag) => vsbs(i:)
649  i=i+size(bsmt)
650  bsirc(1:ngtc,1:ndmag) => vsbs(i:)
651 ! allocate the Kohn-Sham magnetic field on the intersitial grid
652  if (allocated(bsir)) deallocate(bsir)
653  allocate(bsir(ngtot,ndmag))
654 end if
655 ! interstitial Kohn-Sham potential
656 if (allocated(vsir)) deallocate(vsir)
657 allocate(vsir(ngtot))
658 ! interstitial Kohn-Sham potential in G-space
659 if (allocated(vsig)) deallocate(vsig)
660 allocate(vsig(ngvc))
661 ! kinetic energy density and meta-GGA exchange-correlation potential
662 if (allocated(taumt)) deallocate(taumt)
663 if (allocated(tauir)) deallocate(tauir)
664 if (allocated(taucr)) deallocate(taucr)
665 if (allocated(wxcmt)) deallocate(wxcmt)
666 if (allocated(wxcir)) deallocate(wxcir)
667 if (any(xcgrad == [3,4,5,6])) then
669  allocate(taucr(npmtmax,natmtot,nspinor))
670  allocate(wxcmt(npmtmax,natmtot),wxcir(ngtot))
671 ! approximate kinetic energy density functional used to compute the functional
672 ! derivative δτ(r')/δρ(r) for meta-GGA
673  call getxcdata(ktype,kdescr,xcspin_,kgrad,hybrid_,hybridc_)
674 end if
675 ! spin-orbit coupling radial function
676 if (allocated(socfr)) deallocate(socfr)
677 if (spinorb) then
678  allocate(socfr(nrcmtmax,natmtot))
679 end if
680 ! allocate muffin-tin charge and moment arrays
681 if (allocated(chgcrlk)) deallocate(chgcrlk)
682 allocate(chgcrlk(natmtot))
683 if (allocated(chgmt)) deallocate(chgmt)
684 allocate(chgmt(natmtot))
685 if (allocated(mommt)) deallocate(mommt)
686 allocate(mommt(3,natmtot))
687 ! check if scaled spin exchange-correlation should be used
688 tssxc=(abs(sxcscf-1.d0) > 1.d-6)
689 ! spin-spiral phase factors
690 if (ssdph) then
691  if (allocated(zqss)) deallocate(zqss)
692  allocate(zqss(natmtot))
693  do ias=1,natmtot
694  is=idxis(ias)
695  ia=idxia(ias)
696  t1=-0.5d0*dot_product(vqcss(1:3),atposc(1:3,ia,is))
697  zqss(ias)=cmplx(cos(t1),sin(t1),8)
698  end do
699 end if
700 ! mixing vector: either density/magnetisation or potential/magnetic field
701 if (mixrho) then
702  vmixer => rhmg
703 else
704  vmixer => vsbs
705 end if
706 ! zero the mixing vector
707 vmixer(:)=0.d0
708 
709 !-------------------------!
710 ! force variables !
711 !-------------------------!
712 if (tforce) then
713  if (allocated(forcehf)) deallocate(forcehf)
714  allocate(forcehf(3,natmtot))
715  if (allocated(forcetot)) deallocate(forcetot)
716  allocate(forcetot(3,natmtot))
717 end if
718 
719 !-------------------------------------------------!
720 ! DFT+U and fixed tensor moment variables !
721 !-------------------------------------------------!
722 if ((dftu /= 0).or.(ftmtype /= 0)) then
723 ! density matrix elements in each muffin-tin
724  if (allocated(dmatmt)) deallocate(dmatmt)
726 ! potential matrix elements in each muffin-tin
727  if (allocated(vmatmt)) deallocate(vmatmt)
729 ! zero the potential matrix
730  vmatmt(:,:,:,:,:)=0.d0
731 ! matrix elements in spherical coordinates for TDDFT+U
732  if (any(task == [460,461,462,463,478])) then
733  if (allocated(vmatmti)) deallocate(vmatmti)
735  if (allocated(vmatmto)) deallocate(vmatmto)
737  end if
738 ! require the potential matrix elements be calculated
739  tvmatmt=.true.
740 ! flags for non-zero muffin-tin potential matrices
741  if (allocated(tvmmt)) deallocate(tvmmt)
742  allocate(tvmmt(0:lmaxdm,natmtot))
743  tvmmt(:,:)=.false.
744 ! require second-variational eigenvectors
745  tevecsv=.true.
746 end if
747 if (dftu /= 0) then
748  if (any(task == [5,300,600,601,610,620,630,640])) then
749  write(*,*)
750  write(*,'("Error(init0): DFT+U does not work with task ",I8)') task
751  write(*,*)
752  stop
753  end if
754 ! DFT+U energy for each atom
755  if (allocated(engyadu)) deallocate(engyadu)
756  allocate(engyadu(natmmax,ndftu))
757 ! flag the muffin-tin potential matrices which are non-zero
758  do idu=1,ndftu
759  is=isldu(1,idu)
760  if (is > nspecies) then
761  write(*,*)
762  write(*,'("Error(init0): invalid species number : ",I8)') is
763  write(*,*)
764  stop
765  end if
766  l=isldu(2,idu)
767  do ia=1,natoms(is)
768  ias=idxas(ia,is)
769  tvmmt(l,ias)=.true.
770  end do
771  end do
772 ! zero the initial values of screening length
773  lamdu0(:)=0.d0
774 ! write to VARIABLES.OUT
775  if (wrtvars) then
776  call writevars('udufix',nv=ndftu,rva=udufix)
777  end if
778 end if
779 if (ftmtype /= 0) then
780 ! allocate and zero the fixed tensor moment potential array
781  if (allocated(vmftm)) deallocate(vmftm)
783  vmftm(:,:,:,:,:)=0.d0
784 ! flag the muffin-tin potential matrices which are non-zero
785  do i=1,ntmfix
786  is=itmfix(1,i)
787  ia=itmfix(2,i)
788  ias=idxas(ia,is)
789  l=itmfix(3,i)
790  tvmmt(l,ias)=.true.
791  end do
792 end if
793 
794 !-----------------------!
795 ! miscellaneous !
796 !-----------------------!
797 ! determine nuclear radii and volumes
798 call nuclei
799 ! determine the nuclear-nuclear energy
800 call energynn
801 ! get smearing function description
802 call getsdata(stype,sdescr)
803 ! get mixing type description
805 ! generate the spherical harmonic transform (SHT) matrices
806 call genshtmat
807 ! find the maximum size of the spherical Bessel function array over all species
808 call findnjcmax
809 ! allocate 1D plotting arrays
810 if (allocated(dvp1d)) deallocate(dvp1d)
811 allocate(dvp1d(nvp1d))
812 if (allocated(vplp1d)) deallocate(vplp1d)
813 allocate(vplp1d(3,npp1d))
814 if (allocated(dpp1d)) deallocate(dpp1d)
815 allocate(dpp1d(npp1d))
816 ! initial self-consistent loop number
817 iscl=1
818 tlast=.false.
819 ! set the Fermi energy to zero
820 efermi=0.d0
821 ! set the temperature from the smearing width
823 
824 call timesec(ts1)
825 timeinit=timeinit+ts1-ts0
826 
827 end subroutine
828 !EOC
829 
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:907
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:1299
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:935
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:1152
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:923
real(8), dimension(:), allocatable dpp1d
Definition: modmain.f90:1126
integer iscl
Definition: modmain.f90:1051
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:1263
real(8) swidth
Definition: modmain.f90:895
real(8), dimension(:,:), allocatable evalcr
Definition: modmain.f90:937
subroutine gendmftm
Definition: gendmftm.f90:7
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
logical hybrid0
Definition: modmain.f90:1152
real(8) epsocc
Definition: modmain.f90:903
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:941
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:1217
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:993
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:1104
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:1227
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:1223
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:1081
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
complex(4), dimension(:,:,:,:,:), allocatable vmatmto
Definition: moddftu.f90:22
logical tforce
Definition: modmain.f90:989
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:40
real(8), dimension(3, 3) avecref
Definition: modmain.f90:1030
logical, dimension(maxstsp, maxspecies) spcore
Definition: modmain.f90:127
real(8) timeinit
Definition: modmain.f90:1215
integer nrcmtmax
Definition: modmain.f90:175
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
real(8) timepot
Definition: modmain.f90:1225
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:1053
real(8), dimension(:), allocatable dvp1d
Definition: modmain.f90:1122
integer nspncr
Definition: modmain.f90:945
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:901
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:1081
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:991
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:1032
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:1237
real(8), dimension(maxdftu) udufix
Definition: moddftu.f90:60
integer stype
Definition: modmain.f90:891
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:893
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:1219
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:1221
integer, dimension(maxspecies) nstsp
Definition: modmain.f90:113
real(8) gkmax
Definition: modmain.f90:495
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition: modmain.f90:939
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:1234
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:1114
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:1116
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:1154
integer xcgrad
Definition: modmain.f90:602
real(8), dimension(:,:), allocatable vplp1d
Definition: modmain.f90:1124
integer, dimension(2) jspnfv
Definition: modmain.f90:291
logical spincore
Definition: modmain.f90:943
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:1104
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