The Elk Code
eveqnit.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 subroutine eveqnit(nmatp,ngp,igpig,vpl,vgpl,vgpc,apwalm,evalfv,evecfv)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: nmatp,ngp,igpig(ngkmax)
12 real(8), intent(in) :: vpl(3),vgpl(3,ngkmax),vgpc(3,ngkmax)
13 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
14 real(8), intent(out) :: evalfv(nstfv)
15 complex(8), intent(out) :: evecfv(nmatmax,nstfv)
16 ! local variables
17 integer ns,ist,it,i,j
18 integer is,ias,nthd
19 integer nb,lwork,info
20 real(8) t1
21 real(8) ts1,ts0
22 ! allocatable arrays
23 integer, allocatable :: idx(:)
24 real(8), allocatable :: w(:),rwork(:)
25 complex(8), allocatable :: h(:,:),o(:,:),hv(:,:),ov(:,:)
26 complex(8), allocatable :: u(:,:),hu(:,:),ou(:,:)
27 complex(8), allocatable :: hs(:,:),os(:,:),work(:)
28 ! external functions
29 integer, external :: ilaenv
30 real(8), external :: ddot
31 ! compute Hamiltonian and overlap matrices
32 call timesec(ts0)
33 allocate(h(nmatp,nmatp),o(nmatp,nmatp))
34 call holdthd(2,nthd)
35 !$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
36 !$OMP PRIVATE(j,ias,is) &
37 !$OMP NUM_THREADS(nthd)
38 !$OMP SECTION
39 ! Hamiltonian
40 call hmlistl(ngp,igpig,vgpc,nmatp,h)
41 do j=ngp+1,nmatp
42  h(1:j,j)=0.d0
43 end do
44 do ias=1,natmtot
45  is=idxis(ias)
46  call hmlaa(tefvr,is,ias,ngp,apwalm(:,:,:,ias),nmatp,h)
47  call hmlalo(is,ias,ngp,apwalm(:,:,:,ias),nmatp,h)
48  call hmllolo(is,ias,ngp,nmatp,h)
49 end do
50 !$OMP SECTION
51 ! overlap
52 call olpistl(ngp,igpig,nmatp,o)
53 do j=ngp+1,nmatp
54  o(1:j,j)=0.d0
55 end do
56 do ias=1,natmtot
57  is=idxis(ias)
58  call olpaa(tefvr,is,ngp,apwalm(:,:,:,ias),nmatp,o)
59  call olpalo(is,ias,ngp,apwalm(:,:,:,ias),nmatp,o)
60  call olplolo(is,ias,ngp,nmatp,o)
61 end do
62 !$OMP END PARALLEL SECTIONS
63 call freethd(nthd)
64 call timesec(ts1)
65 !$OMP ATOMIC
66 timemat=timemat+ts1-ts0
67 if (trdstate.or.(iscl >= 2)) then
68 ! read in the eigenvectors from file if available
69  call getevecfv(filext,0,vpl,vgpl,evecfv)
70 else
71 ! initialise the eigenvalues/vectors using the diagonal elements of H
72  allocate(idx(nmatp),w(nmatp))
73  do i=1,nmatp
74  w(i)=h(i,i)
75  end do
76  call sortidx(nmatp,w,idx)
77  evecfv(1:nmatp,1:nstfv)=0.d0
78  do ist=1,nstfv
79  evecfv(idx(ist),ist)=1.d0
80  end do
81  deallocate(idx,w)
82 end if
83 call timesec(ts0)
84 ns=2*nstfv
85 allocate(w(ns),rwork(3*ns))
86 allocate(hv(nmatp,nstfv),ov(nmatp,nstfv))
87 allocate(u(nmatp,nstfv),hu(nmatp,nstfv),ou(nmatp,nstfv))
88 allocate(hs(ns,ns),os(ns,ns))
89 ! find the optimal blocksize for allocating the work array
90 nb=ilaenv(1,'ZHETRD','U',nmatp,-1,-1,-1)
91 nb=max(nb,1)
92 lwork=(nb+1)*ns
93 allocate(work(lwork))
94 call holdthd(nstfv,nthd)
95 ! iteration loop
96 !$OMP PARALLEL DEFAULT(SHARED) &
97 !$OMP PRIVATE(it,ist,t1) &
98 !$OMP NUM_THREADS(nthd)
99 do it=1,nefvit
100 !$OMP DO SCHEDULE(DYNAMIC)
101  do ist=1,nstfv
102 ! operate with O on the current eigenvector
103  call zhemv('U',nmatp,zone,o,nmatp,evecfv(:,ist),1,zzero,ov(:,ist),1)
104 ! normalise the eigenvector
105  t1=ddot(2*nmatp,evecfv(:,ist),1,ov(:,ist),1)
106  if (t1 > 0.d0) then
107  t1=1.d0/sqrt(t1)
108  call zdscal(nmatp,t1,evecfv(:,ist),1)
109  call zdscal(nmatp,t1,ov(:,ist),1)
110  end if
111 ! operate with H on the current eigenvector
112  call zhemv('U',nmatp,zone,h,nmatp,evecfv(:,ist),1,zzero,hv(:,ist),1)
113 ! estimate the eigenvalue
114  evalfv(ist)=ddot(2*nmatp,evecfv(:,ist),1,hv(:,ist),1)
115 ! compute the residual |u⟩ = (H - eO)|v⟩
116  call zcopy(nmatp,hv(:,ist),1,u(:,ist),1)
117  t1=-evalfv(ist)
118  u(1:nmatp,ist)=u(1:nmatp,ist)+t1*ov(1:nmatp,ist)
119 ! apply the overlap matrix to the residual
120  call zhemv('U',nmatp,zone,o,nmatp,u(:,ist),1,zzero,ou(:,ist),1)
121 ! apply the Hamiltonian matrix to the residual
122  call zhemv('U',nmatp,zone,h,nmatp,u(:,ist),1,zzero,hu(:,ist),1)
123  end do
124 !$OMP END DO
125 ! compute the Hamiltonian and overlap matrices in the subspace formed by the
126 ! eigenvectors and their residuals
127 !$OMP DO SCHEDULE(DYNAMIC)
128  do ist=1,nstfv
129  call zgemv('C',nmatp,nstfv,zone,evecfv,nmatmax,hv(:,ist),1,zzero, &
130  hs(1,ist),1)
131  call zgemv('C',nmatp,nstfv,zone,evecfv,nmatmax,hu(:,ist),1,zzero, &
132  hs(1,nstfv+ist),1)
133  call zgemv('C',nmatp,nstfv,zone,u,nmatp,hu(:,ist),1,zzero, &
134  hs(nstfv+1,nstfv+ist),1)
135  end do
136 !$OMP END DO NOWAIT
137 !$OMP DO SCHEDULE(DYNAMIC)
138  do ist=1,nstfv
139  call zgemv('C',nmatp,nstfv,zone,evecfv,nmatmax,ov(:,ist),1,zzero, &
140  os(1,ist),1)
141  call zgemv('C',nmatp,nstfv,zone,evecfv,nmatmax,ou(:,ist),1,zzero, &
142  os(1,nstfv+ist),1)
143  call zgemv('C',nmatp,nstfv,zone,u,nmatp,ou(:,ist),1,zzero, &
144  os(nstfv+1,nstfv+ist),1)
145  end do
146 !$OMP END DO
147 ! solve the generalised eigenvalue problem in the subspace (one thread only)
148 !$OMP SINGLE
149  call zhegv(1,'V','U',ns,hs,ns,os,ns,w,work,lwork,rwork,info)
150 !$OMP END SINGLE
151  if (info /= 0) exit
152 ! construct the new eigenvectors
153 !$OMP DO SCHEDULE(DYNAMIC)
154  do ist=1,nstfv
155  call zgemv('N',nmatp,nstfv,zone,evecfv,nmatmax,hs(1,ist),1,zzero, &
156  ov(:,ist),1)
157  call zgemv('N',nmatp,nstfv,zone,u,nmatp,hs(nstfv+1,ist),1,zone,ov(:,ist),1)
158  end do
159 !$OMP END DO
160 !$OMP DO SCHEDULE(DYNAMIC)
161  do ist=1,nstfv
162  call zcopy(nmatp,ov(:,ist),1,evecfv(:,ist),1)
163  end do
164 !$OMP END DO
165 ! end iteration loop
166 end do
167 !$OMP END PARALLEL
168 call freethd(nthd)
169 deallocate(w,rwork,h,o,hv,ov)
170 deallocate(u,hu,ou,hs,os,work)
171 call timesec(ts1)
172 !$OMP ATOMIC
173 timefv=timefv+ts1-ts0
174 end subroutine
175 
pure subroutine hmlistl(ngp, igpig, vgpc, ld, h)
Definition: hmlistl.f90:10
character(256) filext
Definition: modmain.f90:1299
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
integer iscl
Definition: modmain.f90:1049
Definition: modomp.f90:6
subroutine hmlaa(thr, is, ias, ngp, apwalm, ld, h)
Definition: hmlaa.f90:10
complex(8), parameter zone
Definition: modmain.f90:1238
real(8) timemat
Definition: modmain.f90:1215
logical tefvr
Definition: modmain.f90:870
subroutine olpaa(tor, is, ngp, apwalm, ld, o)
Definition: olpaa.f90:7
pure subroutine olpalo(is, ias, ngp, apwalm, ld, o)
Definition: olpalo.f90:7
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
complex(8), parameter zzero
Definition: modmain.f90:1238
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine timesec(ts)
Definition: timesec.f90:10
pure subroutine olpistl(ngp, igpig, ld, o)
Definition: olpistl.f90:10
integer nefvit
Definition: modmain.f90:875
real(8) timefv
Definition: modmain.f90:1217
subroutine eveqnit(nmatp, ngp, igpig, vpl, vgpl, vgpc, apwalm, evalfv, evecfv)
Definition: eveqnit.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
logical trdstate
Definition: modmain.f90:682
pure subroutine olplolo(is, ias, ngp, ld, o)
Definition: olplolo.f90:7
pure subroutine hmllolo(is, ias, ngp, ld, h)
Definition: hmllolo.f90:7
pure subroutine hmlalo(is, ias, ngp, apwalm, ld, h)
Definition: hmlalo.f90:7