The Elk Code
gradzf.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 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 gradzf(zfmt,zfir,gzfmt,gzfir)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 complex(8), intent(in) :: zfmt(npmtmax,natmtot),zfir(ngtot)
12 complex(8), intent(out) :: gzfmt(npmtmax,natmtot,3),gzfir(ngtot,3)
13 ! local variables
14 integer is,ias,ld,i
15 integer ig,ifg,nthd
16 ! allocatable arrays
17 complex(8), allocatable :: zfft(:)
18 ! muffin-tin gradient
19 ld=npmtmax*natmtot
20 call holdthd(natmtot+1,nthd)
21 !$OMP PARALLEL DEFAULT(SHARED) &
22 !$OMP PRIVATE(zfft,is,i,ig,ifg) &
23 !$OMP NUM_THREADS(nthd)
24 !$OMP DO SCHEDULE(DYNAMIC)
25 do ias=1,natmtot
26  is=idxis(ias)
27  call gradzfmt(nrmt(is),nrmti(is),rlmt(:,-1,is),wcrmt(:,:,is),zfmt(:,ias),ld, &
28  gzfmt(1,ias,1))
29 end do
30 !$OMP END DO NOWAIT
31 ! interstitial gradient
32 !$OMP SINGLE
33 allocate(zfft(ngtot))
34 zfft(1:ngtot)=zfir(1:ngtot)
35 call zfftifc(3,ngridg,-1,zfft)
36 do i=1,3
37  do ig=1,ngvec
38  ifg=igfft(ig)
39  gzfir(ifg,i)=vgc(i,ig)*cmplx(-zfft(ifg)%im,zfft(ifg)%re,8)
40  end do
41  gzfir(igfft(ngvec+1:ngtot),i)=0.d0
42  call zfftifc(3,ngridg,1,gzfir(:,i))
43 end do
44 deallocate(zfft)
45 !$OMP END SINGLE
46 !$OMP END PARALLEL
47 call freethd(nthd)
48 end subroutine
49 
integer, dimension(3) ngridg
Definition: modmain.f90:386
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
Definition: modomp.f90:6
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
integer ngvec
Definition: modmain.f90:396
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine gradzf(zfmt, zfir, gzfmt, gzfir)
Definition: gradzf.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150