The Elk Code
genhmlt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 K. Krieger, 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 genhmlt(ik,vmt,vir,bmt,bir,kmat,pmat,h)
7 use modmain
8 use modtddft
9 use modmpi
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik
13 real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtc)
14 real(8), intent(in) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
15 complex(8), intent(in) :: kmat(nstsv,nstsv),pmat(nstsv,nstsv,3)
16 complex(8), intent(out) :: h(nstsv,nstsv)
17 ! local variables
18 integer jst,i
19 real(8) ca,t1
20 ! allocatable arrays
21 complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
22 complex(4), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:)
23 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
24 ! get the ground-state eigenvectors from file for input k-point
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 input k-point
31 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
32 call genwfsv_sp(.false.,.true.,nstsv,[0],ngridg,igfft,ngk(:,ik),igkig(:,:,ik), &
33  apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
34 deallocate(apwalm,evecfv)
35 ! Kohn-Sham potential and magnetic field matrix elements
36 if (spinpol) then
37  call genvbmatk(vmt,vir,bmt,bir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,h)
38 else
39  call genvmatk(vmt,vir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,h)
40 end if
41 deallocate(wfmt,wfgk)
42 ! add the kinetic matrix elements in the second-variational basis
43 do jst=1,nstsv
44  h(1:jst,jst)=h(1:jst,jst)+kmat(1:jst,jst)
45 end do
46 ! coupling constant of the external A-field (-1/c)
47 ca=-1.d0/solsc
48 ! add the A-field matrix elements in the second-variational basis
49 do i=1,3
50  t1=ca*afieldt(i,itimes)
51  if (abs(t1) > 1.d-10) then
52  do jst=1,nstsv
53  h(1:jst,jst)=h(1:jst,jst)+t1*pmat(1:jst,jst,i)
54  end do
55  end if
56 end do
57 ! add the spin-polarised A-field if required
58 if (tafspt) call genhafspt(evecsv,pmat,h)
59 deallocate(evecsv)
60 end subroutine
61 
integer nmatmax
Definition: modmain.f90:853
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
subroutine genhafspt(evecsv, pmat, h)
Definition: genhafspt.f90:7
integer, dimension(3) ngridg
Definition: modmain.f90:386
logical spinpol
Definition: modmain.f90:228
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), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine genvmatk(vmt, vir, ngp, igpig, wfmt, ld, wfgp, vmat)
Definition: genvmatk.f90:7
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 nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8) solsc
Definition: modmain.f90:1253
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
Definition: modmpi.f90:6
logical tafspt
Definition: modtddft.f90:60
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer itimes
Definition: modtddft.f90:46
subroutine genhmlt(ik, vmt, vir, bmt, bir, kmat, pmat, h)
Definition: genhmlt.f90:7
subroutine genvbmatk(vmt, vir, bmt, bir, ngp, igpig, wfmt, ld, wfgp, vbmat)
Definition: genvbmatk.f90:7
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887