File: dhmlistl.f90

package info (click to toggle)
elkcode 10.8.1-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 10,672 kB
  • sloc: f90: 52,747; perl: 964; makefile: 352; sh: 296; python: 105; ansic: 67
file content (35 lines) | stat: -rw-r--r-- 969 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35

! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

pure subroutine dhmlistl(ngp,ngpq,igpig,igpqig,vgpc,vgpqc,ld,dh)
use modmain
use modphonon
implicit none
! arguments
integer, intent(in) :: ngp,ngpq
integer, intent(in) :: igpig(ngkmax),igpqig(ngkmax)
real(8), intent(in) :: vgpc(3,ngkmax),vgpqc(3,ngkmax)
integer, intent(in) :: ld
complex(8), intent(out) :: dh(ld,*)
! local variables
integer j1,j2,j3,ig,i,j
real(8) v1,v2,v3,t1
do j=1,ngp
  ig=igpig(j)
  j1=ivg(1,ig); j2=ivg(2,ig); j3=ivg(3,ig)
  v1=0.5d0*vgpc(1,j); v2=0.5d0*vgpc(2,j); v3=0.5d0*vgpc(3,j)
  do i=1,ngpq
    ig=igpqig(i)
    ig=ivgig(ivg(1,ig)-j1,ivg(2,ig)-j2,ivg(3,ig)-j3)
    if (ig <= ngvc) then
      t1=vgpqc(1,i)*v1+vgpqc(2,i)*v2+vgpqc(3,i)*v3
      dh(i,j)=dvsig(ig)+t1*dcfunig(ig)
    else
      dh(i,j)=0.d0
    end if
  end do
end do
end subroutine