The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine hdbulrk(ik0,hdb)
7use modmain
8use modulr
9implicit none
10! arguments
11integer, intent(in) :: ik0
12complex(8), intent(out) :: hdb(nstsv,nstsv,2:nkpa)
13! local variables
14integer ik,ikk,ikpa,ist
15real(8) t1
16! allocatable arrays
17complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
18complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
19complex(8), allocatable :: wfmtk(:,:,:,:),wfgkk(:,:,:)
20complex(8), allocatable :: ok(:,:),b(:,:)
21! central k-point
22ik=(ik0-1)*nkpa+1
23! get the ground-state eigenvectors from file for central k-point
24allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
25call getevecfv('.OUT',ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
26call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
27! find the matching coefficients
28allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
29call 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
31allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfir(ngtc,nspinor,nstsv))
32call 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
36allocate(wfmtk(npcmtmax,natmtot,nspinor,nstsv),wfgkk(ngkmax,nspinor,nstsv))
37allocate(ok(nstsv,nstsv),b(nstsv,nstsv))
38do 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)
59end do
60deallocate(apwalm,evecfv,evecsv)
61deallocate(wfmt,wfir,wfmtk,wfgkk)
62deallocate(ok,b)
63end subroutine
64
subroutine genolpq(nst, expqmt, ngpq, igpqig, wfmt, wfir, wfmtq, wfgpq, oq)
Definition genolpq.f90:7
subroutine genwfsv(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv.f90:11
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine hdbulrk(ik0, hdb)
Definition hdbulrk.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
complex(8), parameter zzero
Definition modmain.f90:1238
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
integer nspinor
Definition modmain.f90:267
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
complex(8), parameter zone
Definition modmain.f90:1238
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer lmmaxapw
Definition modmain.f90:199
integer apwordmax
Definition modmain.f90:760
integer ngtc
Definition modmain.f90:392
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer npcmtmax
Definition modmain.f90:216
integer nmatmax
Definition modmain.f90:848
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
integer nstfv
Definition modmain.f90:884
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
complex(8), dimension(:,:,:), allocatable expqmt
Definition modulr.f90:46
subroutine unitary(n, a)
Definition unitary.f90:10