The Elk Code
 
Loading...
Searching...
No Matches
init4.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine init4
7use modmain
8use modphonon
9use modpw
10use modvars
11implicit none
12! local variables
13integer ik,jspn,n,i
14integer l1,l2,l3,m1,m2,m3
15integer lm1,lm2,lm3
16real(8) vl(3),vc(3)
17! external functions
18real(8), external :: gaunt
19
20!---------------------------!
21! H+k-vector arrays !
22!---------------------------!
23if (any(task == [135,170,171,172,173])) then
24 if (task == 135) hkmax=0.5d0*gmaxvr-epslat
26! allocate the H+k-vector arrays
27 if (allocated(nhk)) deallocate(nhk)
28 allocate(nhk(nspnfv,nkpt))
29 if (allocated(ihkig)) deallocate(ihkig)
30 allocate(ihkig(nhkmax,nspnfv,nkpt))
31 if (allocated(vhkl)) deallocate(vhkl)
32 allocate(vhkl(3,nhkmax,nspnfv,nkpt))
33 if (allocated(vhkc)) deallocate(vhkc)
34 allocate(vhkc(3,nhkmax,nspnfv,nkpt))
35 if (allocated(hkc)) deallocate(hkc)
36 allocate(hkc(nhkmax,nspnfv,nkpt))
37 if (allocated(sfachk)) deallocate(sfachk)
39! initialise H+k-vectors arrays
40 do ik=1,nkpt
41 do jspn=1,nspnfv
42 vl(:)=vkl(:,ik)
43 vc(:)=vkc(:,ik)
44! spin-spiral case
45 if (spinsprl) then
46 if (jspn == 1) then
47 vl(:)=vl(:)+0.5d0*vqlss(:)
48 vc(:)=vc(:)+0.5d0*vqcss(:)
49 else
50 vl(:)=vl(:)-0.5d0*vqlss(:)
51 vc(:)=vc(:)-0.5d0*vqcss(:)
52 end if
53 end if
54! generate the H+k-vectors
55 call gengkvec(ngvec,ivg,vgc,vl,vc,hkmax,nhkmax,nhk(jspn,ik), &
56 ihkig(:,jspn,ik),vhkl(:,:,jspn,ik),vhkc(:,:,jspn,ik),hkc(:,jspn,ik))
57! generate structure factors for H+k-vectors
58 call gensfacgp(nhk(jspn,ik),vhkc(:,:,jspn,ik),nhkmax,sfachk(:,:,jspn,ik))
59 end do
60 end do
61! write to VARIABLES.OUT
62 if (wrtvars) then
63 call writevars('hkmax',rv=hkmax)
64 call writevars('nhk',nv=nspnfv*nkpt,iva=nhk)
65 do ik=1,nkpt
66 do jspn=1,nspnfv
67 call writevars('ihkig',n1=jspn,n2=ik,nv=nhk(jspn,ik), &
68 iva=ihkig(:,jspn,ik))
69 end do
70 end do
71 end if
72end if
73
74!-----------------------------!
75! G+k+q-vector arrays !
76!-----------------------------!
77if (task == 205) then
78 if (allocated(vkql)) deallocate(vkql)
79 allocate(vkql(3,nkptnr))
80 if (allocated(vkqc)) deallocate(vkqc)
81 allocate(vkqc(3,nkptnr))
82 if (allocated(ngkq)) deallocate(ngkq)
83 allocate(ngkq(nspnfv,nkptnr))
84 if (allocated(igkqig)) deallocate(igkqig)
85 allocate(igkqig(ngkmax,nspnfv,nkptnr))
86 if (allocated(vgkql)) deallocate(vgkql)
87 allocate(vgkql(3,ngkmax,nspnfv,nkptnr))
88 if (allocated(vgkqc)) deallocate(vgkqc)
89 allocate(vgkqc(3,ngkmax,nspnfv,nkptnr))
90 if (allocated(gkqc)) deallocate(gkqc)
91 allocate(gkqc(ngkmax,nspnfv,nkptnr))
92 if (allocated(sfacgkq)) deallocate(sfacgkq)
94end if
95
96!---------------------------!
97! G+q-vector arrays !
98!---------------------------!
99if (task == 205) then
100 if (allocated(vgqc)) deallocate(vgqc)
101 allocate(vgqc(3,ngvec))
102 if (allocated(gqc)) deallocate(gqc)
103 allocate(gqc(ngvec))
104 if (allocated(gclgq)) deallocate(gclgq)
105 allocate(gclgq(ngvec))
106 if (allocated(jlgqrmt)) deallocate(jlgqrmt)
107 allocate(jlgqrmt(0:lnpsd,ngvec,nspecies))
108 if (allocated(ylmgq)) deallocate(ylmgq)
109 allocate(ylmgq(lmmaxo,ngvec))
110 if (allocated(sfacgq)) deallocate(sfacgq)
111 allocate(sfacgq(ngvec,natmtot))
112 if (allocated(ffacgq)) deallocate(ffacgq)
113 allocate(ffacgq(ngvec,nspecies))
114 if (allocated(dcfunig)) deallocate(dcfunig)
115 allocate(dcfunig(ngvec))
116 if (allocated(dcfunir)) deallocate(dcfunir)
117 allocate(dcfunir(ngtot))
118end if
119
120!-----------------------------------------------------------------!
121! phonon density functional perturbation theory variables !
122!-----------------------------------------------------------------!
123if (task == 205) then
124 if (allocated(drhomt)) deallocate(drhomt)
125 allocate(drhomt(npmtmax,natmtot+1))
126 if (allocated(drhoir)) deallocate(drhoir)
127 allocate(drhoir(ngtot))
128 if (allocated(dmagmt)) deallocate(dmagmt)
129 if (allocated(dmagir)) deallocate(dmagir)
130 if (spinpol) then
131 allocate(dmagmt(npmtmax,natmtot,ndmag))
132 allocate(dmagir(ngtot,ndmag))
133 end if
134 if (allocated(dvclmt)) deallocate(dvclmt)
135 allocate(dvclmt(npmtmax,natmtot))
136 if (allocated(dvclir)) deallocate(dvclir)
137 allocate(dvclir(ngtot))
138 if (allocated(gvsmt)) deallocate(gvsmt)
139 allocate(gvsmt(npmtmax))
140! combined target array for Kohn-Sham potential and magnetic field derivative
141 if (allocated(dvsbs)) deallocate(dvsbs)
144 allocate(dvsbs(n))
145! zero the array
146 dvsbs(:)=0.d0
147! associate pointer arrays with target
148 dvsmt(1:npmtmax,1:natmtot) => dvsbs(1:)
149 i=npmtmax*natmtot+1
150 dvsig(1:ngvc) => dvsbs(i:)
151 if (spinpol) then
152 i=i+ngvc
153 dbsmt(1:npcmtmax,1:natmtot,1:ndmag) => dvsbs(i:)
155 dbsir(1:ngtot,1:ndmag) => dvsbs(i:)
156 end if
157 if (allocated(dvsir)) deallocate(dvsir)
158 allocate(dvsir(ngtot))
159 if (allocated(dsocfr)) deallocate(dsocfr)
160 if (spinorb) then
161 allocate(dsocfr(nrcmtmax,natmtot))
162 end if
163 if (allocated(dhaa)) deallocate(dhaa)
165 if (allocated(dhloa)) deallocate(dhloa)
167 if (allocated(dhlolo)) deallocate(dhlolo)
169! allocate and generate real Gaunt coefficient array
170 if (allocated(gntyyy)) deallocate(gntyyy)
171 allocate(gntyyy(lmmaxo,lmmaxapw,lmmaxapw))
172 do l1=0,lmaxapw
173 do m1=-l1,l1
174 lm1=l1*(l1+1)+m1+1
175 do l3=0,lmaxapw
176 do m3=-l3,l3
177 lm3=l3*(l3+1)+m3+1
178 do l2=0,lmaxo
179 do m2=-l2,l2
180 lm2=l2*(l2+1)+m2+1
181 gntyyy(lm2,lm3,lm1)=gaunt(l1,l2,l3,m1,m2,m3)
182 end do
183 end do
184 end do
185 end do
186 end do
187 end do
188 if (allocated(devalfv)) deallocate(devalfv)
189 allocate(devalfv(nstfv,nspnfv,nkptnr))
190 if (allocated(devalsv)) deallocate(devalsv)
191 allocate(devalsv(nstsv,nkptnr))
192 if (allocated(doccsv)) deallocate(doccsv)
193 allocate(doccsv(nstsv,nkptnr))
194end if
195
196end subroutine
197
pure subroutine findngkmax(nkpt, vkc, nspnfv, vqcss, ngv, vgc, gkmax, ngkmax)
real(8) function gaunt(l1, l2, l3, m1, m2, m3)
Definition gaunt.f90:10
pure subroutine gengkvec(ngv, ivg, vgc, vkl, vkc, gkmax, ngkmax, ngk, igkig, vgkl, vgkc, gkc)
Definition gengkvec.f90:11
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition gensfacgp.f90:10
subroutine init4
Definition init4.f90:7
real(8) gmaxvr
Definition modmain.f90:384
integer ngtot
Definition modmain.f90:390
integer nkptnr
Definition modmain.f90:463
logical spinpol
Definition modmain.f90:228
integer ngvec
Definition modmain.f90:396
logical spinorb
Definition modmain.f90:230
integer lmmaxapw
Definition modmain.f90:199
real(8), dimension(3) vqcss
Definition modmain.f90:295
integer apwordmax
Definition modmain.f90:760
integer nlomax
Definition modmain.f90:788
real(8) epslat
Definition modmain.f90:24
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer nspnfv
Definition modmain.f90:289
integer lmaxapw
Definition modmain.f90:197
integer npcmtmax
Definition modmain.f90:216
integer task
Definition modmain.f90:1298
integer nstsv
Definition modmain.f90:886
integer npmtmax
Definition modmain.f90:216
logical spinsprl
Definition modmain.f90:283
integer lmaxo
Definition modmain.f90:201
integer, dimension(:,:), allocatable ivg
Definition modmain.f90:400
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
real(8), dimension(3) vqlss
Definition modmain.f90:293
integer nrcmtmax
Definition modmain.f90:175
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
integer lmmaxo
Definition modmain.f90:203
integer ngvc
Definition modmain.f90:398
integer nstfv
Definition modmain.f90:884
integer nspecies
Definition modmain.f90:34
integer lnpsd
Definition modmain.f90:628
integer ndmag
Definition modmain.f90:238
real(8), dimension(:,:), allocatable vkql
Definition modphonon.f90:54
integer, dimension(:,:), allocatable ngkq
Definition modphonon.f90:84
complex(8), dimension(:,:,:,:), allocatable sfacgkq
Definition modphonon.f90:92
real(8), dimension(:,:), allocatable doccsv
real(8), dimension(:,:,:,:), allocatable vgkqc
Definition modphonon.f90:88
complex(8), dimension(:,:), allocatable ylmgq
Definition modphonon.f90:70
complex(8), dimension(:), allocatable dcfunir
Definition modphonon.f90:78
real(8), dimension(:), allocatable gclgq
Definition modphonon.f90:66
real(8), dimension(:,:,:), allocatable jlgqrmt
Definition modphonon.f90:68
real(8), dimension(:), allocatable gqc
Definition modphonon.f90:64
complex(8), dimension(:), allocatable dvsir
complex(8), dimension(:,:,:,:), allocatable dhlolo
complex(8), dimension(:,:), pointer, contiguous dvsmt
complex(8), dimension(:,:,:), allocatable dmagmt
complex(8), dimension(:,:,:,:,:), allocatable dhloa
complex(8), dimension(:,:), pointer, contiguous dbsir
complex(8), dimension(:,:), allocatable dsocfr
real(8), dimension(:,:), allocatable ffacgq
Definition modphonon.f90:74
complex(8), dimension(:), allocatable dcfunig
Definition modphonon.f90:76
complex(8), dimension(:,:,:,:,:,:), allocatable dhaa
real(8), dimension(:,:), allocatable devalsv
complex(8), dimension(:,:), allocatable dmagir
complex(8), dimension(:), allocatable drhoir
Definition modphonon.f90:98
real(8), dimension(:,:,:,:), allocatable vgkql
Definition modphonon.f90:88
real(8), dimension(:,:,:), allocatable gntyyy
complex(8), dimension(:,:), allocatable drhomt
Definition modphonon.f90:98
real(8), dimension(:,:,:), allocatable devalfv
real(8), dimension(:,:,:), allocatable gkqc
Definition modphonon.f90:90
real(8), dimension(:,:), allocatable vgqc
Definition modphonon.f90:62
real(8), dimension(:,:), allocatable vkqc
Definition modphonon.f90:56
integer, dimension(:,:,:), allocatable igkqig
Definition modphonon.f90:86
complex(8), dimension(:,:), allocatable sfacgq
Definition modphonon.f90:72
complex(8), dimension(:), allocatable, target dvsbs
complex(8), dimension(:), allocatable gvsmt
complex(8), dimension(:,:), allocatable dvclmt
complex(8), dimension(:), allocatable dvclir
complex(8), dimension(:,:,:), pointer, contiguous dbsmt
complex(8), dimension(:), pointer, contiguous dvsig
Definition modpw.f90:6
integer, dimension(:,:,:), allocatable ihkig
Definition modpw.f90:44
real(8), dimension(:,:,:), allocatable hkc
Definition modpw.f90:50
real(8), dimension(:,:,:,:), allocatable vhkc
Definition modpw.f90:48
real(8) hkmax
Definition modpw.f90:38
integer nhkmax
Definition modpw.f90:42
integer, dimension(:,:), allocatable nhk
Definition modpw.f90:40
complex(8), dimension(:,:,:,:), allocatable sfachk
Definition modpw.f90:52
real(8), dimension(:,:,:,:), allocatable vhkl
Definition modpw.f90:46
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