The Elk Code
 
Loading...
Searching...
No Matches
zfirctof.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
6subroutine zfirctof(zfirc,zfir)
7use modmain
8implicit none
9! arguments
10complex(8), intent(in) :: zfirc(ngtc)
11complex(8), intent(out) :: zfir(ngtot)
12! local variables
13integer ig
14! automatic arrays
15complex(8) zfftc(ngtc)
16! Fourier transform function on coarse grid to G-space
17zfftc(1:ngtc)=zfirc(1:ngtc)
18call zfftifc(3,ngdgc,-1,zfftc)
19! Fourier transform to fine real-space grid
20do ig=1,ngvc
21 zfir(igfft(ig))=zfftc(igfc(ig))
22end do
23zfir(igfft(ngvc+1:ngtot))=0.d0
24call zfftifc(3,ngridg,1,zfir)
25end subroutine
26
integer, dimension(3) ngridg
Definition modmain.f90:386
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
integer ngvc
Definition modmain.f90:398
subroutine zfftifc(nd, n, sgn, z)
subroutine zfirctof(zfirc, zfir)
Definition zfirctof.f90:7