fftRenumber.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 Description
25  Multi-dimensional renumbering used in the Numerical Recipes
26  fft routine. This version is recursive, so works in n-d :
27  determined by the length of array nn
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "fftRenumber.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
36 (
37  List<complex>& data,
38  List<complex>& renumData,
39  const labelList& nn,
40  label nnprod,
41  label ii,
42  label l1,
43  label l2
44 )
45 {
46  if (ii == nn.size())
47  {
48  // we've worked out the renumbering scheme. Now copy
49  // the components across
50 
51  data[l1] = complex(renumData[l2].Re(),renumData[l2].Im());
52  }
53  else
54  {
55  // do another level of folding. First work out the
56  // multiplicative value of the index
57 
58  nnprod /= nn[ii];
59  label i_1(0);
60 
61  for (label i=0; i<nn[ii]; i++)
62  {
63  // now evaluate the indices (both from array 1 and to
64  // array 2). These get multiplied by nnprod to (cumulatively)
65  // find the real position in the list corresponding to
66  // this set of indices.
67 
68  if (i<nn[ii]/2)
69  {
70  i_1 = i + nn[ii]/2;
71  }
72  else
73  {
74  i_1 = i - nn[ii]/2;
75  }
76 
77 
78  // go to the next level of recursion.
79 
81  (
82  data,
83  renumData,
84  nn,
85  nnprod,
86  ii+1,
87  l1+i*nnprod,
88  l2+i_1*nnprod
89  );
90  }
91  }
92 }
93 
94 
96 (
97  List<complex>& data,
98  const labelList& nn
99 )
100 {
101  List<complex> renumData(data);
102 
103  label nnprod(1);
104  forAll(nn, i)
105  {
106  nnprod *= nn[i];
107  }
108 
109  label ii(0), l1(0), l2(0);
110 
112  (
113  data,
114  renumData,
115  nn,
116  nnprod,
117  ii,
118  l1,
119  l2
120  );
121 }
122 
123 
124 // ************************************************************************* //
void fftRenumber(List< complex > &data, const labelList &nn)
Definition: fftRenumber.C:96
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Multi-dimensional renumbering used in the Numerical Recipes fft routine.
scalarField Im(const UList< complex > &cf)
Extension to the c++ complex library type.
Definition: complex.H:76
void fftRenumberRecurse(List< complex > &data, List< complex > &renumData, const labelList &nn, label nnprod, label ii, label l1, label l2)
Definition: fftRenumber.C:36
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97