The Elk Code
 
Loading...
Searching...
No Matches
gradrf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2006 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
6subroutine gradrf(rfmt,rfir,grfmt,grfir)
7use modmain
8use modomp
9implicit none
10! arguments
11real(8), intent(in) :: rfmt(npmtmax,natmtot),rfir(ngtot)
12real(8), intent(out) :: grfmt(npmtmax,natmtot,3),grfir(ngtot,3)
13! local variables
14integer is,ias,ld,i
15integer ig,ifg,nthd
16! allocatable arrays
17complex(8), allocatable :: zfft1(:),zfft2(:)
18! muffin-tin gradient
19ld=npmtmax*natmtot
20call holdthd(natmtot+1,nthd)
21!$OMP PARALLEL DEFAULT(SHARED) &
22!$OMP PRIVATE(zfft1,zfft2) &
23!$OMP PRIVATE(is,i,ifg,ig) &
24!$OMP NUM_THREADS(nthd)
25!$OMP DO SCHEDULE(DYNAMIC)
26do ias=1,natmtot
27 is=idxis(ias)
28 call gradrfmt(nrmt(is),nrmti(is),rlmt(:,-1,is),wcrmt(:,:,is),rfmt(:,ias),ld, &
29 grfmt(1,ias,1))
30end do
31!$OMP END DO NOWAIT
32! interstitial gradient
33!$OMP SINGLE
34allocate(zfft1(nfgrz),zfft2(nfgrz))
35call rzfftifc(3,ngridg,-1,rfir,zfft1)
36do i=1,3
37 do ifg=1,nfgrz
38 ig=igrzf(ifg)
39 if (ig <= ngvec) then
40 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
41 else
42 zfft2(ifg)=0.d0
43 end if
44 end do
45 call rzfftifc(3,ngridg,1,grfir(:,i),zfft2)
46end do
47deallocate(zfft1,zfft2)
48!$OMP END SINGLE
49!$OMP END PARALLEL
50call freethd(nthd)
51end subroutine
52
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition gradrf.f90:7
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition gradrfmt.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 nfgrz
Definition modmain.f90:412
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
complex(8), parameter zi
Definition modmain.f90:1239
integer ngvec
Definition modmain.f90:396
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
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 rzfftifc(nd, n, sgn, r, z)