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_common_3d.f90
253 lines (208 loc) · 6.7 KB
/
fft_common_3d.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
!=======================================================================
! 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 file contains 3D c2c/r2c/c2r transform subroutines which are
! identical for several FFT engines
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! 3D FFT - complex to complex
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine fft_3d_c2c(in, out, isign)
implicit none
complex(mytype), dimension(:,:,:), intent(INOUT) :: in
complex(mytype), dimension(:,:,:), intent(OUT) :: out
integer, intent(IN) :: isign
#ifndef OVERWRITE
complex(mytype), allocatable, dimension(:,:,:) :: wk1
#endif
if (format==PHYSICAL_IN_X .AND. isign==DECOMP_2D_FFT_FORWARD .OR. &
format==PHYSICAL_IN_Z .AND. isign==DECOMP_2D_FFT_BACKWARD) then
! ===== 1D FFTs in X =====
#ifdef OVERWRITE
call c2c_1m_x(in,isign,ph)
#else
allocate (wk1(ph%xsz(1),ph%xsz(2),ph%xsz(3)))
wk1 = in
call c2c_1m_x(wk1,isign,ph)
#endif
! ===== Swap X --> Y; 1D FFTs in Y =====
if (dims(1)>1) then
#ifdef OVERWRITE
call transpose_x_to_y(in,wk2_c2c,ph)
#else
call transpose_x_to_y(wk1,wk2_c2c,ph)
#endif
call c2c_1m_y(wk2_c2c,isign,ph)
else
#ifdef OVERWRITE
call c2c_1m_y(in,isign,ph)
#else
call c2c_1m_y(wk1,isign,ph)
#endif
end if
! ===== Swap Y --> Z; 1D FFTs in Z =====
if (dims(1)>1) then
call transpose_y_to_z(wk2_c2c,out,ph)
else
#ifdef OVERWRITE
call transpose_y_to_z(in,out,ph)
#else
call transpose_y_to_z(wk1,out,ph)
#endif
end if
call c2c_1m_z(out,isign,ph)
else if (format==PHYSICAL_IN_X .AND. isign==DECOMP_2D_FFT_BACKWARD &
.OR. &
format==PHYSICAL_IN_Z .AND. isign==DECOMP_2D_FFT_FORWARD) then
! ===== 1D FFTs in Z =====
#ifdef OVERWRITE
call c2c_1m_z(in,isign,ph)
#else
allocate (wk1(ph%zsz(1),ph%zsz(2),ph%zsz(3)))
wk1 = in
call c2c_1m_z(wk1,isign,ph)
#endif
! ===== Swap Z --> Y; 1D FFTs in Y =====
if (dims(1)>1) then
#ifdef OVERWRITE
call transpose_z_to_y(in,wk2_c2c,ph)
#else
call transpose_z_to_y(wk1,wk2_c2c,ph)
#endif
call c2c_1m_y(wk2_c2c,isign,ph)
else ! out==wk2_c2c if 1D decomposition
#ifdef OVERWRITE
call transpose_z_to_y(in,out,ph)
#else
call transpose_z_to_y(wk1,out,ph)
#endif
call c2c_1m_y(out,isign,ph)
end if
! ===== Swap Y --> X; 1D FFTs in X =====
if (dims(1)>1) then
call transpose_y_to_x(wk2_c2c,out,ph)
end if
call c2c_1m_x(out,isign,ph)
end if
return
end subroutine fft_3d_c2c
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! 3D forward FFT - real to complex
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine fft_3d_r2c(in_r, out_c)
implicit none
real(mytype), dimension(:,:,:), intent(IN) :: in_r
complex(mytype), dimension(:,:,:), intent(OUT) :: out_c
if (format==PHYSICAL_IN_X) then
! ===== 1D FFTs in X =====
call r2c_1m_x(in_r,wk13)
! ===== Swap X --> Y; 1D FFTs in Y =====
if (dims(1)>1) then
call transpose_x_to_y(wk13,wk2_r2c,sp)
call c2c_1m_y(wk2_r2c,-1,sp)
else
call c2c_1m_y(wk13,-1,sp)
end if
! ===== Swap Y --> Z; 1D FFTs in Z =====
if (dims(1)>1) then
call transpose_y_to_z(wk2_r2c,out_c,sp)
else
call transpose_y_to_z(wk13,out_c,sp)
end if
call c2c_1m_z(out_c,-1,sp)
else if (format==PHYSICAL_IN_Z) then
! ===== 1D FFTs in Z =====
call r2c_1m_z(in_r,wk13)
! ===== Swap Z --> Y; 1D FFTs in Y =====
if (dims(1)>1) then
call transpose_z_to_y(wk13,wk2_r2c,sp)
call c2c_1m_y(wk2_r2c,-1,sp)
else ! out_c==wk2_r2c if 1D decomposition
call transpose_z_to_y(wk13,out_c,sp)
call c2c_1m_y(out_c,-1,sp)
end if
! ===== Swap Y --> X; 1D FFTs in X =====
if (dims(1)>1) then
call transpose_y_to_x(wk2_r2c,out_c,sp)
end if
call c2c_1m_x(out_c,-1,sp)
end if
return
end subroutine fft_3d_r2c
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! 3D inverse FFT - complex to real
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine fft_3d_c2r(in_c, out_r)
implicit none
complex(mytype), dimension(:,:,:), intent(INOUT) :: in_c
real(mytype), dimension(:,:,:), intent(OUT) :: out_r
#ifndef OVERWRITE
complex(mytype), allocatable, dimension(:,:,:) :: wk1
#endif
if (format==PHYSICAL_IN_X) then
! ===== 1D FFTs in Z =====
#ifdef OVERWRITE
call c2c_1m_z(in_c,1,sp)
#else
allocate(wk1(sp%zsz(1),sp%zsz(2),sp%zsz(3)))
wk1 = in_c
call c2c_1m_z(wk1,1,sp)
#endif
! ===== Swap Z --> Y; 1D FFTs in Y =====
#ifdef OVERWRITE
call transpose_z_to_y(in_c,wk2_r2c,sp)
#else
call transpose_z_to_y(wk1,wk2_r2c,sp)
#endif
call c2c_1m_y(wk2_r2c,1,sp)
! ===== Swap Y --> X; 1D FFTs in X =====
if (dims(1)>1) then
call transpose_y_to_x(wk2_r2c,wk13,sp)
call c2r_1m_x(wk13,out_r)
else
call c2r_1m_x(wk2_r2c,out_r)
end if
else if (format==PHYSICAL_IN_Z) then
! ===== 1D FFTs in X =====
#ifdef OVERWRITE
call c2c_1m_x(in_c,1,sp)
#else
allocate(wk1(sp%xsz(1),sp%xsz(2),sp%xsz(3)))
wk1 = in_c
call c2c_1m_x(wk1,1,sp)
#endif
! ===== Swap X --> Y; 1D FFTs in Y =====
if (dims(1)>1) then
#ifdef OVERWRITE
call transpose_x_to_y(in_c,wk2_r2c,sp)
#else
call transpose_x_to_y(wk1,wk2_r2c,sp)
#endif
call c2c_1m_y(wk2_r2c,1,sp)
else ! in_c==wk2_r2c if 1D decomposition
#ifdef OVERWRITE
call c2c_1m_y(in_c,1,sp)
#else
call c2c_1m_y(wk1,1,sp)
#endif
end if
! ===== Swap Y --> Z; 1D FFTs in Z =====
if (dims(1)>1) then
call transpose_y_to_z(wk2_r2c,wk13,sp)
else
#ifdef OVERWRITE
call transpose_y_to_z(in_c,wk13,sp)
#else
call transpose_y_to_z(wk1,wk13,sp)
#endif
end if
call c2r_1m_z(wk13,out_r)
end if
return
end subroutine fft_3d_c2r