The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine gradzf(zfmt,zfir,gzfmt,gzfir)
7use modmain
8use modomp
9implicit none
10! arguments
11complex(8), intent(in) :: zfmt(npmtmax,natmtot),zfir(ngtot)
12complex(8), intent(out) :: gzfmt(npmtmax,natmtot,3),gzfir(ngtot,3)
13! local variables
14integer is,ias,ld,i
15integer ig,ifg,nthd
16! allocatable arrays
17complex(8), allocatable :: zfft(:)
18! muffin-tin gradient
19ld=npmtmax*natmtot
20call 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)
25do 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))
29end do
30!$OMP END DO NOWAIT
31! interstitial gradient
32!$OMP SINGLE
33allocate(zfft(ngtot))
34zfft(1:ngtot)=zfir(1:ngtot)
35call zfftifc(3,ngridg,-1,zfft)
36do i=1,3
37 do ig=1,ngvec
38 ifg=igfft(ig)
39 gzfir(ifg,i)=vgc(i,ig)*zi*zfft(ifg)
40 end do
41 gzfir(igfft(ngvec+1:ngtot),i)=0.d0
42 call zfftifc(3,ngridg,1,gzfir(:,i))
43end do
44deallocate(zfft)
45!$OMP END SINGLE
46!$OMP END PARALLEL
47call freethd(nthd)
48end subroutine
49
subroutine gradzf(zfmt, zfir, gzfmt, gzfir)
Definition gradzf.f90:7
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer, dimension(3) ngridg
Definition modmain.f90:386
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
complex(8), parameter zi
Definition modmain.f90:1239
integer ngvec
Definition modmain.f90:396
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine zfftifc(nd, n, sgn, z)