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