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