The Elk Code
plotulr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 T. Mueller, 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 plotulr(np,vpl,nf,zfmt,zfir,fp)
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! arguments
12 integer, intent(in) :: np
13 real(8), intent(in) :: vpl(3,np)
14 integer, intent(in) :: nf
15 complex(8), intent(in) :: zfmt(npcmtmax,natmtot,nf,nfqrz)
16 complex(8), intent(in) :: zfir(ngtot,nf,nfqrz)
17 real(8), intent(out) :: fp(np,nf)
18 ! local variables
19 integer iq,ifq0,ifq
20 integer jf,ip
21 real(8) sm,t1
22 complex(8) z1
23 ! allocatable arrays
24 complex(8), allocatable :: fpq(:,:)
25 allocate(fpq(np,nfqrz))
26 ! include or exclude the Q=0 component as required
27 if (tplotq0) then
28  ifq0=1
29 else
30  ifq0=2
31 end if
32 ! loop over the number of functions
33 do jf=1,nf
34 ! loop over real-complex FFT points
35  do ifq=ifq0,nfqrz
36 ! evaluate the complex function at all the plot points
37  call zfpts(np,vpl,zfmt(:,:,jf,ifq),zfir(:,jf,ifq),fpq(:,ifq))
38  end do
39  do ip=1,np
40  sm=0.d0
41  do ifq=ifq0,nfqrz
42  iq=iqrzf(ifq)
43 ! multiply complex function by phase factor exp(iQ.r)
44  t1=twopi*dot_product(vql(:,iq),vpl(:,ip))
45  z1=cmplx(cos(t1),sin(t1),8)
46  t1=dble(fpq(ip,ifq)*z1)
47  if (ifq > 1) t1=t1*2.d0
48  sm=sm+t1
49  end do
50  fp(ip,jf)=sm
51  end do
52 end do
53 deallocate(fpq)
54 end subroutine
55 
real(8), parameter twopi
Definition: modmain.f90:1231
subroutine plotulr(np, vpl, nf, zfmt, zfir, fp)
Definition: plotulr.f90:7
Definition: modomp.f90:6
integer, dimension(:), allocatable iqrzf
Definition: modmain.f90:543
subroutine zfpts(np, vrl, zfmt, zfir, fp)
Definition: zfpts.f90:7
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
Definition: modulr.f90:6
logical tplotq0
Definition: modulr.f90:87