fft.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fft.H"
27 #include "fftRenumber.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 #define SWAP(a,b) tempr=(a); (a)=(b); (b)=tempr
37 #define TWOPI 6.28318530717959
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 (
43  complexField& field,
44  const labelList& nn,
45  transformDirection isign
46 )
47 {
48  forAll(nn, idim)
49  {
50  // Check for power of two
51  unsigned int dimCount = nn[idim];
52  if (!dimCount || (dimCount & (dimCount - 1)))
53  {
55  << "number of elements in direction " << idim
56  << " is not a power of 2" << endl
57  << " Number of elements in each direction = " << nn
58  << abort(FatalError);
59  }
60  }
61 
62  const label ndim = nn.size();
63 
64  label i1, i2, i3, i2rev, i3rev, ip1, ip2, ip3, ifp1, ifp2;
65  label ibit, k1, k2, n, nprev, nrem, idim;
66  scalar tempi, tempr;
67  scalar theta, wi, wpi, wpr, wr, wtemp;
68  scalar* data = reinterpret_cast<scalar*>(field.begin()) - 1;
69 
70 
71  // if inverse transform : renumber before transform
72 
73  if (isign == REVERSE_TRANSFORM)
74  {
75  fftRenumber(field, nn);
76  }
77 
78 
79  label ntot = 1;
80  forAll(nn, idim)
81  {
82  ntot *= nn[idim];
83  }
84 
85 
86  nprev = 1;
87 
88  for (idim=ndim; idim>=1; idim--)
89  {
90  n = nn[idim-1];
91  nrem = ntot/(n*nprev);
92  ip1 = nprev << 1;
93  ip2 = ip1*n;
94  ip3 = ip2*nrem;
95  i2rev = 1;
96 
97  for (i2=1; i2<=ip2; i2+=ip1)
98  {
99  if (i2 < i2rev)
100  {
101  for (i1=i2; i1<=i2 + ip1 - 2; i1+=2)
102  {
103  for (i3=i1; i3<=ip3; i3+=ip2)
104  {
105  i3rev = i2rev + i3 - i2;
106  SWAP(data[i3], data[i3rev]);
107  SWAP(data[i3 + 1], data[i3rev + 1]);
108  }
109  }
110  }
111 
112  ibit = ip2 >> 1;
113  while (ibit >= ip1 && i2rev > ibit)
114  {
115  i2rev -= ibit;
116  ibit >>= 1;
117  }
118 
119  i2rev += ibit;
120  }
121 
122  ifp1 = ip1;
123 
124  while (ifp1 < ip2)
125  {
126  ifp2 = ifp1 << 1;
127  theta = isign*TWOPI/(ifp2/ip1);
128  wtemp = sin(0.5*theta);
129  wpr = -2.0*wtemp*wtemp;
130  wpi = sin(theta);
131  wr = 1.0;
132  wi = 0.0;
133 
134  for (i3 = 1; i3 <= ifp1; i3 += ip1)
135  {
136  for (i1 = i3; i1 <= i3 + ip1 - 2; i1 += 2)
137  {
138  for (i2 = i1; i2 <= ip3; i2 += ifp2)
139  {
140  k1 = i2;
141  k2 = k1 + ifp1;
142  tempr = scalar(wr*data[k2]) - scalar(wi*data[k2 + 1]);
143  tempi = scalar(wr*data[k2 + 1]) + scalar(wi*data[k2]);
144  data[k2] = data[k1] - tempr;
145  data[k2 + 1] = data[k1 + 1] - tempi;
146  data[k1] += tempr;
147  data[k1 + 1] += tempi;
148  }
149  }
150 
151  wr = (wtemp = wr)*wpr - wi*wpi + wr;
152  wi = wi*wpr + wtemp*wpi + wi;
153  }
154  ifp1 = ifp2;
155  }
156  nprev *= n;
157  }
158 
159 
160  // if forward transform : renumber after transform
161 
162  if (isign == FORWARD_TRANSFORM)
163  {
164  fftRenumber(field, nn);
165  }
166 
167 
168  // scale result (symmetric scale both forward and inverse transform)
169 
170  scalar recRootN = 1.0/sqrt(scalar(ntot));
171 
172  forAll(field, i)
173  {
174  field[i] *= recRootN;
175  }
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #undef SWAP
182 #undef TWOPI
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
187 (
188  const tmp<complexField>& tfield,
189  const labelList& nn
190 )
191 {
192  tmp<complexField> tfftField(new complexField(tfield));
193 
194  transform(tfftField.ref(), nn, FORWARD_TRANSFORM);
195 
196  tfield.clear();
197 
198  return tfftField;
199 }
200 
201 
203 (
204  const tmp<complexField>& tfield,
205  const labelList& nn
206 )
207 {
208  tmp<complexField> tifftField(new complexField(tfield));
209 
210  transform(tifftField.ref(), nn, REVERSE_TRANSFORM);
211 
212  tfield.clear();
213 
214  return tifftField;
215 }
216 
217 
219 (
220  const tmp<complexVectorField>& tfield,
221  const labelList& nn
222 )
223 {
224  tmp<complexVectorField> tfftVectorField
225  (
227  (
228  tfield().size()
229  )
230  );
231 
232  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
233  {
234  tfftVectorField.ref().replace
235  (
236  cmpt,
237  forwardTransform(tfield().component(cmpt), nn)
238  );
239  }
240 
241  tfield.clear();
242 
243  return tfftVectorField;
244 }
245 
246 
248 (
249  const tmp<complexVectorField>& tfield,
250  const labelList& nn
251 )
252 {
253  tmp<complexVectorField> tifftVectorField
254  (
256  (
257  tfield().size()
258  )
259  );
260 
261  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
262  {
263  tifftVectorField.ref().replace
264  (
265  cmpt,
266  reverseTransform(tfield().component(cmpt), nn)
267  );
268  }
269 
270  tfield.clear();
271 
272  return tifftVectorField;
273 }
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Pre-declare SubField and related Field type.
Definition: Field.H:83
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
iterator begin()
Return an iterator to begin traversing the UList.
Definition: UListI.H:216
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:105
static tmp< complexField > reverseTransform(const tmp< complexField > &field, const labelList &nn)
Definition: fft.C:203
static void transform(complexField &field, const labelList &nn, transformDirection fftDirection)
Definition: fft.C:42
transformDirection
Definition: fft.H:57
@ REVERSE_TRANSFORM
Definition: fft.H:59
@ FORWARD_TRANSFORM
Definition: fft.H:58
static tmp< complexField > forwardTransform(const tmp< complexField > &field, const labelList &nn)
Definition: fft.C:187
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Multi-dimensional renumbering used in the Numerical Recipes fft routine.
#define TWOPI
Definition: fft.C:37
#define SWAP(a, b)
Definition: fft.C:36
Namespace for OpenFOAM.
void fftRenumber(List< complex > &data, const labelList &nn)
Definition: fftRenumber.C:96
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
Field< complex > complexField
Specialisation of Field<T> for complex.
Definition: complexFields.H:55
dimensionedScalar sqrt(const dimensionedScalar &ds)
error FatalError
uint8_t direction
Definition: direction.H:45