-
Notifications
You must be signed in to change notification settings - Fork 0
/
fft.F90
272 lines (229 loc) · 9.65 KB
/
fft.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
!C******************************************************************************|
!C fft.f, the FFT package for diablo. VERSION 0.9
!C
!C This file isolates all calls to the FFTW package (available at: www.fftw.org)
!C These wrapper routines were written by T. Bewley (spring 2001).
!C******************************************************************************
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C The arrangement of the significant real numbers in the arrays (denoted by +)
!C in physical space, in Fourier space, and in Fourier space after packing are
!C shown below for the 2D (X-Z) plane. The third direction (Y) is handled in
!C an identical matter as the Z direction shown here.
!C
!C oooooooooooooooooo oooooooooooooooooo oooooooooooooooooo
!C oooooooooooooooooo oooooooooooooooooo oooooooooooooooooo
!C NZ-1 ++++++++++++++++oo -1 ++++++++++++oooooo oooooooooooooooooo
!C ++++++++++++++++oo -2 ++++++++++++oooooo oooooooooooooooooo
!C ++++++++++++++++oo -3 ++++++++++++oooooo oooooooooooooooooo
!C ++++++++++++++++oo ++++++++++++oooooo oooooooooooooooooo
!C ++++++++++++++++oo -NKZ ++++++++++++oooooo oooooooooooooooooo
!C ++++++++++++++++oo oooooooooooooooooo -1 ++++++++++++oooooo
!C ++++++++++++++++oo oooooooooooooooooo -2 ++++++++++++oooooo
!C ++++++++++++++++oo oooooooooooooooooo -3 ++++++++++++oooooo
!C ++++++++++++++++oo oooooooooooooooooo ++++++++++++oooooo
!C ++++++++++++++++oo oooooooooooooooooo -NKZ ++++++++++++oooooo
!C ++++++++++++++++oo NKZ ++++++++++++oooooo NKZ ++++++++++++oooooo
!C ++++++++++++++++oo ++++++++++++oooooo ++++++++++++oooooo
!C 3 ++++++++++++++++oo 3 ++++++++++++oooooo 3 ++++++++++++oooooo
!C 2 ++++++++++++++++oo 2 ++++++++++++oooooo 2 ++++++++++++oooooo
!C 1 ++++++++++++++++oo 1 ++++++++++++oooooo 1 ++++++++++++oooooo
!C 0 ++++++++++++++++oo 0 +o++++++++++oooooo 0 +o++++++++++oooooo
!C ^^^^ ^ ^ ^ ^ ^ ^ ^ ^ ^
!C 0123 NX-1 0 1 2 NKX 0 1 2 NKX
!C
!C PHYSICAL SPACE FOURIER SPACE FOURIER SPACE (PACKED)
!C
!C After the Real->Fourier transform, the significant coefficients are put next
!C to each other in the array, so a loop such as
!C
!C DO K=0,TNKZ [where TNKZ = 2*NKZ = 2*(NZ/3) ]
!C DO I=0,NKX [where NKX = NX/3 ]
!C CP(I,K,J)= ...
!C END DO
!C END DO
!C
!C includes all the Fourier coefficients of interest. The subsequent loops in
!C Fourier space just work on these coefficients in the matrix.
!C
!C Before a Fourier->Real transform, the significant coefficients are unpacked
!C and the higher wavenumbers are SET TO ZERO before the inverse transform.
!C This has the effect of doing the required dealiasing.
!C
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
SUBROUTINE INIT_FFT
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
! INCLUDE 'header_duct'
use ntypes
use Domain
use Grid
use Fft_var
use mpi_var
use TIME_STEP_VAR
implicit none
INTEGER I,J,K
INTEGER FFTW_FORWARD, FFTW_BACKWARD, &
FFTW_ESTIMATE, FFTW_MEASURE, &
FFTW_OUT_OF_PLACE, FFTW_IN_PLACE, &
FFTW_USE_WISDOM, FFTW_THREADSAFE
PARAMETER( FFTW_FORWARD=-1, FFTW_BACKWARD=1, &
FFTW_ESTIMATE=0, FFTW_MEASURE=1, &
FFTW_OUT_OF_PLACE=0, FFTW_IN_PLACE=8, &
FFTW_USE_WISDOM=16, FFTW_THREADSAFE=128 )
if (rank .eq. 0) then
WRITE(6,*) 'Initializing FFTW package.'
endif
PI = 4. * ATAN(1.0)
CI = CMPLX(0.0,1.0)
EPS= 0.000000001
IF (NUM_PER_DIR .GT. 0) THEN
CALL RFFTWND_F77_CREATE_PLAN(FFTW_X_TO_F_PLAN, 1, NX, &
FFTW_FORWARD, FFTW_MEASURE + FFTW_IN_PLACE )
CALL RFFTWND_F77_CREATE_PLAN(FFTW_X_TO_P_PLAN, 1, NX, &
FFTW_BACKWARD, FFTW_MEASURE + FFTW_IN_PLACE )
! NKX=NX/3
RNX=1.0*NX
DO I=0,NKX
KX(I)=I*(2.*PI)/LX
KX2(I)=KX(I)*KX(I)
CIKX(I)=CI*KX(I)
END DO
DO I=0,NX2P
KXP(I)=(I+rank*(NX2P+1))*(2.*PI)/LX
KX2P(I)=KXP(I)*KXP(I)
CIKXP(I)=CI*KXP(I)
END DO
END IF
! DO I=0,NKX
! DO K=0,NZ
! CZX_PLANE(K,I)=CMPLX(0.0,0.0)
! END DO
! END DO
! DO K=0,TNKZ
! DO J=0,NY
! CYZ_PLANE(J,K)=CMPLX(0.0,0.0)
! END DO
! END DO
if (rank .eq. 0) then
WRITE(6,*) 'FFTW package initialized.'
endif
RETURN
END
!C******************************************************************************|
!C-------------> The transform routines for the duct flow follow. <-------------|
!C******************************************************************************|
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
SUBROUTINE FFT_X_TO_FOURIER(U,CU,JMIN,JMAX,KMIN,KMAX)
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C This routine transforms (in 1 direction) planes JMIN-JMAX to Fourier space.
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
! INCLUDE 'header_duct'
use ntypes
use Domain
use Grid
use Fft_var
use TIME_STEP_VAR
implicit none
INTEGER JMIN, JMAX, KMIN, KMAX, I, J, K
REAL*8 U (0:NX+1,0:NZ+1,0:NY+1)
COMPLEX*16 CU(0:NX/2,0:NZ+1,0:NY+1)
!C Looping over the planes of interest, simply perform a real -> complex
!C transform in place in the big storage array, scaling appropriately.
DO J=JMIN,JMAX
CALL RFFTWND_F77_REAL_TO_COMPLEX(FFTW_X_TO_F_PLAN,(KMAX-KMIN+1), &
U(0,KMIN,J), 1, NX+2, CU(0,KMIN,J), 1, NX/2+1)
DO K=KMIN,KMAX
DO I=0,NKX
CU(I,K,J)=CU(I,K,J)/RNX
END DO
END DO
END DO
RETURN
END
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
SUBROUTINE FFT_X_TO_PHYSICAL(CU,U,JMIN,JMAX,KMIN,KMAX)
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C This routine transforms (in 1 direction) planes JMIN-JMAX to physical space.
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
! INCLUDE 'header_duct'
use ntypes
use Domain
use Grid
use Fft_var
use TIME_STEP_VAR
implicit none
INTEGER JMIN, JMAX, KMIN, KMAX, I, J, K
REAL*8 U (0:NX+1,0:NZ+1,0:NY+1)
COMPLEX*16 CU(0:NX/2,0:NZ+1,0:NY+1)
!C Looping over the planes of interest, simply set the higher wavenumbers to
!C zero and then perform a complex -> real transform in place in the big
!C storage array.
DO J=JMIN,JMAX
DO K=KMIN,KMAX
DO I=NKX+1,NX/2
CU(I,K,J)=0.
END DO
END DO
CALL RFFTWND_F77_COMPLEX_TO_REAL(FFTW_X_TO_P_PLAN,(KMAX-KMIN+1), &
CU(0,KMIN,J), 1, NX/2+1, U(0,KMIN,J), 1, NX+2)
END DO
RETURN
END
!C******************************************************************************|
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
SUBROUTINE FFT_X_TO_FOURIER_OP(U,CU,JMIN,JMAX,KMIN,KMAX)
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C This routine transforms (in 1 direction) planes JMIN-JMAX to Fourier space.
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
use ntypes
use Domain
use Grid
use Fft_var
use TIME_STEP_VAR
!use run_variable
implicit none
INTEGER JMIN, JMAX, KMIN, KMAX, I, J, K
REAL*8 U (0:NX+1,0:NZP,0:NY+1)
COMPLEX*16 CU(0:NX/2,0:NZP,0:NY+1)
!C Looping over the planes of interest, simply perform a real -> complex
!C transform in place in the big storage array, scaling appropriately.
DO J=JMIN,JMAX
CALL RFFTWND_F77_REAL_TO_COMPLEX(FFTW_X_TO_F_PLAN,(KMAX-KMIN+1),&
U(0,KMIN,J), 1, NX+2, CU(0,KMIN,J), 1, NX/2+1)
DO K=KMIN,KMAX
DO I=0,NKX
CU(I,K,J)=CU(I,K,J)/RNX
END DO
END DO
END DO
RETURN
END
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
SUBROUTINE FFT_X_TO_PHYSICAL_OP(CU,U,JMIN,JMAX,KMIN,KMAX)
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
!C This routine transforms (in 1 direction) planes JMIN-JMAX to physical space.
!C----*|--.---------.---------.---------.---------.---------.---------.-|-------|
use ntypes
use Domain
use Grid
use Fft_var
use TIME_STEP_VAR
!use run_variable
implicit none
INTEGER JMIN, JMAX, KMIN, KMAX, I, J, K
REAL*8 U (0:NX+1,0:NZP,0:NY+1)
COMPLEX*16 CU(0:NX/2,0:NZP,0:NY+1)
!C Looping over the planes of interest, simply set the higher wavenumbers to
!C zero and then perform a complex -> real transform in place in the big
!C storage array.
DO J=JMIN,JMAX
DO K=KMIN,KMAX
DO I=NKX+1,NX/2
CU(I,K,J)=0.
END DO
END DO
CALL RFFTWND_F77_COMPLEX_TO_REAL(FFTW_X_TO_P_PLAN,(KMAX-KMIN+1),&
CU(0,KMIN,J), 1, NX/2+1, U(0,KMIN,J), 1, NX+2)
END DO
RETURN
END