The Elk Code
 
Loading...
Searching...
No Matches
deveqnfv.f90
Go to the documentation of this file.
1
2! Copyright (C) 2013 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 deveqnfv(ngp,ngpq,igpig,igpqig,vgpc,vgpqc,evalfv,apwalm,apwalmq, &
7 dapwalm,dapwalmq,evecfv,devalfvp,devecfv)
8use modmain
9use modphonon
10use modomp
11implicit none
12! arguments
13integer, intent(in) :: ngp,ngpq
14integer, intent(in) :: igpig(ngkmax),igpqig(ngkmax)
15real(8), intent(in) :: vgpc(3,ngkmax),vgpqc(3,ngkmax)
16real(8), intent(in) :: evalfv(nstfv)
17complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
18complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw,natmtot)
19complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
20complex(8), intent(in) :: dapwalmq(ngkmax,apwordmax,lmmaxapw)
21complex(8), intent(in) :: evecfv(nmatmax,nstfv)
22real(8), intent(out) :: devalfvp(nstfv)
23complex(8), intent(out) :: devecfv(nmatmax,nstfv)
24! local variables
25integer nm,nmq,is,ias,jst,i
26integer lwork,info,nthd
27real(8) t1
28complex(8) z1
29! allocatable arrays
30real(8), allocatable :: w(:),rwork(:)
31complex(8), allocatable :: h(:,:),o(:,:),dh(:,:),od(:,:)
32complex(8), allocatable :: x(:),y(:),work(:)
33! external functions
34real(8), external :: ddot
35! matrix sizes for k and k+q
36nm=ngp+nlotot
37nmq=ngpq+nlotot
38allocate(h(nmq,nmq),o(nmq,nmq))
39! compute the Hamiltonian and overlap matrices at p+q
40call holdthd(2,nthd)
41!$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
42!$OMP NUM_THREADS(nthd)
43!$OMP SECTION
44call hmlfv(nmq,ngpq,igpqig,vgpqc,apwalmq,h)
45!$OMP SECTION
46call olpfv(nmq,ngpq,igpqig,apwalmq,o)
47!$OMP END PARALLEL SECTIONS
48call freethd(nthd)
49! solve the generalised eigenvalue problem (H - eⱼ O)|vⱼ> = 0
50! (note: these are also the eigenvalues/vectors of O⁻¹ H )
51lwork=2*nmq
52allocate(w(nmq),rwork(3*nmq),work(lwork))
53call zhegv(1,'V','U',nmq,h,nmq,o,nmq,w,work,lwork,rwork,info)
54if (info /= 0) then
55 write(*,*)
56 write(*,'("Error(deveqnfv): diagonalisation failed")')
57 write(*,'(" ZHEGV returned INFO = ",I8)') info
58 write(*,*)
59 stop
60end if
61deallocate(rwork,o,work)
62! compute the Hamiltonian and overlap matrix derivatives
63allocate(dh(nmq,nm),od(nmq,nm))
64call holdthd(2,nthd)
65!$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
66!$OMP PRIVATE(ias,is) &
67!$OMP NUM_THREADS(nthd)
68!$OMP SECTION
69call dhmlistl(ngp,ngpq,igpig,igpqig,vgpc,vgpqc,nmq,dh)
70dh(ngpq+1:nmq,1:ngp)=0.d0
71dh(1:nmq,ngp+1:nm)=0.d0
72do ias=1,natmtot
73 is=idxis(ias)
74 call dhmlaa(is,ias,ngp,ngpq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),dapwalm, &
75 dapwalmq,nmq,dh)
76 call dhmlalo(is,ias,ngp,ngpq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),dapwalm, &
77 dapwalmq,nmq,dh)
78 call dhmllolo(is,ias,ngp,ngpq,nmq,dh)
79end do
80!$OMP SECTION
81call dolpistl(ngp,ngpq,igpig,igpqig,nmq,od)
82od(ngpq+1:nmq,1:ngp)=0.d0
83od(1:nmq,ngp+1:nm)=0.d0
84do ias=1,natmtot
85 is=idxis(ias)
86 call dolpaa(is,ias,ngp,ngpq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),dapwalm, &
87 dapwalmq,nmq,od)
88 call dolpalo(is,ias,ngp,ngpq,dapwalm,dapwalmq,nmq,od)
89end do
90!$OMP END PARALLEL SECTIONS
91call freethd(nthd)
92allocate(x(nmq),y(nmq))
93! loop over states
94do jst=1,nstfv
95! compute |dvⱼ> = V (eⱼ - D)⁻¹ V† (dH - eⱼ dO)|vⱼ>
96 z1=-evalfv(jst)
97 call zgemv('N',nmq,nm,z1,od,nmq,evecfv(:,jst),1,zzero,x,1)
98 call zgemv('N',nmq,nm,zone,dh,nmq,evecfv(:,jst),1,zone,x,1)
99! compute the first-order change in eigenvalue
100 if (tphq0) then
101 devalfvp(jst)=ddot(2*nmq,evecfv(:,jst),1,x,1)
102 else
103 devalfvp(jst)=0.d0
104 end if
105 call zgemv('C',nmq,nmq,zone,h,nmq,x,1,zzero,y,1)
106 do i=1,nmq
107 t1=evalfv(jst)-w(i)
108 if (abs(t1) > epsdev) then
109 y(i)=y(i)/t1
110 else
111 y(i)=0.d0
112 end if
113 end do
114 call zgemv('N',nmq,nmq,zone,h,nmq,y,1,zzero,devecfv(:,jst),1)
115end do
116deallocate(w,h,dh,od,x,y)
117end subroutine
118
subroutine deveqnfv(ngp, ngpq, igpig, igpqig, vgpc, vgpqc, evalfv, apwalm, apwalmq, dapwalm, dapwalmq, evecfv, devalfvp, devecfv)
Definition deveqnfv.f90:8
subroutine dhmlaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, dh)
Definition dhmlaa.f90:7
pure subroutine dhmlalo(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, dh)
Definition dhmlalo.f90:7
pure subroutine dhmlistl(ngp, ngpq, igpig, igpqig, vgpc, vgpqc, ld, dh)
Definition dhmlistl.f90:7
pure subroutine dhmllolo(is, ias, ngp, ngpq, ld, dh)
Definition dhmllolo.f90:7
subroutine dolpaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, od)
Definition dolpaa.f90:7
pure subroutine dolpalo(is, ias, ngp, ngpq, dapwalm, dapwalmq, ld, od)
Definition dolpalo.f90:7
pure subroutine dolpistl(ngp, ngpq, igpig, igpqig, ld, od)
Definition dolpistl.f90:7
subroutine hmlfv(nmatp, ngp, igpig, vgpc, apwalm, h)
Definition hmlfv.f90:7
complex(8), parameter zzero
Definition modmain.f90:1238
complex(8), parameter zone
Definition modmain.f90:1238
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer nlotot
Definition modmain.f90:790
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
logical tphq0
Definition modphonon.f90:17
real(8) epsdev
subroutine olpfv(nmatp, ngp, igpig, apwalm, o)
Definition olpfv.f90:7