The Elk Code
 
Loading...
Searching...
No Matches
rfpts.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 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 rfpts(np,vpl,rfmt,rfir,fp)
7use modmain
8use modomp
9implicit none
10! arguments
11integer, intent(in) :: np
12real(8), intent(in) :: vpl(3,np),rfmt(npmtmax,natmtot),rfir(ngtot)
13real(8), intent(out) :: fp(np)
14! local variables
15integer ip,nthd
16! allocatable arrays
17complex(8), allocatable :: zfft(:)
18! Fourier transform rfir to G-space
19allocate(zfft(ngtot))
20zfft(1:ngtot)=rfir(1:ngtot)
21call zfftifc(3,ngridg,-1,zfft)
22! begin loop over all points
23call holdthd(np,nthd)
24!$OMP PARALLEL DO DEFAULT(SHARED) &
25!$OMP SCHEDULE(DYNAMIC) &
26!$OMP NUM_THREADS(nthd)
27do ip=1,np
28 call rfpt(ip)
29end do
30!$OMP END PARALLEL DO
31call freethd(nthd)
32deallocate(zfft)
33return
34
35contains
36
37subroutine rfpt(ip)
38implicit none
39! arguments
40integer, intent(in) :: ip
41! local variables
42integer is,ias,nr,nri
43integer ir,ir0,i0,ig,j
44integer lmax,lmmax,lm
45real(8) v(3),r,sm,t1
46! automatic arrays
47real(8) ya(4,lmmaxo),rlm(lmmaxo)
48! check if point is in a muffin-tin
49call findmtpt(vpl(:,ip),ias,ir,v,r)
50if (ir > 0) then
51! point is inside a muffin-tin
52 is=idxis(ias)
53 nr=nrmt(is)
54 nri=nrmti(is)
55 if (ir <= 3) then
56 ir0=1
57 else if (ir > nr-2) then
58 ir0=nr-3
59 else
60 ir0=ir-2
61 end if
62 r=max(r,rsp(1,is))
63 if (ir0 <= nri) then
64 lmax=lmaxi
65 lmmax=lmmaxi
66 else
67 lmax=lmaxo
68 lmmax=lmmaxo
69 end if
70 do j=1,4
71 ir=ir0+j-1
72 if (ir <= nri) then
73 i0=lmmaxi*(ir-1)
74 else
75 i0=lmmaxi*nri+lmmaxo*(ir-nri-1)
76 end if
77 ya(j,1:lmmax)=rfmt(i0+1:i0+lmmax,ias)
78 end do
79! generate real spherical harmonics
80 call genrlmv(lmax,v,rlm)
81 sm=0.d0
82 do lm=1,lmmax
83 t1=poly4(rsp(ir0,is),ya(:,lm),r)
84 sm=sm+t1*rlm(lm)
85 end do
86else
87! otherwise use direct Fourier transform of the interstitial function
88 v(1:3)=vpl(1,ip)*avec(1:3,1)+vpl(2,ip)*avec(1:3,2)+vpl(3,ip)*avec(1:3,3)
89 sm=0.d0
90!$OMP SIMD PRIVATE(t1) REDUCTION(+:sm)
91 do ig=1,ngvec
92 t1=vgc(1,ig)*v(1)+vgc(2,ig)*v(2)+vgc(3,ig)*v(3)
93 sm=sm+dble(zfft(igfft(ig))*cmplx(cos(t1),sin(t1),8))
94 end do
95end if
96fp(ip)=sm
97end subroutine
98
99pure real(8) function poly4(xa,ya,x)
100implicit none
101! arguments
102real(8), intent(in) :: xa(4),ya(4),x
103! local variables
104real(8) x0,x1,x2,x3,y0,y1,y2,y3
105real(8) c1,c2,c3,t0,t1,t2,t3,t4,t5,t6
106! evaluate the polynomial coefficients
107x0=xa(1)
108x1=xa(2)-x0; x2=xa(3)-x0; x3=xa(4)-x0
109y0=ya(1)
110y1=ya(2)-y0; y2=ya(3)-y0; y3=ya(4)-y0
111t4=x1-x2; t5=x1-x3; t6=x2-x3
112t1=x1*x2*y3; t2=x2*x3*y1; t3=x1*x3
113t0=1.d0/(x2*t3*t4*t5*t6)
114t3=t3*y2
115c3=t1*t4+t2*t6-t3*t5
116t4=x1**2; t5=x2**2; t6=x3**2
117c2=t1*(t5-t4)+t2*(t6-t5)+t3*(t4-t6)
118c1=t1*(x2*t4-x1*t5)+t2*(x3*t5-x2*t6)+t3*(x1*t6-x3*t4)
119t1=x-x0
120! evaluate the polynomial
121poly4=y0+t0*t1*(c1+t1*(c2+c3*t1))
122end function
123
124end subroutine
125
subroutine findmtpt(vrl, ias, ir, v, r)
Definition findmtpt.f90:7
subroutine genrlmv(lmax, v, rlm)
Definition genrlmv.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer, dimension(3) ngridg
Definition modmain.f90:386
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer lmmaxi
Definition modmain.f90:207
integer ngvec
Definition modmain.f90:396
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
integer lmaxo
Definition modmain.f90:201
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer lmaxi
Definition modmain.f90:205
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfpt(ip)
Definition rfpts.f90:38
subroutine rfpts(np, vpl, rfmt, rfir, fp)
Definition rfpts.f90:7
pure real(8) function poly4(xa, ya, x)
Definition rfpts.f90:100
subroutine zfftifc(nd, n, sgn, z)