Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TFFTRealComplex.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///
14/// \class TFFTRealComplex
15///
16/// One of the interface classes to the FFTW package, can be used directly
17/// or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
18///
19/// Computes a real input/complex output discrete Fourier transform in 1 or more
20/// dimensions. However, only out-of-place transforms are now supported for transforms
21/// in more than 1 dimension. For detailed information about the computed transforms,
22/// please refer to the FFTW manual
23///
24/// How to use it:
25/// 1. Create an instance of TFFTRealComplex - this will allocate input and output
26/// arrays (unless an in-place transform is specified)
27/// 2. Run the Init() function with the desired flags and settings (see function
28/// comments for possible kind parameters)
29/// 3. Set the data (via SetPoints()or SetPoint() functions)
30/// 4. Run the Transform() function
31/// 5. Get the output (via GetPoints() or GetPoint() functions)
32/// 6. Repeat steps 3)-5) as needed
33/// For a transform of the same size, but with different flags,
34/// rerun the Init() function and continue with steps 3)-5)
35///
36/// NOTE:
37/// 1. running Init() function will overwrite the input array! Don't set any data
38/// before running the Init() function
39/// 2. FFTW computes unnormalized transform, so doing a transform followed by
40/// its inverse will lead to the original array scaled by the transform size
41///
42///
43/////////////////////////////////////////////////////////////////////////////////
44
45#include "TFFTRealComplex.h"
46#include "fftw3.h"
47#include "TComplex.h"
48
49
50
51////////////////////////////////////////////////////////////////////////////////
52///default
53
55{
56 fIn = nullptr;
57 fOut = nullptr;
58 fPlan = nullptr;
59 fN = nullptr;
60 fNdim = 0;
61 fTotalSize = 0;
62}
63
64////////////////////////////////////////////////////////////////////////////////
65///For 1d transforms
66///Allocates memory for the input array, and, if inPlace = kFALSE, for the output array
67
69{
70
71 if (!inPlace){
72 fIn = fftw_malloc(sizeof(Double_t)*n);
73 fOut = fftw_malloc(sizeof(fftw_complex)*(n/2+1));
74 } else {
75 fIn = fftw_malloc(sizeof(Double_t)*(2*(n/2+1)));
76 fOut = nullptr;
77 }
78 fN = new Int_t[1];
79 fN[0] = n;
80 fTotalSize = n;
81 fNdim = 1;
82 fPlan = nullptr;
83}
84
85////////////////////////////////////////////////////////////////////////////////
86///For ndim-dimensional transforms
87///Second argument contains sizes of the transform in each dimension
88
90{
91 if (ndim>1 && inPlace==kTRUE){
92 Error("TFFTRealComplex", "multidimensional in-place r2c transforms are not implemented yet");
93 return;
94 }
95 fNdim = ndim;
96 fTotalSize = 1;
97 fN = new Int_t[fNdim];
98 for (Int_t i=0; i<fNdim; i++){
99 fN[i] = n[i];
100 fTotalSize*=n[i];
101 }
102 Int_t sizeout = Int_t(Double_t(fTotalSize)*(n[ndim-1]/2+1)/n[ndim-1]);
103 if (!inPlace){
106 } else {
107 fIn = fftw_malloc(sizeof(Double_t)*(2*sizeout));
108 fOut = nullptr;
109 }
110 fPlan = nullptr;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114///Destroys the data arrays and the plan. However, some plan information stays around
115///until the root session is over, and is reused if other plans of the same size are
116///created
117
119{
121 fPlan = nullptr;
122 fftw_free(fIn);
123 fIn = nullptr;
124 if (fOut)
126 fOut = nullptr;
127 if (fN)
128 delete [] fN;
129 fN = nullptr;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133///Creates the fftw-plan
134///
135///NOTE: input and output arrays are overwritten during initialisation,
136/// so don't set any points, before running this function!!!!!
137///
138///Arguments sign and kind are dummy and not need to be specified
139///Possible flag_options:
140///
141/// - "ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
142/// performance
143/// - "M" (from "measure") - some time spend in finding the optimal way to do the transform
144/// - "P" (from "patient") - more time spend in finding the optimal way to do the transform
145/// - "EX" (from "exhaustive") - the most optimal way is found
146///
147///This option should be chosen depending on how many transforms of the same size and
148///type are going to be done. Planning is only done once, for the first transform of this
149///size and type.
150
151void TFFTRealComplex::Init(Option_t *flags,Int_t /*sign*/, const Int_t* /*kind*/)
152{
153 fFlags = flags;
154
155 if (fPlan)
157 fPlan = nullptr;
158
159 if (fOut)
161 else
163}
164
165////////////////////////////////////////////////////////////////////////////////
166///Computes the transform, specified in Init() function
167
169{
170
171 if (fPlan){
173 }
174 else {
175 Error("Transform", "transform hasn't been initialised");
176 return;
177 }
178}
179
180////////////////////////////////////////////////////////////////////////////////
181///Fills the array data with the computed transform.
182///Only (roughly) a half of the transform is copied (exactly the output of FFTW),
183///the rest being Hermitian symmetric with the first half
184
186{
187 if (fromInput){
188 for (Int_t i=0; i<fTotalSize; i++)
189 data[i] = ((Double_t*)fIn)[i];
190 } else {
191 Int_t realN = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
192 if (fOut){
193 for (Int_t i=0; i<realN; i+=2){
194 data[i] = ((fftw_complex*)fOut)[i/2][0];
195 data[i+1] = ((fftw_complex*)fOut)[i/2][1];
196 }
197 }
198 else {
199 for (Int_t i=0; i<realN; i++)
200 data[i] = ((Double_t*)fIn)[i];
201 }
202 }
203}
204
205////////////////////////////////////////////////////////////////////////////////
206///Returns the real part of the point `#ipoint` from the output or the point `#ipoint`
207///from the input
208
210{
211 if (fOut && !fromInput){
212 Warning("GetPointReal", "Output is complex. Only real part returned");
213 return ((fftw_complex*)fOut)[ipoint][0];
214 }
215 else
216 return ((Double_t*)fIn)[ipoint];
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// Returns the real part of the point `#ipoint` from the output or the point
221/// `#ipoint` from the input
222
224{
225 Int_t ireal = ipoint[0];
226 for (Int_t i=0; i<fNdim-1; i++)
227 ireal=fN[i+1]*ireal + ipoint[i+1];
228
229 if (fOut && !fromInput){
230 Warning("GetPointReal", "Output is complex. Only real part returned");
231 return ((fftw_complex*)fOut)[ireal][0];
232 }
233 else
234 return ((Double_t*)fIn)[ireal];
235}
236
237
238////////////////////////////////////////////////////////////////////////////////
239///Returns the point `#ipoint`.
240///For 1d, if ipoint > fN/2+1 (the point is in the Hermitian symmetric part), it is still
241///returned. For >1d, only the first (roughly)half of points can be returned
242///For 2d, see function GetPointComplex(Int_t *ipoint,...)
243
245{
246 if (fromInput){
247 re = ((Double_t*)fIn)[ipoint];
248 } else {
249 if (fNdim==1){
250 if (fOut){
251 if (ipoint<fN[0]/2+1){
252 re = ((fftw_complex*)fOut)[ipoint][0];
253 im = ((fftw_complex*)fOut)[ipoint][1];
254 } else {
255 re = ((fftw_complex*)fOut)[fN[0]-ipoint][0];
256 im = -((fftw_complex*)fOut)[fN[0]-ipoint][1];
257 }
258 } else {
259 if (ipoint<fN[0]/2+1){
260 re = ((Double_t*)fIn)[2*ipoint];
261 im = ((Double_t*)fIn)[2*ipoint+1];
262 } else {
263 re = ((Double_t*)fIn)[2*(fN[0]-ipoint)];
264 im = ((Double_t*)fIn)[2*(fN[0]-ipoint)+1];
265 }
266 }
267 }
268 else {
269 Int_t realN = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
270 if (ipoint>realN/2){
271 Error("GetPointComplex", "Illegal index value");
272 return;
273 }
274 if (fOut){
275 re = ((fftw_complex*)fOut)[ipoint][0];
276 im = ((fftw_complex*)fOut)[ipoint][1];
277 } else {
278 re = ((Double_t*)fIn)[2*ipoint];
279 im = ((Double_t*)fIn)[2*ipoint+1];
280 }
281 }
282 }
283}
284////////////////////////////////////////////////////////////////////////////////
285///For multidimensional transforms. Returns the point `#ipoint`.
286///In case of transforms of more than 2 dimensions,
287///only points from the first (roughly)half are returned, the rest being Hermitian symmetric
288
290{
291
292 Int_t ireal = ipoint[0];
293 for (Int_t i=0; i<fNdim-2; i++)
294 ireal=fN[i+1]*ireal + ipoint[i+1];
295 //special treatment of the last dimension
296 ireal = (fN[fNdim-1]/2+1)*ireal + ipoint[fNdim-1];
297
298 if (fromInput){
299 re = ((Double_t*)fIn)[ireal];
300 return;
301 }
302 if (fNdim==1){
303 if (fOut){
304 if (ipoint[0] <fN[0]/2+1){
305 re = ((fftw_complex*)fOut)[ipoint[0]][0];
306 im = ((fftw_complex*)fOut)[ipoint[0]][1];
307 } else {
308 re = ((fftw_complex*)fOut)[fN[0]-ipoint[0]][0];
309 im = -((fftw_complex*)fOut)[fN[0]-ipoint[0]][1];
310 }
311 } else {
312 if (ireal <fN[0]/2+1){
313 re = ((Double_t*)fIn)[2*ipoint[0]];
314 im = ((Double_t*)fIn)[2*ipoint[0]+1];
315 } else {
316 re = ((Double_t*)fIn)[2*(fN[0]-ipoint[0])];
317 im = ((Double_t*)fIn)[2*(fN[0]-ipoint[0])+1];
318 }
319 }
320 }
321 else if (fNdim==2){
322 if (fOut){
323 if (ipoint[1]<fN[1]/2+1){
324 re = ((fftw_complex*)fOut)[ipoint[1]+(fN[1]/2+1)*ipoint[0]][0];
325 im = ((fftw_complex*)fOut)[ipoint[1]+(fN[1]/2+1)*ipoint[0]][1];
326 } else {
327 if (ipoint[0]==0){
328 re = ((fftw_complex*)fOut)[fN[1]-ipoint[1]][0];
329 im = -((fftw_complex*)fOut)[fN[1]-ipoint[1]][1];
330 } else {
331 re = ((fftw_complex*)fOut)[(fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0])][0];
332 im = -((fftw_complex*)fOut)[(fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0])][1];
333 }
334 }
335 } else {
336 if (ipoint[1]<fN[1]/2+1){
337 re = ((Double_t*)fIn)[2*(ipoint[1]+(fN[1]/2+1)*ipoint[0])];
338 im = ((Double_t*)fIn)[2*(ipoint[1]+(fN[1]/2+1)*ipoint[0])+1];
339 } else {
340 if (ipoint[0]==0){
341 re = ((Double_t*)fIn)[2*(fN[1]-ipoint[1])];
342 im = -((Double_t*)fIn)[2*(fN[1]-ipoint[1])+1];
343 } else {
344 re = ((Double_t*)fIn)[2*((fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0]))];
345 im = -((Double_t*)fIn)[2*((fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0]))+1];
346 }
347 }
348 }
349 }
350 else {
351
352 if (fOut){
353 re = ((fftw_complex*)fOut)[ireal][0];
354 im = ((fftw_complex*)fOut)[ireal][1];
355 } else {
356 re = ((Double_t*)fIn)[2*ireal];
357 im = ((Double_t*)fIn)[2*ireal+1];
358 }
359 }
360}
361
362
363////////////////////////////////////////////////////////////////////////////////
364///Returns the input array// One of the interface classes to the FFTW package, can be used directly
365/// or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
366
368{
369 if (!fromInput){
370 Error("GetPointsReal", "Output array is complex");
371 return nullptr;
372 }
373 return (Double_t*)fIn;
374}
375
376////////////////////////////////////////////////////////////////////////////////
377///Fills the argument arrays with the real and imaginary parts of the computed transform.
378///Only (roughly) a half of the transform is copied, the rest being Hermitian
379///symmetric with the first half
380
382{
383 Int_t nreal;
384 if (fOut && !fromInput){
385 nreal = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
386 for (Int_t i=0; i<nreal; i++){
387 re[i] = ((fftw_complex*)fOut)[i][0];
388 im[i] = ((fftw_complex*)fOut)[i][1];
389 }
390 } else if (fromInput) {
391 for (Int_t i=0; i<fTotalSize; i++){
392 re[i] = ((Double_t *)fIn)[i];
393 im[i] = 0;
394 }
395 }
396 else {
397 nreal = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
398 for (Int_t i=0; i<nreal; i+=2){
399 re[i/2] = ((Double_t*)fIn)[i];
400 im[i/2] = ((Double_t*)fIn)[i+1];
401 }
402 }
403}
404
405////////////////////////////////////////////////////////////////////////////////
406///Fills the argument arrays with the real and imaginary parts of the computed transform.
407///Only (roughly) a half of the transform is copied, the rest being Hermitian
408///symmetric with the first half
409
411{
412 Int_t nreal;
413
414 if (fOut && !fromInput){
415 nreal = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
416
417 for (Int_t i=0; i<nreal; i+=2){
418 data[i] = ((fftw_complex*)fOut)[i/2][0];
419 data[i+1] = ((fftw_complex*)fOut)[i/2][1];
420 }
421 } else if (fromInput){
422 for (Int_t i=0; i<fTotalSize; i+=2){
423 data[i] = ((Double_t*)fIn)[i/2];
424 data[i+1] = 0;
425 }
426 }
427 else {
428
429 nreal = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
430 for (Int_t i=0; i<nreal; i++)
431 data[i] = ((Double_t*)fIn)[i];
432 }
433}
434
435////////////////////////////////////////////////////////////////////////////////
436///Set the point `#ipoint`
437
439{
440 ((Double_t *)fIn)[ipoint] = re;
441}
442
443////////////////////////////////////////////////////////////////////////////////
444///For multidimensional transforms. Set the point `#ipoint`
445
447{
448 Int_t ireal = ipoint[0];
449 for (Int_t i=0; i<fNdim-1; i++)
450 ireal=fN[i+1]*ireal + ipoint[i+1];
451 ((Double_t*)fIn)[ireal]=re;
452}
453
454////////////////////////////////////////////////////////////////////////////////
455///Set all input points
456
458{
459 for (Int_t i=0; i<fTotalSize; i++){
460 ((Double_t*)fIn)[i]=data[i];
461 }
462}
463
464////////////////////////////////////////////////////////////////////////////////
465///Sets the point `#ipoint` (only the real part of the argument is taken)
466
471
472////////////////////////////////////////////////////////////////////////////////
473///Set all points. Only the real array is used
474
476{
477 for (Int_t i=0; i<fTotalSize; i++)
478 ((Double_t *)fIn)[i] = re[i];
479}
480
481////////////////////////////////////////////////////////////////////////////////
482///allowed options:
483///"ES"
484///"M"
485///"P"
486///"EX"
487
489{
490
491 TString opt = flag;
492 opt.ToUpper();
493 if (opt.Contains("ES"))
494 return FFTW_ESTIMATE;
495 if (opt.Contains("M"))
496 return FFTW_MEASURE;
497 if (opt.Contains("P"))
498 return FFTW_PATIENT;
499 if (opt.Contains("EX"))
500 return FFTW_EXHAUSTIVE;
501 return FFTW_ESTIMATE;
502}
#define c(i)
Definition RSha256.hxx:101
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
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
UInt_t MapFlag(Option_t *flag)
allowed options: "ES" "M" "P" "EX"
void Init(Option_t *flags, Int_t, const Int_t *) override
Creates the fftw-plan.
void GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const override
Returns the point #ipoint.
void SetPoints(const Double_t *data) override
Set all input points.
Double_t * GetPointsReal(Bool_t fromInput=kFALSE) const override
Returns the input array// One of the interface classes to the FFTW package, can be used directly or v...
Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const override
Returns the real part of the point #ipoint from the output or the point #ipoint from the input.
void SetPointComplex(Int_t ipoint, TComplex &c) override
Sets the point #ipoint (only the real part of the argument is taken)
void GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput=kFALSE) const override
Fills the argument arrays with the real and imaginary parts of the computed transform.
void SetPointsComplex(const Double_t *re, const Double_t *im) override
Set all points. Only the real array is used.
void SetPoint(Int_t ipoint, Double_t re, Double_t im=0) override
Set the point #ipoint
void GetPoints(Double_t *data, Bool_t fromInput=kFALSE) const override
Fills the array data with the computed transform.
void Transform() override
Computes the transform, specified in Init() function.
~TFFTRealComplex() override
Destroys the data arrays and the plan.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1057
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