The Elk Code
genhmlu.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 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 genhmlu(ik0,h)
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik0
13 complex(8), intent(out) :: h(nstulr,nstulr)
14 ! local variables
15 integer ik,ikk,ist,jst,ispn,idm
16 integer ikpa,jkpa,iq,ifq,ngk0,igk
17 integer i1,i2,i3,j1,j2,j3,i,j,nthd
18 ! automatic arrays
19 complex(8) zvir(ngtc),zbir(ngtc,ndmag),vmat(nstsv,nstsv)
20 ! allocatable arrays
21 complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
22 complex(4), allocatable :: wfmt(:,:,:,:),wfir(:,:,:),wfgk(:,:,:)
23 ! central k-point
24 ik=(ik0-1)*nkpa+1
25 ! number of G+k-vectors for central k-point
26 ngk0=ngk(1,ik)
27 ! get the ground-state eigenvectors from file for central k-point
28 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
29 call getevecfv('.OUT',ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
30 call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
31 ! find the matching coefficients
32 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
33 call match(ngk0,vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
34 ! calculate the wavefunctions for all states of the central k-point
35 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngk0,nspinor,nstsv))
36 call genwfsv_sp(.false.,.true.,nstsv,[0],ngridg,igfft,ngk0,igkig(:,1,ik), &
37  apwalm,evecfv,evecsv,wfmt,ngk0,wfgk)
38 deallocate(apwalm,evecfv,evecsv)
39 ! determine the interstitial wavefunctions in real-space (without 1/sqrt(omega))
40 allocate(wfir(ngtc,nspinor,nstsv))
41 call holdthd(nstsv,nthd)
42 !$OMP PARALLEL DO DEFAULT(SHARED) &
43 !$OMP PRIVATE(ispn,igk) &
44 !$OMP SCHEDULE(DYNAMIC) &
45 !$OMP NUM_THREADS(nthd)
46 do ist=1,nstsv
47  do ispn=1,nspinor
48  wfir(1:ngtc,ispn,ist)=0.e0
49  do igk=1,ngk0
50  wfir(igfc(igkig(igk,1,ik)),ispn,ist)=wfgk(igk,ispn,ist)
51  end do
52  call cfftifc(3,ngdgc,1,wfir(:,ispn,ist))
53  end do
54 end do
55 !$OMP END PARALLEL DO
56 call freethd(nthd)
57 ! generate the matrix elements for all Q-vectors
58 call holdthd(nfqrz,nthd)
59 !$OMP PARALLEL DO DEFAULT(SHARED) &
60 !$OMP PRIVATE(zvir,zbir,vmat) &
61 !$OMP PRIVATE(iq,idm,ikpa,jkpa,jst) &
62 !$OMP PRIVATE(i1,i2,i3,j1,j2,j3,i,j) &
63 !$OMP SCHEDULE(DYNAMIC) &
64 !$OMP NUM_THREADS(nthd)
65 do ifq=1,nfqrz
66  iq=iqrzf(ifq)
67 ! multiply long-range interstitial potential by characteristic function and
68 ! convert to coarse grid
69  call zfirftoc(vsqir(:,ifq),zvir)
70  if (spinpol) then
71 ! convert interstitial magnetic field to coarse grid
72  do idm=1,ndmag
73  call zfirftoc(bsqir(:,idm,ifq),zbir(:,idm))
74  end do
75 ! calculate matrix elements for this Q-vector
76  call genzvbmatk(vsqmt(:,:,ifq),zvir,bsqmt(:,:,:,ifq),zbir,ngk0, &
77  igkig(:,1,ik),wfmt,wfir,wfgk,vmat)
78  else
79  call genzvmatk(vsqmt(:,:,ifq),zvir,ngk0,igkig(:,1,ik),wfmt,wfir,wfgk,vmat)
80  end if
81  do jkpa=1,nkpa
82  j1=ivq(1,jkpa); j2=ivq(2,jkpa); j3=ivq(3,jkpa)
83  do jst=1,nstsv
84  j=(jkpa-1)*nstsv+jst
85  do ikpa=1,jkpa-1
86  i=(ikpa-1)*nstsv+1
87  i1=ivq(1,ikpa)-j1; i2=ivq(2,ikpa)-j2; i3=ivq(3,ikpa)-j3
88  if (ivqiq(i1,i2,i3) == iq) then
89 ! copy matrix elements for κ_i - κ_j in Q-point set
90  h(i:i+nstsv-1,j)=vmat(1:nstsv,jst)
91  else if (ivqiq(-i1,-i2,-i3) == iq) then
92 ! otherwise use conjugate transpose
93  h(i:i+nstsv-1,j)=conjg(vmat(jst,1:nstsv))
94  end if
95  end do
96 ! copy only the upper triangular part for Q = 0
97  if (ifq == 1) then
98  i=(jkpa-1)*nstsv+1
99  h(i:i+jst-1,j)=vmat(1:jst,jst)
100  end if
101  end do
102  end do
103 end do
104 !$OMP END PARALLEL DO
105 call freethd(nthd)
106 ! add the second-variational eigenvalues of k+κ to the diagonal
107 do ikpa=1,nkpa
108  ikk=(ik0-1)*nkpa+ikpa
109  i=(ikpa-1)*nstsv
110  do ist=1,nstsv
111  j=i+ist
112  h(j,j)=h(j,j)+evalsv(ist,ikk)
113  end do
114 end do
115 deallocate(wfmt,wfir,wfgk)
116 end subroutine
117 
integer nmatmax
Definition: modmain.f90:853
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition: modulr.f90:82
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer npcmtmax
Definition: modmain.f90:216
integer, dimension(3) ngridg
Definition: modmain.f90:386
subroutine zfirftoc(zfir, zfirc)
Definition: zfirftoc.f90:7
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:919
logical spinpol
Definition: modmain.f90:228
integer lmmaxapw
Definition: modmain.f90:199
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition: modulr.f90:82
integer, dimension(:,:), allocatable ivq
Definition: modmain.f90:529
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
Definition: modomp.f90:6
integer, dimension(:), allocatable iqrzf
Definition: modmain.f90:543
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition: modulr.f90:83
subroutine genzvbmatk(zvmt, zvir, zbmt, zbir, ngp, igpig, wfmt, wfir, wfgp, vbmat)
Definition: genzvbmatk.f90:7
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine cfftifc(nd, n, sgn, c)
Definition: cfftifc_fftw.f90:7
subroutine genzvmatk(zvmt, zvir, ngp, igpig, wfmt, wfir, wfgp, vmat)
Definition: genzvmatk.f90:7
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv_sp.f90:8
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer nspinor
Definition: modmain.f90:267
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition: modulr.f90:83
integer apwordmax
Definition: modmain.f90:760
integer nfqrz
Definition: modmain.f90:539
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
Definition: modulr.f90:6
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
subroutine genhmlu(ik0, h)
Definition: genhmlu.f90:7
integer nkpa
Definition: modulr.f90:24
integer nstfv
Definition: modmain.f90:885