FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

fft.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include <randomProcesses/fft.H>
00027 #include <randomProcesses/fftRenumber.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 
00034 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00035 
00036 #define SWAP(a,b) tempr=(a); (a)=(b); (b)=tempr
00037 #define TWOPI 6.28318530717959
00038 
00039 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00040 
00041 void fft::transform
00042 (
00043     complexField& field,
00044     const labelList& nn,
00045     transformDirection isign
00046 )
00047 {
00048     forAll(nn, idim)
00049     {
00050         // Check for power of two
00051         unsigned int dimCount = nn[idim];
00052         if (!dimCount || (dimCount & (dimCount - 1)))
00053         {
00054             FatalErrorIn
00055             (
00056                  "fft::transform(complexField&, const labelList&, "
00057                  "transformDirection)"
00058             )   << "number of elements in direction " << idim
00059                 << " is not a power of 2" << endl
00060                 << "    Number of elements in each direction = " << nn
00061                 << abort(FatalError);
00062         }
00063     }
00064 
00065     const label ndim = nn.size();
00066 
00067     label i1, i2, i3, i2rev, i3rev, ip1, ip2, ip3, ifp1, ifp2;
00068     label ibit, k1, k2, n, nprev, nrem, idim;
00069     scalar tempi, tempr;
00070     scalar theta, wi, wpi, wpr, wr, wtemp;
00071     scalar* data = reinterpret_cast<scalar*>(field.begin()) - 1;
00072 
00073 
00074     // if inverse transform : renumber before transform
00075 
00076     if (isign == REVERSE_TRANSFORM)
00077     {
00078         fftRenumber(field, nn);
00079     }
00080 
00081 
00082     label ntot = 1;
00083     forAll(nn, idim)
00084     {
00085         ntot *= nn[idim];
00086     }
00087 
00088 
00089     nprev = 1;
00090 
00091     for (idim=ndim; idim>=1; idim--)
00092     {
00093         n = nn[idim-1];
00094         nrem = ntot/(n*nprev);
00095         ip1 = nprev << 1;
00096         ip2 = ip1*n;
00097         ip3 = ip2*nrem;
00098         i2rev = 1;
00099 
00100         for (i2=1; i2<=ip2; i2+=ip1)
00101         {
00102             if (i2 < i2rev)
00103             {
00104                 for (i1=i2; i1<=i2 + ip1 - 2; i1+=2)
00105                 {
00106                     for (i3=i1; i3<=ip3; i3+=ip2)
00107                     {
00108                         i3rev = i2rev + i3 - i2;
00109                         SWAP(data[i3], data[i3rev]);
00110                         SWAP(data[i3 + 1], data[i3rev + 1]);
00111                     }
00112                 }
00113             }
00114 
00115             ibit = ip2 >> 1;
00116             while (ibit >= ip1 && i2rev > ibit)
00117             {
00118                 i2rev -= ibit;
00119                 ibit >>= 1;
00120             }
00121 
00122             i2rev += ibit;
00123         }
00124 
00125         ifp1 = ip1;
00126 
00127         while (ifp1 < ip2)
00128         {
00129             ifp2 = ifp1 << 1;
00130             theta = isign*TWOPI/(ifp2/ip1);
00131             wtemp = sin(0.5*theta);
00132             wpr = -2.0*wtemp*wtemp;
00133             wpi = sin(theta);
00134             wr = 1.0;
00135             wi = 0.0;
00136 
00137             for (i3 = 1; i3 <= ifp1; i3 += ip1)
00138             {
00139                 for (i1 = i3; i1 <= i3 + ip1 - 2; i1 += 2)
00140                 {
00141                     for (i2 = i1; i2 <= ip3; i2 += ifp2)
00142                     {
00143                         k1 = i2;
00144                         k2 = k1 + ifp1;
00145                         tempr = scalar(wr*data[k2]) - scalar(wi*data[k2 + 1]);
00146                         tempi = scalar(wr*data[k2 + 1]) + scalar(wi*data[k2]);
00147                         data[k2] = data[k1] - tempr;
00148                         data[k2 + 1] = data[k1 + 1] - tempi;
00149                         data[k1] += tempr;
00150                         data[k1 + 1] += tempi;
00151                     }
00152                 }
00153 
00154                 wr = (wtemp = wr)*wpr - wi*wpi + wr;
00155                 wi = wi*wpr + wtemp*wpi + wi;
00156             }
00157             ifp1 = ifp2;
00158         }
00159         nprev *= n;
00160     }
00161 
00162 
00163     // if forward transform : renumber after transform
00164 
00165     if (isign == FORWARD_TRANSFORM)
00166     {
00167         fftRenumber(field, nn);
00168     }
00169 
00170 
00171     // scale result (symmetric scale both forward and inverse transform)
00172 
00173     scalar recRootN = 1.0/sqrt(scalar(ntot));
00174 
00175     forAll(field, i)
00176     {
00177         field[i] *= recRootN;
00178     }
00179 }
00180 
00181 
00182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00183 
00184 #undef SWAP
00185 #undef TWOPI
00186 
00187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00188 
00189 tmp<complexField> fft::forwardTransform
00190 (
00191     const tmp<complexField>& tfield,
00192     const labelList& nn
00193 )
00194 {
00195     tmp<complexField> tfftField(new complexField(tfield));
00196 
00197     transform(tfftField(), nn, FORWARD_TRANSFORM);
00198 
00199     tfield.clear();
00200 
00201     return tfftField;
00202 }
00203 
00204 
00205 tmp<complexField> fft::reverseTransform
00206 (
00207     const tmp<complexField>& tfield,
00208     const labelList& nn
00209 )
00210 {
00211     tmp<complexField> tifftField(new complexField(tfield));
00212 
00213     transform(tifftField(), nn, REVERSE_TRANSFORM);
00214 
00215     tfield.clear();
00216 
00217     return tifftField;
00218 }
00219 
00220 
00221 tmp<complexVectorField> fft::forwardTransform
00222 (
00223     const tmp<complexVectorField>& tfield,
00224     const labelList& nn
00225 )
00226 {
00227     tmp<complexVectorField> tfftVectorField
00228     (
00229         new complexVectorField
00230         (
00231             tfield().size()
00232         )
00233     );
00234 
00235     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
00236     {
00237         tfftVectorField().replace
00238         (
00239             cmpt,
00240             forwardTransform(tfield().component(cmpt), nn)
00241         );
00242     }
00243 
00244     tfield.clear();
00245 
00246     return tfftVectorField;
00247 }
00248 
00249 
00250 tmp<complexVectorField> fft::reverseTransform
00251 (
00252     const tmp<complexVectorField>& tfield,
00253     const labelList& nn
00254 )
00255 {
00256     tmp<complexVectorField> tifftVectorField
00257     (
00258         new complexVectorField
00259         (
00260             tfield().size()
00261         )
00262     );
00263 
00264     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
00265     {
00266         tifftVectorField().replace
00267         (
00268             cmpt,
00269             reverseTransform(tfield().component(cmpt), nn)
00270         );
00271     }
00272 
00273     tfield.clear();
00274 
00275     return tifftVectorField;
00276 }
00277 
00278 
00279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00280 
00281 } // End namespace Foam
00282 
00283 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines