The Elk Code
 
Loading...
Searching...
No Matches
rfirctof.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!BOP
7! !ROUTINE: rfirctof
8! !INTERFACE:
9subroutine rfirctof(rfirc,rfir)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! rfirc : real input function on coarse grid (in,real(ngtc))
14! rfir : real output function on fine grid (out,real(ngtot))
15! !DESCRIPTION:
16! Converts a real function on a coarse grid given by sizes {\tt ngdgc} to a
17! function on a fine grid given by {\tt ngridg}. This is done by first
18! Fourier transforming {\tt rfirc} to $G$-space, zeroing the missing values
19! and then transforming back to {\tt rfir}.
20!
21! !REVISION HISTORY:
22! Created March 2020 (JKD)
23!EOP
24!BOC
25implicit none
26! arguments
27real(8), intent(in) :: rfirc(ngtc)
28real(8), intent(out) :: rfir(ngtot)
29! local variables
30integer ig,ifg
31! automatic arrays
32complex(8) zfftc(ngtc)
33! allocatable arrays
34complex(8), allocatable :: zfft(:)
35! Fourier transform function on coarse grid to G-space
36zfftc(1:ngtc)=rfirc(1:ngtc)
37call zfftifc(3,ngdgc,-1,zfftc)
38! Fourier transform to fine real-space grid
39allocate(zfft(nfgrz))
40do ifg=1,nfgrz
41 ig=igrzf(ifg)
42 if (ig <= ngvc) then
43 zfft(ifg)=zfftc(igfc(ig))
44 else
45 zfft(ifg)=0.d0
46 end if
47end do
48call rzfftifc(3,ngridg,1,rfir,zfft)
49deallocate(zfft)
50end subroutine
51!EOC
52
integer, dimension(3) ngridg
Definition modmain.f90:386
integer nfgrz
Definition modmain.f90:412
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer ngvc
Definition modmain.f90:398
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine zfftifc(nd, n, sgn, z)
subroutine rzfftifc(nd, n, sgn, r, z)