This repository has been archived by the owner on Nov 7, 2019. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 2
/
fft_generic.f90
303 lines (226 loc) · 7.52 KB
/
fft_generic.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
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
!=======================================================================
! This is part of the 2DECOMP&FFT library
!
! 2DECOMP&FFT is a software framework for general-purpose 2D (pencil)
! decomposition. It also implements a highly scalable distributed
! three-dimensional Fast Fourier Transform (FFT).
!
! Copyright (C) 2009-2011 Ning Li, the Numerical Algorithms Group (NAG)
!
!=======================================================================
! This is the 'generic' implementation of the FFT library
module decomp_2d_fft
use decomp_2d ! 2D decomposition module
use glassman
implicit none
private ! Make everything private unless declared public
! engine-specific global variables
complex(mytype), allocatable, dimension(:) :: buf, scratch
! common code used for all engines, including global variables,
! generic interface definitions and several subroutines
#include "fft_common.f90"
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! This routine performs one-time initialisations for the FFT engine
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine init_fft_engine
implicit none
integer :: cbuf_size
if (nrank==0) then
write(*,*) ' '
write(*,*) '***** Using the generic FFT engine *****'
write(*,*) ' '
end if
cbuf_size = max(ph%xsz(1), ph%ysz(2))
cbuf_size = max(cbuf_size, ph%zsz(3))
allocate(buf(cbuf_size))
allocate(scratch(cbuf_size))
return
end subroutine init_fft_engine
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! This routine performs one-time finalisations for the FFT engine
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine finalize_fft_engine
implicit none
deallocate(buf,scratch)
return
end subroutine finalize_fft_engine
! Following routines calculate multiple one-dimensional FFTs to form
! the basis of three-dimensional FFTs.
! c2c transform, multiple 1D FFTs in x direction
subroutine c2c_1m_x(inout, isign, decomp)
implicit none
complex(mytype), dimension(:,:,:), intent(INOUT) :: inout
integer, intent(IN) :: isign
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k
do k=1,decomp%xsz(3)
do j=1,decomp%xsz(2)
do i=1,decomp%xsz(1)
buf(i) = inout(i,j,k)
end do
call spcfft(buf,decomp%xsz(1),isign,scratch)
do i=1,decomp%xsz(1)
inout(i,j,k) = buf(i)
end do
end do
end do
return
end subroutine c2c_1m_x
! c2c transform, multiple 1D FFTs in y direction
subroutine c2c_1m_y(inout, isign, decomp)
implicit none
complex(mytype), dimension(:,:,:), intent(INOUT) :: inout
integer, intent(IN) :: isign
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k
do k=1,decomp%ysz(3)
do i=1,decomp%ysz(1)
do j=1,decomp%ysz(2)
buf(j) = inout(i,j,k)
end do
call spcfft(buf,decomp%ysz(2),isign,scratch)
do j=1,decomp%ysz(2)
inout(i,j,k) = buf(j)
end do
end do
end do
return
end subroutine c2c_1m_y
! c2c transform, multiple 1D FFTs in z direction
subroutine c2c_1m_z(inout, isign, decomp)
implicit none
complex(mytype), dimension(:,:,:), intent(INOUT) :: inout
integer, intent(IN) :: isign
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k
do j=1,decomp%zsz(2)
do i=1,decomp%zsz(1)
do k=1,decomp%zsz(3)
buf(k) = inout(i,j,k)
end do
call spcfft(buf,decomp%zsz(3),isign,scratch)
do k=1,decomp%zsz(3)
inout(i,j,k) = buf(k)
end do
end do
end do
return
end subroutine c2c_1m_z
! r2c transform, multiple 1D FFTs in x direction
subroutine r2c_1m_x(input, output)
implicit none
real(mytype), dimension(:,:,:), intent(IN) :: input
complex(mytype), dimension(:,:,:), intent(OUT) :: output
integer :: i,j,k, s1,s2,s3, d1
s1 = size(input,1)
s2 = size(input,2)
s3 = size(input,3)
d1 = size(output,1)
do k=1,s3
do j=1,s2
! Glassman's FFT is c2c only,
! needing some pre- and post-processing for r2c
! pack real input in complex storage
do i=1,s1
buf(i) = cmplx(input(i,j,k),0._mytype, kind=mytype)
end do
call spcfft(buf,s1,-1,scratch)
! note d1 ~ s1/2+1
! simply drop the redundant part of the complex output
do i=1,d1
output(i,j,k) = buf(i)
end do
end do
end do
return
end subroutine r2c_1m_x
! r2c transform, multiple 1D FFTs in z direction
subroutine r2c_1m_z(input, output)
implicit none
real(mytype), dimension(:,:,:), intent(IN) :: input
complex(mytype), dimension(:,:,:), intent(OUT) :: output
integer :: i,j,k, s1,s2,s3, d3
s1 = size(input,1)
s2 = size(input,2)
s3 = size(input,3)
d3 = size(output,3)
do j=1,s2
do i=1,s1
! Glassman's FFT is c2c only,
! needing some pre- and post-processing for r2c
! pack real input in complex storage
do k=1,s3
buf(k) = cmplx(input(i,j,k),0._mytype, kind=mytype)
end do
call spcfft(buf,s3,-1,scratch)
! note d3 ~ s3/2+1
! simply drop the redundant part of the complex output
do k=1,d3
output(i,j,k) = buf(k)
end do
end do
end do
return
end subroutine r2c_1m_z
! c2r transform, multiple 1D FFTs in x direction
subroutine c2r_1m_x(input, output)
implicit none
complex(mytype), dimension(:,:,:), intent(IN) :: input
real(mytype), dimension(:,:,:), intent(OUT) :: output
integer :: i,j,k, d1,d2,d3
d1 = size(output,1)
d2 = size(output,2)
d3 = size(output,3)
do k=1,d3
do j=1,d2
! Glassman's FFT is c2c only,
! needing some pre- and post-processing for c2r
do i=1,d1/2+1
buf(i) = input(i,j,k)
end do
! expanding to a full-size complex array
! For odd N, the storage is:
! 1, 2, ...... N/2+1 integer division rounded down
! N, ...... N/2+2 => a(i) is conjugate of a(N+2-i)
! For even N, the storage is:
! 1, 2, ...... N/2 , N/2+1
! N, ...... N/2+2 again a(i) conjugate of a(N+2-i)
do i=d1/2+2,d1
buf(i) = conjg(buf(d1+2-i))
end do
call spcfft(buf,d1,1,scratch)
do i=1,d1
! simply drop imaginary part
output(i,j,k) = real(buf(i), kind=mytype)
end do
end do
end do
return
end subroutine c2r_1m_x
! c2r transform, multiple 1D FFTs in z direction
subroutine c2r_1m_z(input, output)
implicit none
complex(mytype), dimension(:,:,:), intent(IN) :: input
real(mytype), dimension(:,:,:), intent(OUT) :: output
integer :: i,j,k, d1,d2,d3
d1 = size(output,1)
d2 = size(output,2)
d3 = size(output,3)
do j=1,d2
do i=1,d1
do k=1,d3/2+1
buf(k) = input(i,j,k)
end do
do k=d3/2+2,d3
buf(k) = conjg(buf(d3+2-k))
end do
call spcfft(buf,d3,1,scratch)
do k=1,d3
output(i,j,k) = real(buf(k), kind=mytype)
end do
end do
end do
return
end subroutine c2r_1m_z
#include "fft_common_3d.f90"
end module decomp_2d_fft