The Elk Code
cfftifc_fftw.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2022 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 cfftifc(nd,n,sgn,c)
7 implicit none
8 ! arguments
9 integer, intent(in) :: nd,n(nd),sgn
10 complex(4), intent(inout) :: c(*)
11 ! local variables
12 integer, parameter :: FFTW_ESTIMATE=64
13 integer p
14 integer(8) plan
15 real(4) t1
16 ! interface to FFTW version 3
17 !$OMP CRITICAL(fftw_)
18 call sfftw_plan_dft(plan,nd,n,c,c,sgn,fftw_estimate)
19 !$OMP END CRITICAL(fftw_)
20 call sfftw_execute(plan)
21 !$OMP CRITICAL(fftw_)
22 call sfftw_destroy_plan(plan)
23 !$OMP END CRITICAL(fftw_)
24 if (sgn == -1) then
25  p=product(n(1:nd))
26  t1=1.e0/real(p)
27  call csscal(p,t1,c,1)
28 end if
29 end subroutine
30 
31 subroutine rcfftifc(nd,n,sgn,r,c)
32 implicit none
33 ! arguments
34 integer, intent(in) :: nd,n(nd),sgn
35 real(4), intent(inout) :: r(*)
36 complex(4), intent(inout) :: c(*)
37 ! local variables
38 integer, parameter :: FFTW_ESTIMATE=64
39 integer p
40 integer(8) plan
41 real(4) t1
42 !$OMP CRITICAL(fftw_)
43 if (sgn == -1) then
44  call sfftw_plan_dft_r2c(plan,nd,n,r,c,fftw_estimate)
45 else
46  call sfftw_plan_dft_c2r(plan,nd,n,c,r,fftw_estimate)
47 end if
48 !$OMP END CRITICAL(fftw_)
49 call sfftw_execute(plan)
50 !$OMP CRITICAL(fftw_)
51 call sfftw_destroy_plan(plan)
52 !$OMP END CRITICAL(fftw_)
53 if (sgn == -1) then
54  p=product(n(1:nd))
55  t1=1.e0/real(p)
56  p=p/n(1)
57  p=p*(n(1)/2+1)
58  call csscal(p,t1,c,1)
59 end if
60 end subroutine
61 
subroutine cfftifc(nd, n, sgn, c)
Definition: cfftifc_fftw.f90:7
subroutine rcfftifc(nd, n, sgn, r, c)