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
ifq0=merge(1,2,
tplotq0
)
28
! loop over the number of functions
29
do
jf=1,nf
30
! loop over real-complex FFT points
31
do
ifq=ifq0,nfqrz
32
! evaluate the complex function at all the plot points
33
call
zfpts
(np,vpl,zfmt(:,:,jf,ifq),zfir(:,jf,ifq),fpq(:,ifq))
34
end do
35
do
ip=1,np
36
sm=0.d0
37
do
ifq=ifq0,nfqrz
38
iq=
iqrzf
(ifq)
39
! multiply complex function by phase factor exp(iQ.r)
40
t1=
twopi
*dot_product(
vql
(:,iq),vpl(:,ip))
41
z1=cmplx(cos(t1),sin(t1),8)
42
t1=dble(fpq(ip,ifq)*z1)
43
if
(ifq > 1) t1=t1*2.d0
44
sm=sm+t1
45
end do
46
fp(ip,jf)=sm
47
end do
48
end do
49
deallocate
(fpq)
50
end subroutine
51
modmain::twopi
real(8), parameter twopi
Definition:
modmain.f90:1220
plotulr
subroutine plotulr(np, vpl, nf, zfmt, zfir, fp)
Definition:
plotulr.f90:7
modomp
Definition:
modomp.f90:6
modmain::iqrzf
integer, dimension(:), allocatable iqrzf
Definition:
modmain.f90:543
modmain
Definition:
modmain.f90:6
zfpts
subroutine zfpts(np, vrl, zfmt, zfir, fp)
Definition:
zfpts.f90:7
modmain::vql
real(8), dimension(:,:), allocatable vql
Definition:
modmain.f90:545
modulr
Definition:
modulr.f90:6
modulr::tplotq0
logical tplotq0
Definition:
modulr.f90:89
plotulr.f90
Generated by
1.8.14