Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TFFTReal.cxx
Go to the documentation of this file.
1// @(#)root/fft:$Id$
2// Author: Anna Kreshuk 07/4/2006
3
4/*************************************************************************
5 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12////////////////////////////////////////////////////////////////////////////////
13/// \class TFFTReal
14/// One of the interface classes to the FFTW package, can be used directly
15/// or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
16///
17/// Computes transforms called r2r in FFTW manual:
18/// - transforms of real input and output in "halfcomplex" format i.e.
19/// real and imaginary parts for a transform of size n stored as
20/// (r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1)
21/// - discrete Hartley transform
22/// - sine and cosine transforms (DCT-I,II,III,IV and DST-I,II,III,IV)
23/// For the detailed information on the computed
24/// transforms please refer to the FFTW manual, chapter "What FFTW really computes".
25///
26/// How to use it:
27/// 1. Create an instance of TFFTReal - this will allocate input and output
28/// arrays (unless an in-place transform is specified)
29/// 2. Run the Init() function with the desired flags and settings (see function
30/// comments for possible kind parameters)
31/// 3. Set the data (via SetPoints()or SetPoint() functions)
32/// 4. Run the Transform() function
33/// 5. Get the output (via GetPoints() or GetPoint() functions)
34/// 6. Repeat steps 3)-5) as needed
35///
36/// For a transform of the same size, but of different kind (or with different flags),
37/// rerun the Init() function and continue with steps 3)-5)
38///
39/// NOTE:
40/// 1. running Init() function will overwrite the input array! Don't set any data
41/// before running the Init() function!
42/// 2. FFTW computes unnormalized transform, so doing a transform followed by
43/// its inverse will lead to the original array scaled BY:
44/// - transform size (N) for R2HC, HC2R, DHT transforms
45/// - 2*(N-1) for DCT-I (REDFT00)
46/// - 2*(N+1) for DST-I (RODFT00)
47/// - 2*N for the remaining transforms
48///
49/// Transform inverses:
50/// - R2HC<-->HC2R
51/// - DHT<-->DHT
52/// - DCT-I<-->DCT-I
53/// - DCT-II<-->DCT-III
54/// - DCT-IV<-->DCT-IV
55/// - DST-I<-->DST-I
56/// - DST-II<-->DST-III
57/// - DST-IV<-->DST-IV
58///
59////////////////////////////////////////////////////////////////////////////////
60
61#include "TFFTReal.h"
62#include "fftw3.h"
63
64
65////////////////////////////////////////////////////////////////////////////////
66///default
67
69{
70 fIn = nullptr;
71 fOut = nullptr;
72 fPlan = nullptr;
73 fN = nullptr;
74 fKind = nullptr;
75 fNdim = 0;
76 fTotalSize = 0;
77}
78
79////////////////////////////////////////////////////////////////////////////////
80///For 1d transforms
81///n here is the physical size of the transform (see FFTW manual for more details)
82
84{
85 fIn = fftw_malloc(sizeof(Double_t)*n);
86 if (inPlace) fOut = nullptr;
87 else fOut = fftw_malloc(sizeof(Double_t)*n);
88
89 fPlan = nullptr;
90 fNdim = 1;
91 fN = new Int_t[1];
92 fN[0] = n;
93 fKind = nullptr;
94 fTotalSize = n;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98///For multidimensional transforms
99///1st parameter is the # of dimensions,
100///2nd is the sizes (physical) of the transform in each dimension
101
103{
104 fTotalSize = 1;
105 fNdim = ndim;
106 fN = new Int_t[ndim];
107 fKind = nullptr;
108 fPlan = nullptr;
109 for (Int_t i=0; i<ndim; i++){
110 fTotalSize*=n[i];
111 fN[i] = n[i];
112 }
114 if (!inPlace)
116 else
117 fOut = nullptr;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121///clean-up
122
124{
126 fPlan = nullptr;
127 fftw_free(fIn);
128 fIn = nullptr;
129 if (fOut){
131 }
132 fOut = nullptr;
133 if (fN)
134 delete [] fN;
135 fN = nullptr;
136 if (fKind)
138 fKind = nullptr;
139}
140
141////////////////////////////////////////////////////////////////////////////////
142///Creates the fftw-plan
143///
144///NOTE: input and output arrays are overwritten during initialisation,
145/// so don't set any points, before running this function!!!!!
146///
147/// #### 1st parameter:
148/// Possible flag_options:
149///
150/// - "ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
151/// performance
152/// - "M" (from "measure") - some time spend in finding the optimal way to do the transform
153/// - "P" (from "patient") - more time spend in finding the optimal way to do the transform
154/// - "EX" (from "exhaustive") - the most optimal way is found
155///
156/// This option should be chosen depending on how many transforms of the same size and
157/// type are going to be done. Planning is only done once, for the first transform of this
158/// size and type.
159///
160/// #### 2nd parameter:
161/// is dummy and doesn't need to be specified
162///
163/// #### 3rd parameter:
164/// transform kind for each dimension
165/// 4 different kinds of sine and cosine transforms are available
166/// - DCT-I - kind=0
167/// - DCT-II - kind=1
168/// - DCT-III - kind=2
169/// - DCT-IV - kind=3
170/// - DST-I - kind=4
171/// - DST-II - kind=5
172/// - DSTIII - kind=6
173/// - DSTIV - kind=7
174
175void TFFTReal::Init( Option_t* flags,Int_t /*sign*/, const Int_t *kind)
176{
177 if (fPlan)
179 fPlan = nullptr;
180
181 if (!fKind)
183
184 if (MapOptions(kind)){
185 if (fOut)
187 else
189 fFlags = flags;
190 }
191}
192
193////////////////////////////////////////////////////////////////////////////////
194///Computes the transform, specified in Init() function
195
197{
198 if (fPlan)
200 else {
201 Error("Transform", "transform hasn't been initialised");
202 return;
203 }
204}
205
206////////////////////////////////////////////////////////////////////////////////
207///Returns the type of the transform
208
210{
211 if (!fKind) {
212 Error("GetType", "Type not defined yet (kind not set)");
213 return "";
214 }
215 if (((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC) return "R2HC";
216 if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R) return "HC2R";
217 if (((fftw_r2r_kind*)fKind)[0]==FFTW_DHT) return "DHT";
218 else return "R2R";
219}
220
221////////////////////////////////////////////////////////////////////////////////
222///Copies the output (or input) points into the provided array, that should
223///be big enough
224
226{
227 const Double_t * array = GetPointsReal(fromInput);
228 if (!array) return;
229 std::copy(array, array+fTotalSize, data);
230}
231
232////////////////////////////////////////////////////////////////////////////////
233///For 1d transformations. Returns point `#ipoint`
234
236{
238 Error("GetPointReal", "No such point");
239 return 0;
240 }
241 const Double_t * array = GetPointsReal(fromInput);
242 return ( array ) ? array[ipoint] : 0;
243}
244
245////////////////////////////////////////////////////////////////////////////////
246///For multidim.transforms. Returns point `#ipoint`
247
249{
250 Int_t ireal = ipoint[0];
251 for (Int_t i=0; i<fNdim-1; i++)
252 ireal=fN[i+1]*ireal + ipoint[i+1];
253
254 const Double_t * array = GetPointsReal(fromInput);
255 return ( array ) ? array[ireal] : 0;
256}
257
258////////////////////////////////////////////////////////////////////////////////
259///Only for input of HC2R and output of R2HC
260
262{
263 const Double_t * array = GetPointsReal(fromInput);
264 if (!array) return;
265 if ( ( ((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC && !fromInput ) ||
266 ( ((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R && fromInput ) )
267 {
268 if (ipoint<fN[0]/2+1){
269 re = array[ipoint];
270 im = array[fN[0]-ipoint];
271 } else {
272 re = array[fN[0]-ipoint];
273 im = -array[ipoint];
274 }
275 if ((fN[0]%2)==0 && ipoint==fN[0]/2) im = 0;
276 }
277}
278////////////////////////////////////////////////////////////////////////////////
279///Only for input of HC2R and output of R2HC and for 1d
280
285
286////////////////////////////////////////////////////////////////////////////////
287///Returns the output (or input) array
288/// we have 4 different cases:
289///
290/// - fromInput = false; fOut = !NULL (transformed is not in place) : return fOut
291/// - fromInput = false; fOut = NULL (transformed is in place) : return fIn
292/// - fromInput = true; fOut = !NULL : return fIn
293/// - fromInput = true; fOut = NULL return an error since input array is overwritten
294
296{
297
298 if (!fromInput && fOut)
299 return (Double_t*)fOut;
300 else if (fromInput && !fOut) {
301 Error("GetPointsReal","Input array was destroyed");
302 return nullptr;
303 }
304 //R__ASSERT(fIn);
305 return (Double_t*)fIn;
306}
307
308////////////////////////////////////////////////////////////////////////////////
309
311{
313 Error("SetPoint", "illegal point index");
314 return;
315 }
316 if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R){
317 if ((fN[0]%2)==0 && ipoint==fN[0]/2)
318 ((Double_t*)fIn)[ipoint] = re;
319 else {
320 ((Double_t*)fIn)[ipoint] = re;
321 ((Double_t*)fIn)[fN[0]-ipoint]=im;
322 }
323 }
324 else
325 ((Double_t*)fIn)[ipoint]=re;
326}
327
328////////////////////////////////////////////////////////////////////////////////
329///Since multidimensional R2HC and HC2R transforms are not supported,
330///third parameter is dummy
331
333{
334 Int_t ireal = ipoint[0];
335 for (Int_t i=0; i<fNdim-1; i++)
336 ireal=fN[i+1]*ireal + ipoint[i+1];
338 Error("SetPoint", "illegal point index");
339 return;
340 }
341 ((Double_t*)fIn)[ireal]=re;
342}
343
344////////////////////////////////////////////////////////////////////////////////
345///Sets all points
346
348{
349 for (Int_t i=0; i<fTotalSize; i++)
350 ((Double_t*)fIn)[i] = data[i];
351}
352
353////////////////////////////////////////////////////////////////////////////////
354///transfers the r2r_kind parameters to fftw type
355
357{
358 if (kind[0] == 10){
359 if (fNdim>1){
360 Error("Init", "Multidimensional R2HC transforms are not supported, use R2C interface instead");
361 return 0;
362 }
364 }
365 else if (kind[0] == 11) {
366 if (fNdim>1){
367 Error("Init", "Multidimensional HC2R transforms are not supported, use C2R interface instead");
368 return 0;
369 }
371 }
372 else if (kind[0] == 12) {
373 for (Int_t i=0; i<fNdim; i++)
375 }
376 else {
377 for (Int_t i=0; i<fNdim; i++){
378 switch (kind[i]) {
379 case 0: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT00; break;
380 case 1: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT01; break;
381 case 2: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT10; break;
382 case 3: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT11; break;
383 case 4: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT00; break;
384 case 5: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT01; break;
385 case 6: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT10; break;
386 case 7: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT11; break;
387 default:
388 ((fftw_r2r_kind*)fKind)[i] = FFTW_R2HC; break;
389 }
390 }
391 }
392 return 1;
393}
394
395////////////////////////////////////////////////////////////////////////////////
396///allowed options:
397/// - "ES" - FFTW_ESTIMATE
398/// - "M" - FFTW_MEASURE
399/// - "P" - FFTW_PATIENT
400/// - "EX" - FFTW_EXHAUSTIVE
401
403{
404 TString opt = flag;
405 opt.ToUpper();
406 if (opt.Contains("ES"))
407 return FFTW_ESTIMATE;
408 if (opt.Contains("M"))
409 return FFTW_MEASURE;
410 if (opt.Contains("P"))
411 return FFTW_PATIENT;
412 if (opt.Contains("EX"))
413 return FFTW_EXHAUSTIVE;
414 return FFTW_ESTIMATE;
415}
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Int_t fNdim
Definition TFFTReal.h:25
void Transform() override
Computes the transform, specified in Init() function.
Definition TFFTReal.cxx:196
void SetPoint(Int_t ipoint, Double_t re, Double_t im=0) override
Definition TFFTReal.cxx:310
void * fPlan
Definition TFFTReal.h:24
UInt_t MapFlag(Option_t *flag)
allowed options:
Definition TFFTReal.cxx:402
Double_t * GetPointsReal(Bool_t fromInput=kFALSE) const override
Returns the output (or input) array we have 4 different cases:
Definition TFFTReal.cxx:295
Int_t fTotalSize
Definition TFFTReal.h:26
void SetPoints(const Double_t *data) override
Sets all points.
Definition TFFTReal.cxx:347
void GetPoints(Double_t *data, Bool_t fromInput=kFALSE) const override
Copies the output (or input) points into the provided array, that should be big enough.
Definition TFFTReal.cxx:225
Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const override
For 1d transformations. Returns point #ipoint
Definition TFFTReal.cxx:235
Int_t * fN
Definition TFFTReal.h:27
Option_t * GetType() const override
Returns the type of the transform.
Definition TFFTReal.cxx:209
TString fFlags
Definition TFFTReal.h:29
void * fOut
Definition TFFTReal.h:23
TFFTReal()
default
Definition TFFTReal.cxx:68
void * fIn
Definition TFFTReal.h:22
Int_t MapOptions(const Int_t *kind)
transfers the r2r_kind parameters to fftw type
Definition TFFTReal.cxx:356
void Init(Option_t *flags, Int_t sign, const Int_t *kind) override
Creates the fftw-plan.
Definition TFFTReal.cxx:175
void * fKind
Definition TFFTReal.h:28
void GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const override
Only for input of HC2R and output of R2HC and for 1d.
Definition TFFTReal.cxx:281
~TFFTReal() override
clean-up
Definition TFFTReal.cxx:123
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
Basic string class.
Definition TString.h:138
void ToUpper()
Change string to upper case.
Definition TString.cxx:1202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
const Int_t n
Definition legend1.C:16