The Elk Code
phonon.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 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 
6 subroutine phonon
7 use modmain
8 use modphonon
9 use modpw
10 use modmpi
11 use modomp
12 use moddelf
13 implicit none
14 ! local variables
15 integer ik,jk,jspn,idm
16 integer is,ia,ias,ip,nthd
17 integer nmix,nwork,n,lp
18 ! use Broyden mixing only
19 integer, parameter :: mtype=3
20 real(8) ddv,a,b
21 character(256) fext
22 ! allocatable arrays
23 real(8), allocatable :: evalfv(:,:),work(:)
24 complex(8), allocatable :: dyn(:,:)
25 complex(8), allocatable :: apwalm(:,:,:,:,:),apwalmq(:,:,:,:,:)
26 complex(8), allocatable :: dapwalm(:,:,:,:),dapwalmq(:,:,:,:)
27 complex(8), allocatable :: evecfv(:,:,:),devecfv(:,:,:)
28 complex(8), allocatable :: evecsv(:,:),devecsv(:,:)
29 ! initialise universal variables
30 call init0
31 call init1
32 call init2
33 call init4
34 if (lmaxi < 2) then
35  write(*,*)
36  write(*,'("Error(phonon): lmaxi too small for calculating DFPT phonons : ",&
37  &I0)') lmaxi
38  write(*,'(" Run the ground-state calculation again with lmaxi >= 2")')
39  write(*,*)
40  stop
41 end if
42 if (spinpol) then
43  write(*,*)
44  write(*,'("Error(phonon): spin-polarised phonons not yet available")')
45  write(*,*)
46  stop
47 end if
48 ! check spin-spiral de-phasing is not used
49 if (ssdph) then
50  write(*,*)
51  write(*,'("Error(phonon): ssdph should be .false. for DFPT phonons")')
52  write(*,*)
53  stop
54 end if
55 ! check for zero atoms
56 if (natmtot == 0) return
57 ! read in the density and potentials
58 call readstate
59 ! Fourier transform Kohn-Sham potential to G-space
60 call genvsig
61 ! find the new linearisation energies
62 call linengy
63 ! generate the APW and local-orbital radial functions and integrals
64 call genapwlofr
65 ! generate the spin-orbit coupling radial functions
66 call gensocfr
67 ! get the eigenvalues and occupation numbers from file
68 do ik=1,nkpt
69  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
70  call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
71 end do
72 ! size of mixing vector (complex array)
73 nmix=2*size(dvsbs)
74 ! determine the size of the mixer work array
75 nwork=-1
76 call mixerifc(mtype,nmix,dvsbs,ddv,nwork,work)
77 allocate(work(nwork))
78 ! allocate dynamical matrix column
79 allocate(dyn(3,natmtot))
80 ! begin new dynamical matrix task
81 10 continue
82 call dyntask(80,fext)
83 ! if nothing more to do then return
84 if (iqph == 0) return
85 if (mp_mpi) then
86  write(*,'("Info(phonon): working on ",A)') 'DYN'//trim(fext)
87 ! open RMSDDVS.OUT
88  open(65,file='RMSDDVS'//trim(fext),form='FORMATTED')
89 end if
90 ! zero the dynamical matrix row
91 dyn(:,:)=0.d0
92 ! check to see if mass is considered infinite
93 if (spmass(isph) <= 0.d0) goto 20
94 ! generate the G+k+q-vectors and related quantities
95 call gengkqvec(iqph)
96 ! generate the G+q-vectors and related quantities
97 call gengqvec(iqph)
98 ! generate the Coulomb Green's function in G+q-space
99 call gengclgq(.false.,iqph,ngvec,gqc,gclgq)
100 ! generate the characteristic function derivative
101 call gendcfun
102 ! generate the gradient of the Kohn-Sham potential
103 call gengvsmt
104 ! zero the density derivative
105 drhmg(:)=0.d0
106 ! initialise the potential derivative
107 call dpotks
108 call gendvsig
109 ! initialise the mixer
110 iscl=0
111 call mixerifc(mtype,nmix,dvsbs,ddv,nwork,work)
112 ! initialise the Fermi energy and occupancy derivatives
113 defermi=0.d0
114 doccsv(:,:)=0.d0
115 ! set last self-consistent loop flag
116 tlast=.false.
117 ! begin the self-consistent loop
118 do iscl=1,maxscl
119 ! compute the Hamiltonian radial integral derivatives
120  call dhmlrad
121 ! zero the density and magnetisation derivatives
122  drhmg(:)=0.d0
123 ! parallel loop over k-points
124  call holdthd(nkptnr/np_mpi,nthd)
125 !$OMP PARALLEL DEFAULT(SHARED) &
126 !$OMP PRIVATE(evalfv,apwalm,apwalmq,dapwalm,dapwalmq) &
127 !$OMP PRIVATE(evecfv,devecfv,evecsv,devecsv,jk,jspn) &
128 !$OMP NUM_THREADS(nthd)
129  allocate(evalfv(nstfv,nspnfv))
130  allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
131  allocate(apwalmq(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
132  allocate(dapwalm(ngkmax,apwordmax,lmmaxapw,nspnfv))
133  allocate(dapwalmq(ngkmax,apwordmax,lmmaxapw,nspnfv))
134  allocate(evecfv(nmatmax,nstfv,nspnfv),devecfv(nmatmax,nstfv,nspnfv))
135  allocate(evecsv(nstsv,nstsv),devecsv(nstsv,nstsv))
136 !$OMP DO SCHEDULE(DYNAMIC)
137  do ik=1,nkptnr
138 ! distribute among MPI processes
139  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
140 ! equivalent reduced k-point
141  jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
142 ! compute the matching coefficients and derivatives
143  do jspn=1,nspnfv
144  call match(ngk(jspn,ik),vgkc(:,:,jspn,ik),gkc(:,jspn,ik), &
145  sfacgk(:,:,jspn,ik),apwalm(:,:,:,:,jspn))
146  call dmatch(iasph,ipph,ngk(jspn,ik),vgkc(:,:,jspn,ik), &
147  apwalm(:,:,:,:,jspn),dapwalm(:,:,:,jspn))
148  call match(ngkq(jspn,ik),vgkqc(:,:,jspn,ik),gkqc(:,jspn,ik), &
149  sfacgkq(:,:,jspn,ik),apwalmq(:,:,:,:,jspn))
150  call dmatch(iasph,ipph,ngkq(jspn,ik),vgkqc(:,:,jspn,ik), &
151  apwalmq(:,:,:,:,jspn),dapwalmq(:,:,:,jspn))
152  end do
153 ! get the first- and second-variational eigenvalues and eigenvectors from file
154  call getevalfv(filext,jk,vkl(:,ik),evalfv)
155  call getevecfv(filext,0,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
156  call getevecsv(filext,0,vkl(:,ik),evecsv)
157 ! solve the first-variational eigenvalue equation derivative
158  do jspn=1,nspnfv
159  call deveqnfv(ngk(jspn,ik),ngkq(jspn,ik),igkig(:,jspn,ik), &
160  igkqig(:,jspn,ik),vgkc(:,:,jspn,ik),vgkqc(:,:,jspn,ik),evalfv(:,jspn), &
161  apwalm(:,:,:,:,jspn),apwalmq(:,:,:,:,jspn),dapwalm(:,:,:,jspn), &
162  dapwalmq(:,:,:,jspn),evecfv(:,:,jspn),devalfv(:,jspn,ik), &
163  devecfv(:,:,jspn))
164  end do
165  if (spinsprl) then
166 ! solve the spin-spiral second-variational eigenvalue equation derivative
167 ! call deveqnss(ngk(:,ik),ngkq(:,ik),igkig(:,:,ik),igkqig(:,:,ik),apwalm, &
168 ! dapwalm,devalfv,evecfv,evecfvq,devecfv,evalsv(:,jk),evecsv,devecsv)
169  else
170 ! solve the second-variational eigenvalue equation derivative
171 ! call deveqnsv(ngk(1,ik),ngkq(1,ik),igkig(:,1,ik),igkqig(:,1,ik), &
172 ! vgkqc(:,:,1,ik),apwalm,dapwalm,devalfv,evecfv,evecfvq,devecfv, &
173 ! evalsv(:,jk),evecsv,devecsv)
174  end if
175 
176 !*******
177 devalsv(:,ik)=devalfv(:,1,ik)
178 !*******
179 
180 ! write the eigenvalue/vector derivatives to file
181  call putdevecfv(ik,devecfv)
182  if (tevecsv) call putdevecsv(ik,devecsv)
183 ! add to the density and magnetisation derivatives
184  call drhomagk(ngk(:,ik),ngkq(:,ik),igkig(:,:,ik),igkqig(:,:,ik), &
185  occsv(:,jk),doccsv(:,ik),apwalm,apwalmq,dapwalm,evecfv,devecfv,evecsv, &
186  devecsv)
187  end do
188 !$OMP END DO
189  deallocate(evalfv,apwalm,apwalmq,dapwalm,dapwalmq)
190  deallocate(evecfv,devecfv,evecsv,devecsv)
191 !$OMP END PARALLEL
192  call freethd(nthd)
193 ! broadcast eigenvalue derivative arrays to every process
194  n=nstfv*nspnfv
195  do ik=1,nkptnr
196  lp=mod(ik-1,np_mpi)
197  call mpi_bcast(devalfv(:,:,ik),n,mpi_double_precision,lp,mpicom,ierror)
198  call mpi_bcast(devalsv(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
199  end do
200 ! convert muffin-tin density/magnetisation derivatives to spherical harmonics
201  call drhomagsh
202 ! convert from a coarse to a fine radial mesh
203  call zfmtctof(drhomt)
204  do idm=1,ndmag
205  call zfmtctof(dmagmt(:,:,idm))
206  end do
207 ! convert interstitial density/magnetisation derivatives to fine grid
208  call zfirctof(drhoir,drhoir)
209  do idm=1,ndmag
210  call zfirctof(dmagir(:,idm),dmagir(:,idm))
211  end do
212 ! add densities from each process and redistribute
213  if (np_mpi > 1) then
214  n=size(drhmg)
215  call mpi_allreduce(mpi_in_place,drhmg,n,mpi_double_complex,mpi_sum,mpicom, &
216  ierror)
217  end if
218 ! synchronise MPI processes
219  call mpi_barrier(mpicom,ierror)
220 ! add gradient contribution to density derivative
221  call gradrhomt
222 ! compute the Kohn-Sham potential derivative
223  call dpotks
224 ! Fourier transform Kohn-Sham potential derivative to G+q-space
225  call gendvsig
226 ! mix the old potential and field with the new
227  call mixerifc(mtype,nmix,dvsbs,ddv,nwork,work)
228  if (mp_mpi) then
229  write(65,'(G18.10)') ddv
230  flush(65)
231  end if
232 ! exit self-consistent loop if required
233  if (tlast) goto 20
234 ! check for convergence
235  if (iscl >= 2) then
236  if (ddv < epspot) tlast=.true.
237  end if
238 ! broadcast tlast from master process to all other processes
239  call mpi_bcast(tlast,1,mpi_logical,0,mpicom,ierror)
240 ! compute the occupation number derivatives
241  call doccupy
242 ! end the self-consistent loop
243 end do
244 write(*,*)
245 write(*,'("Warning(phonon): failed to reach self-consistency after ",I0,&
246  &" loops")') maxscl
247 20 continue
248 ! close the RMSDDVS.OUT file
249 if (mp_mpi) close(65)
250 ! synchronise MPI processes
251 call mpi_barrier(mpicom,ierror)
252 ! generate the dynamical matrix row from force derivatives
253 call dforce(dyn)
254 ! synchronise MPI processes
255 call mpi_barrier(mpicom,ierror)
256 ! write dynamical matrix row to file
257 if (mp_mpi) then
258  do ias=1,natmtot
259  is=idxis(ias)
260  ia=idxia(ias)
261  do ip=1,3
262  a=dble(dyn(ip,ias))
263  b=aimag(dyn(ip,ias))
264  if (abs(a) < 1.d-12) a=0.d0
265  if (abs(b) < 1.d-12) b=0.d0
266  write(80,'(2G18.10," : is = ",I4,", ia = ",I4,", ip = ",I4)') a,b,is,ia,ip
267  end do
268  end do
269  close(80)
270 ! write the complex Kohn-Sham potential derivative to file
271  call writedvs(fext)
272 end if
273 ! delete the non-essential files
274 call delfiles(devec=.true.)
275 ! synchronise MPI processes
276 call mpi_barrier(mpicom,ierror)
277 goto 10
278 end subroutine
279 
integer nmatmax
Definition: modmain.f90:851
subroutine doccupy
Definition: doccupy.f90:7
subroutine zfmtctof(zfmt)
Definition: zfmtctof.f90:7
subroutine dpotks
Definition: dpotks.f90:7
subroutine drhomagk(ngp, ngpq, igpig, igpqig, occsvp, doccsvp, apwalm, apwalmq, dapwalm, evecfv, devecfv, evecsv, devecsv)
Definition: drhomagk.f90:8
complex(8), dimension(:,:), pointer, contiguous dmagir
Definition: modphonon.f90:102
subroutine dforce(dyn)
Definition: dforce.f90:7
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1295
subroutine dhmlrad
Definition: dhmlrad.f90:7
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:915
logical mp_mpi
Definition: modmpi.f90:17
logical spinpol
Definition: modmain.f90:228
complex(8), dimension(:,:), pointer, contiguous drhomt
Definition: modphonon.f90:100
integer lmmaxapw
Definition: modmain.f90:199
integer nkpt
Definition: modmain.f90:461
real(8), dimension(:,:,:), allocatable gkqc
Definition: modphonon.f90:90
integer ndmag
Definition: modmain.f90:238
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
logical tevecsv
Definition: modmain.f90:917
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
subroutine gendvsig
Definition: gendvsig.f90:7
integer iscl
Definition: modmain.f90:1045
logical spinsprl
Definition: modmain.f90:283
Definition: modomp.f90:6
integer iasph
Definition: modphonon.f90:15
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
subroutine dyntask(fnum, fext)
Definition: dyntask.f90:7
real(8), dimension(:,:), allocatable devalsv
Definition: modphonon.f90:134
subroutine genvsig
Definition: genvsig.f90:10
integer nkptnr
Definition: modmain.f90:463
subroutine phonon
Definition: phonon.f90:7
subroutine genapwlofr
Definition: genapwlofr.f90:7
subroutine gensocfr
Definition: gensocfr.f90:7
integer np_mpi
Definition: modmpi.f90:13
subroutine gengkqvec(iq)
Definition: gengkqvec.f90:7
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine linengy
Definition: linengy.f90:10
subroutine mixerifc(mtype, n, v, dv, nwork, work)
Definition: mixerifc.f90:7
subroutine putdevecfv(ik, devecfv)
Definition: putdevecfv.f90:7
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition: gengclgq.f90:7
integer nstsv
Definition: modmain.f90:885
complex(8), dimension(:,:,:), pointer, contiguous dmagmt
Definition: modphonon.f90:102
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
logical tlast
Definition: modmain.f90:1047
complex(8), dimension(:,:,:,:), allocatable sfacgkq
Definition: modphonon.f90:92
subroutine putdevecsv(ik, devecsv)
Definition: putdevecsv.f90:7
real(8), dimension(:,:), allocatable doccsv
Definition: modphonon.f90:136
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
complex(8), dimension(:), allocatable, target drhmg
Definition: modphonon.f90:98
real(8), dimension(:,:,:), allocatable devalfv
Definition: modphonon.f90:132
integer ipph
Definition: modphonon.f90:15
real(8), dimension(maxspecies) spmass
Definition: modmain.f90:101
subroutine init2
Definition: init2.f90:7
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:901
subroutine getevalfv(fext, ikp, vpl, evalfv)
Definition: getevalfv.f90:7
integer isph
Definition: modphonon.f90:15
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
subroutine init1
Definition: init1.f90:10
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition: moddelf.f90:25
subroutine drhomagsh
Definition: drhomagsh.f90:7
logical ssdph
Definition: modmain.f90:285
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
Definition: modmpi.f90:6
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition: getoccsv.f90:7
real(8) epspot
Definition: modmain.f90:1055
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
subroutine init4
Definition: init4.f90:7
subroutine deveqnfv(ngp, ngpq, igpig, igpqig, vgpc, vgpqc, evalfv, apwalm, apwalmq, dapwalm, dapwalmq, evecfv, devalfvp, devecfv)
Definition: deveqnfv.f90:8
complex(8), dimension(:), pointer, contiguous drhoir
Definition: modphonon.f90:100
subroutine zfirctof(zfirc, zfir)
Definition: zfirctof.f90:7
Definition: modpw.f90:6
complex(8), dimension(:), allocatable, target dvsbs
Definition: modphonon.f90:108
subroutine readstate
Definition: readstate.f90:10
integer iqph
Definition: modphonon.f90:15
subroutine gradrhomt
Definition: gradrhomt.f90:7
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(:,:), allocatable ngkq
Definition: modphonon.f90:84
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
subroutine gengvsmt
Definition: gengvsmt.f90:7
real(8), dimension(:,:,:,:), allocatable vgkqc
Definition: modphonon.f90:88
subroutine gendcfun
Definition: gendcfun.f90:7
subroutine init0
Definition: init0.f90:10
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:), allocatable gqc
Definition: modphonon.f90:64
pure subroutine dmatch(ias, ip, ngp, vgpc, apwalm, dapwalm)
Definition: dmatch.f90:7
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
real(8) defermi
Definition: modphonon.f90:130
integer nstfv
Definition: modmain.f90:883
subroutine writedvs(fext)
Definition: writedvs.f90:7
integer maxscl
Definition: modmain.f90:1043
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465
integer mpicom
Definition: modmpi.f90:11
integer nspnfv
Definition: modmain.f90:289
subroutine gengqvec(iq)
Definition: gengqvec.f90:7
integer, dimension(:,:,:), allocatable igkqig
Definition: modphonon.f90:86
real(8), dimension(:), allocatable gclgq
Definition: modphonon.f90:66
integer ierror
Definition: modmpi.f90:19
integer lmaxi
Definition: modmain.f90:205