The Elk Code
hdbulrk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 T. Mueller, 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 hdbulrk(ik0,hdb)
7 use modmain
8 use modulr
9 implicit none
10 ! arguments
11 integer, intent(in) :: ik0
12 complex(8), intent(out) :: hdb(nstsv,nstsv,2:nkpa)
13 ! local variables
14 integer ik,ikk,ikpa,ist
15 real(8) t1
16 ! allocatable arrays
17 complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
18 complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
19 complex(8), allocatable :: wfmtk(:,:,:,:),wfgkk(:,:,:)
20 complex(8), allocatable :: ok(:,:),b(:,:)
21 ! central k-point
22 ik=(ik0-1)*nkpa+1
23 ! get the ground-state eigenvectors from file for central k-point
24 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
25 call getevecfv('.OUT',ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
26 call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
27 ! find the matching coefficients
28 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
29 call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
30 ! calculate the wavefunctions for all states of the central k-point
31 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfir(ngtc,nspinor,nstsv))
32 call genwfsv(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ik),igkig(:,1,ik), &
33  apwalm,evecfv,evecsv,wfmt,ngtc,wfir)
34 ! compute the diagonal blocks of the ultra long-range Hamiltonian in the basis
35 ! of the states at the central k-point
36 allocate(wfmtk(npcmtmax,natmtot,nspinor,nstsv),wfgkk(ngkmax,nspinor,nstsv))
37 allocate(ok(nstsv,nstsv),b(nstsv,nstsv))
38 do ikpa=2,nkpa
39  ikk=(ik0-1)*nkpa+ikpa
40  call getevecfv('.OUT',ikk,vkl(:,ikk),vgkl(:,:,:,ikk),evecfv)
41  call getevecsv('.OUT',ikk,vkl(:,ikk),evecsv)
42  call match(ngk(1,ikk),vgkc(:,:,1,ikk),gkc(:,1,ikk),sfacgk(:,:,1,ikk),apwalm)
43  call genwfsv(.false.,.true.,nstsv,[0],ngdgc,igfc,ngk(:,ikk),igkig(:,:,ikk), &
44  apwalm,evecfv,evecsv,wfmtk,ngkmax,wfgkk)
45 ! compute the overlap matrix between the states at k and k+κ
46  call genolpq(nstsv,expqmt(:,:,ikpa),ngk(:,ikk),igkig(:,:,ikk),wfmt,wfir, &
47  wfmtk,wfgkk,ok)
48 ! use singular value decompostion to make the matrix strictly unitary
49  call unitary(nstsv,ok)
50 ! apply the overlap matrix from the right to the eigenvalues at k+κ
51  do ist=1,nstsv
52  t1=evalsv(ist,ikk)
53  b(ist,1:nstsv)=t1*ok(ist,1:nstsv)
54  end do
55 ! apply the conjugate transpose of the overlap matrix from the left to form the
56 ! Hamiltonian matrix in the basis of states at k
57  call zgemm('C','N',nstsv,nstsv,nstsv,zone,ok,nstsv,b,nstsv,zzero, &
58  hdb(:,:,ikpa),nstsv)
59 end do
60 deallocate(apwalm,evecfv,evecsv)
61 deallocate(wfmt,wfir,wfmtk,wfgkk)
62 deallocate(ok,b)
63 end subroutine
64 
integer nmatmax
Definition: modmain.f90:853
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer npcmtmax
Definition: modmain.f90:216
subroutine genwfsv(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv.f90:11
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
subroutine genolpq(nst, expqmt, ngpq, igpqig, wfmt, wfir, wfmtq, wfgpq, oq)
Definition: genolpq.f90:7
integer ngtc
Definition: modmain.f90:392
integer lmmaxapw
Definition: modmain.f90:199
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
complex(8), parameter zone
Definition: modmain.f90:1240
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine hdbulrk(ik0, hdb)
Definition: hdbulrk.f90:7
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
integer natmtot
Definition: modmain.f90:40
subroutine unitary(n, a)
Definition: unitary.f90:10
Definition: modulr.f90:6
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887
complex(8), dimension(:,:,:), allocatable expqmt
Definition: modulr.f90:46