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