File: gradzf.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (45 lines) | stat: -rw-r--r-- 1,122 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
36
37
38
39
40
41
42
43
44
45

! Copyright (C) 2012 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.

subroutine gradzf(zfmt,zfir,gzfmt,gzfir)
use modmain
implicit none
! arguments
complex(8), intent(in) :: zfmt(npmtmax,natmtot),zfir(ngtot)
complex(8), intent(out) :: gzfmt(npmtmax,natmtot,3),gzfir(ngtot,3)
! local variables
integer is,ias,np
integer ig,ifg,i
complex(8) z1
! allocatable arrays
complex(8), allocatable :: gzfmt1(:,:),zfft(:)
! muffin-tin gradient
allocate(gzfmt1(npmtmax,3))
do ias=1,natmtot
  is=idxis(ias)
  np=npmt(is)
  call gradzfmt(nrmt(is),nrmti(is),rsp(:,is),zfmt(:,ias),npmtmax,gzfmt1)
  do i=1,3
    gzfmt(1:np,ias,i)=gzfmt1(1:np,i)
  end do
end do
deallocate(gzfmt1)
! interstitial gradient
allocate(zfft(ngtot))
call zcopy(ngtot,zfir,1,zfft,1)
call zfftifc(3,ngridg,-1,zfft)
do i=1,3
  gzfir(:,i)=0.d0
  do ig=1,ngvec
    ifg=igfft(ig)
    z1=zfft(ifg)
    gzfir(ifg,i)=vgc(i,ig)*cmplx(-aimag(z1),dble(z1),8)
  end do
  call zfftifc(3,ngridg,1,gzfir(:,i))
end do
deallocate(zfft)
return
end subroutine