The Elk Code
 
Loading...
Searching...
No Matches
hmlistl.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: hmlistl
8! !INTERFACE:
9pure subroutine hmlistl(ngp,igpig,vgpc,ld,h)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ngp : number of G+p-vectors (in,integer)
14! igpig : index from G+p-vectors to G-vectors (in,integer(ngkmax))
15! vgpc : G+p-vectors in Cartesian coordinates (in,real(3,ngkmax))
16! ld : leading dimension of h (in,integer)
17! h : Hamiltonian matrix (inout,complex(*))
18! !DESCRIPTION:
19! Computes the interstitial contribution to the Hamiltonian matrix for the APW
20! basis functions. The Hamiltonian is given by
21! $$ H^{\rm I}({\bf G+k,G'+k})=\frac{1}{2}({\bf G+k})\cdot({\bf G'+k})
22! \tilde{\Theta}({\bf G-G'})+V_s({\bf G-G'}), $$
23! where $V_s$ is the interstitial Kohn-Sham potential and $\tilde{\Theta}$ is
24! the characteristic function. See routine {\tt gencfun}.
25!
26! !REVISION HISTORY:
27! Created April 2003 (JKD)
28!EOP
29!BOC
30implicit none
31! arguments
32integer, intent(in) :: ngp,igpig(ngkmax)
33real(8), intent(in) :: vgpc(3,ngkmax)
34integer, intent(in) :: ld
35complex(8), intent(out) :: h(ld,*)
36! local variables
37integer j1,j2,j3,ig,i,j
38real(8) v1,v2,v3
39do j=1,ngp
40 ig=igpig(j)
41 j1=ivg(1,ig); j2=ivg(2,ig); j3=ivg(3,ig)
42 v1=0.5d0*vgpc(1,j); v2=0.5d0*vgpc(2,j); v3=0.5d0*vgpc(3,j)
43 do i=1,j
44 ig=igpig(i)
45 ig=ivgig(ivg(1,ig)-j1,ivg(2,ig)-j2,ivg(3,ig)-j3)
46 h(i,j)=vsig(ig)+(vgpc(1,i)*v1+vgpc(2,i)*v2+vgpc(3,i)*v3)*cfunig(ig)
47 end do
48end do
49end subroutine
50!EOC
51
pure subroutine hmlistl(ngp, igpig, vgpc, ld, h)
Definition hmlistl.f90:10
complex(8), dimension(:), allocatable cfunig
Definition modmain.f90:434
complex(8), dimension(:), allocatable vsig
Definition modmain.f90:662
integer, dimension(:,:), allocatable ivg
Definition modmain.f90:400
integer ngkmax
Definition modmain.f90:499
integer, dimension(:,:,:), allocatable ivgig
Definition modmain.f90:402