OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
fftpack_FFT.h
Go to the documentation of this file.
1//
2// IPPL FFT
3//
4// Copyright (c) 2008-2018
5// Paul Scherrer Institut, Villigen PSI, Switzerland
6// All rights reserved.
7//
8// OPAL is licensed under GNU GPL version 3.
9//
10
11/*
12 Prototypes for accessing Fortran 1D FFT routines from
13 Netlib, and definitions for templated class FFTPACK, which acts as an
14 FFT engine for the FFT class, providing storage for trigonometric
15 information and performing the 1D FFTs as needed.
16*/
17
18#ifndef IPPL_FFT_FFTPACK_FFT_H
19#define IPPL_FFT_FFTPACK_FFT_H
20
21#include "Utility/PAssert.h"
22#include "Utility/IpplInfo.h"
23
24// FFTPACK function prototypes for Fortran routines
25extern "C" {
26 // double-precision CC FFT
27 void cffti (size_t n, double& wsave);
28 void cfftf (size_t n, double& r, double& wsave);
29 void cfftb (size_t n, double& r, double& wsave);
30 // double-precision RC FFT
31 void rffti (size_t n, double& wsave);
32 void rfftf (size_t n, double& r, double& wsave);
33 void rfftb (size_t n, double& r, double& wsave);
34 // double-precision sine transform
35 void sinti (size_t n, double& wsave);
36 void sint (size_t n, double& r, double& wsave);
37 // single-precision CC FFT
38 void fcffti (size_t n, float& wsave);
39 void fcfftf (size_t n, float& r, float& wsave);
40 void fcfftb (size_t n, float& r, float& wsave);
41 // single-precision RC FFT
42 void frffti (size_t n, float& wsave);
43 void frfftf (size_t n, float& r, float& wsave);
44 void frfftb (size_t n, float& r, float& wsave);
45 // single-precision sine transform
46 void fsinti (size_t n, float& wsave);
47 void fsint (size_t n, float& r, float& wsave);
48}
49
50
51// FFTPACK_wrap provides static functions that wrap the Fortran functions
52// in a common interface. We specialize this class on precision type.
53template <class T>
54class FFTPACK_wrap {};
55
56// Specialization for float
57template <>
58class FFTPACK_wrap<float> {
59
60public:
61 // interface functions used by class FFTPACK
62
63 // initialization functions for CC FFT, RC FFT, and sine transform
64 static void ccffti(size_t n, float* wsave) { fcffti (n, *wsave); }
65 static void rcffti(size_t n, float* wsave) { frffti (n, *wsave); }
66 static void rrffti(size_t n, float* wsave) { fsinti (n, *wsave); }
67 // forward and backward CC FFT
68 static void ccfftf(size_t n, float* r, float* wsave) { fcfftf (n, *r, *wsave); }
69 static void ccfftb(size_t n, float* r, float* wsave) { fcfftb (n, *r, *wsave); }
70 // forward and backward RC FFT
71 static void rcfftf(size_t n, float* r, float* wsave) { frfftf (n, *r, *wsave); }
72 static void rcfftb(size_t n, float* r, float* wsave) { frfftb (n, *r, *wsave); }
73 // sine transform
74 static void rrfft(size_t n, float* r, float* wsave) { fsint (n, *r, *wsave); }
75};
76
77// Specialization for double
78template <>
79class FFTPACK_wrap<double> {
80
81public:
82 // interface functions used by class FFTPACK
83
84 // initialization functions for CC FFT, RC FFT, and sine transform
85 static void ccffti(size_t n, double* wsave) { cffti (n, *wsave); }
86 static void rcffti(size_t n, double* wsave) { rffti (n, *wsave); }
87 static void rrffti(size_t n, double* wsave) { sinti (n, *wsave); }
88 // forward and backward CC FFT
89 static void ccfftf(size_t n, double* r, double* wsave) {cfftf (n, *r, *wsave);}
90 static void ccfftb(size_t n, double* r, double* wsave) {cfftb (n, *r, *wsave);}
91 // forward and backward RC FFT
92 static void rcfftf(size_t n, double* r, double* wsave) {rfftf (n, *r, *wsave);}
93 static void rcfftb(size_t n, double* r, double* wsave) {rfftb (n, *r, *wsave);}
94 // sine transform
95 static void rrfft(size_t n, double* r, double* wsave) { sint (n, *r, *wsave); }
96};
97
98
99// Definition of FFT engine class FFTPACK
100template <class T>
101class FFTPACK {
102
103public:
104
105 // definition of complex type
106 typedef std::complex<T> Complex_t;
107
108 // Trivial constructor. Do the real work in setup function.
109 FFTPACK(void) {}
110
111 // destructor
112 ~FFTPACK(void);
113
114 // setup internal storage and prepare to perform FFTs
115 // inputs are number of dimensions to transform, the transform types,
116 // and the lengths of these dimensions.
117 void setup(unsigned numTransformDims, const int* transformTypes,
118 const int* axisLengths);
119
120 // invoke FFT on complex data for given dimension and direction
121 void callFFT(unsigned transformDim, int direction, Complex_t* data);
122
123 // invoke FFT on real data for given dimension and direction
124 void callFFT(unsigned transformDim, int direction, T* data);
125
126private:
127
128 unsigned numTransformDims_m; // number of dimensions to transform
129 int* transformType_m; // transform type for each dimension
130 int* axisLength_m; // length of each transform dimension
131 T** trig_m; // trigonometric tables
132
133};
134
135
136// Inline member function definitions
137
138// destructor
139template <class T>
140inline
142 // delete storage
143 for (unsigned d=0; d<numTransformDims_m; ++d)
144 delete [] trig_m[d];
145 delete [] trig_m;
146 delete [] transformType_m;
147 delete [] axisLength_m;
148}
149
150// setup internal storage to prepare for FFTs
151template <class T>
152inline void
153FFTPACK<T>::setup(unsigned numTransformDims, const int* transformTypes,
154 const int* axisLengths) {
155
156 // store transform types and lengths for each transform dim
157 numTransformDims_m = numTransformDims;
158 transformType_m = new int[numTransformDims_m];
159 axisLength_m = new int[numTransformDims_m];
160 unsigned d;
161 for (d=0; d<numTransformDims_m; ++d) {
162 transformType_m[d] = transformTypes[d];
163 axisLength_m[d] = axisLengths[d];
164 }
165
166 // allocate and initialize trig table
167 trig_m = new T*[numTransformDims_m];
168 for (d=0; d<numTransformDims_m; ++d) {
169 switch (transformType_m[d]) {
170 case 0: // CC FFT
171 trig_m[d] = new T[4 * axisLength_m[d] + 15];
172 FFTPACK_wrap<T>::ccffti(axisLength_m[d], trig_m[d]);
173 break;
174 case 1: // RC FFT
175 trig_m[d] = new T[2 * axisLength_m[d] + 15];
176 FFTPACK_wrap<T>::rcffti(axisLength_m[d], trig_m[d]);
177 break;
178 case 2: // Sine transform
179 trig_m[d] = new T[static_cast<int>(2.5 * axisLength_m[d] + 0.5) + 15];
180 FFTPACK_wrap<T>::rrffti(axisLength_m[d], trig_m[d]);
181 break;
182 default:
183 ERRORMSG("Unknown transform type requested!!" << endl);
184 break;
185 }
186 }
187
188 return;
189}
190
191// invoke FFT on complex data for given dimension and direction
192template <class T>
193inline void
194FFTPACK<T>::callFFT(unsigned transformDim, int direction,
195 FFTPACK<T>::Complex_t* data) {
196
197 // check transform dimension and direction arguments
198 PAssert_LT(transformDim, numTransformDims_m);
199 PAssert_EQ(std::abs(direction), 1);
200
201 // cast complex number pointer to T* for calling Fortran routines
202 T* rdata = reinterpret_cast<T*>(data);
203
204 // branch on transform type for this dimension
205 switch (transformType_m[transformDim]) {
206 case 0: // CC FFT
207 if (direction == +1) {
208 // call forward complex-to-complex FFT
209 FFTPACK_wrap<T>::ccfftf(axisLength_m[transformDim], rdata,
210 trig_m[transformDim]);
211 }
212 else {
213 // call backward complex-to-complex FFT
214 FFTPACK_wrap<T>::ccfftb(axisLength_m[transformDim], rdata,
215 trig_m[transformDim]);
216 }
217 break;
218 case 1: // RC FFT
219 if (direction == +1) {
220 // call forward real-to-complex FFT
221 FFTPACK_wrap<T>::rcfftf(axisLength_m[transformDim], rdata,
222 trig_m[transformDim]);
223 // rearrange output to conform with SGI format for complex result
224 int clen = axisLength_m[transformDim]/2 + 1;
225 data[clen-1] = Complex_t(imag(data[clen-2]),0.0);
226 for (int i = clen-2; i > 0; --i)
227 data[i] = Complex_t(imag(data[i-1]),real(data[i]));
228 data[0] = Complex_t(real(data[0]),0.0);
229 }
230 else {
231 // rearrange input to conform with Netlib format for complex modes
232 int clen = axisLength_m[transformDim]/2 + 1;
233 data[0] = Complex_t(real(data[0]),real(data[1]));
234 for (int i = 1; i < clen-1; ++i)
235 data[i] = Complex_t(imag(data[i]),real(data[i+1]));
236 // call backward complex-to-real FFT
237 FFTPACK_wrap<T>::rcfftb(axisLength_m[transformDim], rdata,
238 trig_m[transformDim]);
239 }
240 break;
241 case 2: // Sine transform
242 ERRORMSG("Input for real-to-real FFT should be real!!" << endl);
243 break;
244 default:
245 ERRORMSG("Unknown transform type requested!!" << endl);
246 break;
247 }
248
249 return;
250}
251
252// invoke FFT on real data for given dimension and direction
253template <class T>
254inline void
255FFTPACK<T>::callFFT(unsigned transformDim, int direction, T* data) {
256
257 // check transform dimension and direction arguments
258 PAssert_LT(transformDim, numTransformDims_m);
259 // avoid unused variable warning if we compile with Release
260 // :FIXME: remove direction
261 (void)direction;
262 PAssert_EQ(std::abs(direction), 1);
263
264 // branch on transform type for this dimension
265 switch (transformType_m[transformDim]) {
266 case 0: // CC FFT
267 ERRORMSG("Input for complex-to-complex FFT should be complex!!" << endl);
268 break;
269 case 1: // RC FFT
270 ERRORMSG("real-to-complex FFT uses complex input!!" << endl);
271 break;
272 case 2: // Sine transform
273 // invoke the real-to-real transform on the data
274 FFTPACK_wrap<T>::rrfft(axisLength_m[transformDim], data,
275 trig_m[transformDim]);
276 break;
277 default:
278 ERRORMSG("Unknown transform type requested!!" << endl);
279 break;
280 }
281
282 return;
283}
284#endif // IPPL_FFT_FFTPACK_FFT_H
285
FLieGenerator< T, N > real(const FLieGenerator< std::complex< T >, N > &)
Take real part of a complex generator.
FLieGenerator< T, N > imag(const FLieGenerator< std::complex< T >, N > &)
Take imaginary part of a complex generator.
void frfftf(size_t n, float &r, float &wsave)
void frffti(size_t n, float &wsave)
void rfftf(size_t n, double &r, double &wsave)
void fsint(size_t n, float &r, float &wsave)
void cffti(size_t n, double &wsave)
void fsinti(size_t n, float &wsave)
void cfftf(size_t n, double &r, double &wsave)
void sint(size_t n, double &r, double &wsave)
void sinti(size_t n, double &wsave)
void cfftb(size_t n, double &r, double &wsave)
void rfftb(size_t n, double &r, double &wsave)
void fcffti(size_t n, float &wsave)
void fcfftf(size_t n, float &r, float &wsave)
void rffti(size_t n, double &wsave)
void fcfftb(size_t n, float &r, float &wsave)
void frfftb(size_t n, float &r, float &wsave)
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
#define PAssert_LT(a, b)
Definition: PAssert.h:106
#define PAssert_EQ(a, b)
Definition: PAssert.h:104
#define ERRORMSG(msg)
Definition: IpplInfo.h:350
static void ccfftb(size_t n, float *r, float *wsave)
Definition: fftpack_FFT.h:69
static void rcfftb(size_t n, float *r, float *wsave)
Definition: fftpack_FFT.h:72
static void rcffti(size_t n, float *wsave)
Definition: fftpack_FFT.h:65
static void ccfftf(size_t n, float *r, float *wsave)
Definition: fftpack_FFT.h:68
static void rcfftf(size_t n, float *r, float *wsave)
Definition: fftpack_FFT.h:71
static void rrfft(size_t n, float *r, float *wsave)
Definition: fftpack_FFT.h:74
static void ccffti(size_t n, float *wsave)
Definition: fftpack_FFT.h:64
static void rrffti(size_t n, float *wsave)
Definition: fftpack_FFT.h:66
static void ccfftb(size_t n, double *r, double *wsave)
Definition: fftpack_FFT.h:90
static void rrfft(size_t n, double *r, double *wsave)
Definition: fftpack_FFT.h:95
static void rrffti(size_t n, double *wsave)
Definition: fftpack_FFT.h:87
static void rcffti(size_t n, double *wsave)
Definition: fftpack_FFT.h:86
static void rcfftb(size_t n, double *r, double *wsave)
Definition: fftpack_FFT.h:93
static void rcfftf(size_t n, double *r, double *wsave)
Definition: fftpack_FFT.h:92
static void ccfftf(size_t n, double *r, double *wsave)
Definition: fftpack_FFT.h:89
static void ccffti(size_t n, double *wsave)
Definition: fftpack_FFT.h:85
int * axisLength_m
Definition: fftpack_FFT.h:130
~FFTPACK(void)
Definition: fftpack_FFT.h:141
T ** trig_m
Definition: fftpack_FFT.h:131
std::complex< T > Complex_t
Definition: fftpack_FFT.h:106
int * transformType_m
Definition: fftpack_FFT.h:129
void callFFT(unsigned transformDim, int direction, Complex_t *data)
unsigned numTransformDims_m
Definition: fftpack_FFT.h:128
FFTPACK(void)
Definition: fftpack_FFT.h:109
void setup(unsigned numTransformDims, const int *transformTypes, const int *axisLengths)
Definition: fftpack_FFT.h:153