The Elk Code
Loading...
Searching...
No Matches
rfirftoc.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
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
rfirftoc
(rfir,rfirc)
7
use
modmain
8
implicit none
9
! arguments
10
real
(8),
intent(in)
:: rfir(ngtot)
11
real
(8),
intent(out)
:: rfirc(ngtc)
12
! local variables
13
integer
ig,ifg
14
! allocatable arrays
15
complex(8)
,
allocatable
:: zfft(:)
16
! automatic arrays
17
complex(8)
zfftc(nfgrzc)
18
allocate
(zfft(ngtot))
19
! multiply by characteristic function and Fourier transform on fine grid
20
zfft(1:ngtot)=rfir(1:ngtot)*
cfunir
(1:ngtot)
21
call
zfftifc
(3,
ngridg
,-1,zfft)
22
! Fourier transform to coarse real-space grid
23
do
ifg=1,nfgrzc
24
ig=
igrzfc
(ifg)
25
if
(ig <=
ngvc
)
then
26
zfftc(ifg)=zfft(
igfft
(ig))
27
else
28
zfftc(ifg)=0.d0
29
end if
30
end do
31
deallocate
(zfft)
32
! Fourier transform to real-space coarse grid
33
call
rzfftifc
(3,
ngdgc
,1,rfirc,zfftc)
34
end subroutine
35
modmain
Definition
modmain.f90:6
modmain::ngridg
integer, dimension(3) ngridg
Definition
modmain.f90:386
modmain::ngdgc
integer, dimension(3) ngdgc
Definition
modmain.f90:388
modmain::cfunir
real(8), dimension(:), allocatable cfunir
Definition
modmain.f90:436
modmain::igfft
integer, dimension(:), allocatable igfft
Definition
modmain.f90:406
modmain::igrzfc
integer, dimension(:), allocatable igrzfc
Definition
modmain.f90:418
modmain::ngvc
integer ngvc
Definition
modmain.f90:398
rfirftoc
subroutine rfirftoc(rfir, rfirc)
Definition
rfirftoc.f90:7
zfftifc
subroutine zfftifc(nd, n, sgn, z)
Definition
zfftifc_fftw.f90:7
rzfftifc
subroutine rzfftifc(nd, n, sgn, r, z)
Definition
zfftifc_fftw.f90:32
rfirftoc.f90
Generated by
1.9.8