Simplified the pow-of-2 check. See
[OpenFOAM-1.5.x.git] / src / randomProcesses / fft / fft.C
blob9509e97ab129a00c37a0f920d86039e9bb8e9502
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "fft.H"
28 #include "fftRenumber.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 #define SWAP(a,b) tempr=(a); (a)=(b); (b)=tempr
38 #define TWOPI 6.28318530717959
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 void fft::transform
44     complexField& field,
45     const labelList& nn,
46     transformDirection isign
49     forAll(nn, idim)
50     {
51         scalar pow2 = log(scalar(nn[idim]))/log(scalar(2));
52         if ((pow2 - int(pow2 + 0.5)) > SMALL)
53         {
54             FatalErrorIn
55             (
56                  "fft::transform(complexField&, const labelList&, "
57                  "transformDirection)"
58             )   << "number of elements in direction " << idim
59                 << " is not a power of 2" << endl
60                 << "    Number of elements in each direction = " << nn
61                 << abort(FatalError);
62         }
63     }
65     const label ndim = nn.size();
67     label i1, i2, i3, i2rev, i3rev, ip1, ip2, ip3, ifp1, ifp2;
68     label ibit, k1, k2, n, nprev, nrem, idim;
69     scalar tempi, tempr;
70     scalar theta, wi, wpi, wpr, wr, wtemp;
71     scalar* data = reinterpret_cast<scalar*>(field.begin()) - 1;
74     // if inverse transform : renumber before transform
76     if (isign == REVERSE_TRANSFORM)
77     {
78         fftRenumber(field, nn);
79     }
82     label ntot = 1;
83     forAll(nn, idim)
84     {
85         ntot *= nn[idim];
86     }
89     nprev = 1;
91     for (idim=ndim; idim>=1; idim--)
92     {
93         n = nn[idim-1];
94         nrem = ntot/(n*nprev);
95         ip1 = nprev << 1;
96         ip2 = ip1*n;
97         ip3 = ip2*nrem;
98         i2rev = 1;
100         for (i2=1; i2<=ip2; i2+=ip1)
101         {
102             if (i2 < i2rev)
103             {
104                 for (i1=i2; i1<=i2 + ip1 - 2; i1+=2)
105                 {
106                     for (i3=i1; i3<=ip3; i3+=ip2)
107                     {
108                         i3rev = i2rev + i3 - i2;
109                         SWAP(data[i3], data[i3rev]);
110                         SWAP(data[i3 + 1], data[i3rev + 1]);
111                     }
112                 }
113             }
115             ibit = ip2 >> 1;
116             while (ibit >= ip1 && i2rev > ibit)
117             {
118                 i2rev -= ibit;
119                 ibit >>= 1;
120             }
122             i2rev += ibit;
123         }
125         ifp1 = ip1;
127         while (ifp1 < ip2)
128         {
129             ifp2 = ifp1 << 1;
130             theta = isign*TWOPI/(ifp2/ip1);
131             wtemp = sin(0.5*theta);
132             wpr = -2.0*wtemp*wtemp;
133             wpi = sin(theta);
134             wr = 1.0;
135             wi = 0.0;
137             for (i3 = 1; i3 <= ifp1; i3 += ip1)
138             {
139                 for (i1 = i3; i1 <= i3 + ip1 - 2; i1 += 2)
140                 {
141                     for (i2 = i1; i2 <= ip3; i2 += ifp2)
142                     {
143                         k1 = i2;
144                         k2 = k1 + ifp1;
145                         tempr = scalar(wr*data[k2]) - scalar(wi*data[k2 + 1]);
146                         tempi = scalar(wr*data[k2 + 1]) + scalar(wi*data[k2]);
147                         data[k2] = data[k1] - tempr;
148                         data[k2 + 1] = data[k1 + 1] - tempi;
149                         data[k1] += tempr;
150                         data[k1 + 1] += tempi;
151                     }
152                 }
154                 wr = (wtemp = wr)*wpr - wi*wpi + wr;
155                 wi = wi*wpr + wtemp*wpi + wi;
156             }
157             ifp1 = ifp2;
158         }
159         nprev *= n;
160     }
163     // if forward transform : renumber after transform
165     if (isign == FORWARD_TRANSFORM)
166     {
167         fftRenumber(field, nn);
168     }
171     // scale result (symmetric scale both forward and inverse transform)
173     scalar recRootN = 1.0/sqrt(scalar(ntot));
175     forAll(field, i)
176     {
177         field[i] *= recRootN;
178     }
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 #undef SWAP
185 #undef TWOPI
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 tmp<complexField> fft::forwardTransform
191     const tmp<complexField>& tfield,
192     const labelList& nn
195     tmp<complexField> tfftField(new complexField(tfield));
197     transform(tfftField(), nn, FORWARD_TRANSFORM);
199     tfield.clear();
201     return tfftField;
205 tmp<complexField> fft::reverseTransform
207     const tmp<complexField>& tfield,
208     const labelList& nn
211     tmp<complexField> tifftField(new complexField(tfield));
213     transform(tifftField(), nn, REVERSE_TRANSFORM);
215     tfield.clear();
217     return tifftField;
221 tmp<complexVectorField> fft::forwardTransform
223     const tmp<complexVectorField>& tfield,
224     const labelList& nn
227     tmp<complexVectorField> tfftVectorField
228     (
229         new complexVectorField
230         (
231             tfield().size()
232         )
233     );
235     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
236     {
237         tfftVectorField().replace
238         (
239             cmpt,
240             forwardTransform(tfield().component(cmpt), nn)
241         );
242     }
244     tfield.clear();
246     return tfftVectorField;
250 tmp<complexVectorField> fft::reverseTransform
252     const tmp<complexVectorField>& tfield,
253     const labelList& nn
256     tmp<complexVectorField> tifftVectorField
257     (
258         new complexVectorField
259         (
260             tfield().size()
261         )
262     );
264     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
265     {
266         tifftVectorField().replace
267         (
268             cmpt,
269             reverseTransform(tfield().component(cmpt), nn)
270         );
271     }
273     tfield.clear();
275     return tifftVectorField;
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 } // End namespace Foam
283 // ************************************************************************* //