The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine genhmlt(ik,vmt,vir,bmt,bir,kmat,pmat,h)
7use modmain
8use modtddft
9use modmpi
10implicit none
11! arguments
12integer, intent(in) :: ik
13real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtc)
14real(8), intent(in) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
15complex(8), intent(in) :: kmat(nstsv,nstsv),pmat(nstsv,nstsv,3)
16complex(8), intent(out) :: h(nstsv,nstsv)
17! local variables
18integer jst,i
19real(8) ca,t1
20! allocatable arrays
21complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
22complex(4), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:)
23allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
24! get the ground-state eigenvectors from file for input k-point
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 input k-point
31allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
32call genwfsv_sp(.false.,.true.,nstsv,[0],ngridg,igfft,ngk(:,ik),igkig(:,:,ik), &
33 apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
34deallocate(apwalm,evecfv)
35! Kohn-Sham potential and magnetic field matrix elements
36if (spinpol) then
37 call genvbmatk(vmt,vir,bmt,bir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,h)
38else
39 call genvmatk(vmt,vir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,h)
40end if
41deallocate(wfmt,wfgk)
42! add the kinetic matrix elements in the second-variational basis
43do jst=1,nstsv
44 h(1:jst,jst)=h(1:jst,jst)+kmat(1:jst,jst)
45end do
46! coupling constant of the external A-field (-1/c)
47ca=-1.d0/solsc
48! add the A-field matrix elements in the second-variational basis
49do 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
56end do
57! add the spin-polarised A-field if required
58if (tafspt) call genhafspt(evecsv,pmat,h)
59deallocate(evecsv)
60end subroutine
61
subroutine genhafspt(evecsv, pmat, h)
Definition genhafspt.f90:7
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
subroutine genvmatk(vmt, vir, ngp, igpig, wfmt, ld, wfgp, vmat)
Definition genvmatk.f90:7
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv_sp.f90:8
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
integer, dimension(3) ngridg
Definition modmain.f90:386
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
integer nspinor
Definition modmain.f90:267
logical spinpol
Definition modmain.f90:228
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
integer lmmaxapw
Definition modmain.f90:199
integer apwordmax
Definition modmain.f90:760
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
integer nmatmax
Definition modmain.f90:848
real(8) solsc
Definition modmain.f90:1252
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
logical tafspt
Definition modtddft.f90:60
real(8), dimension(:,:), allocatable afieldt
Definition modtddft.f90:58
integer itimes
Definition modtddft.f90:46