00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059 #include "TFFTReal.h"
00060 #include "fftw3.h"
00061
00062 ClassImp(TFFTReal)
00063
00064
00065 TFFTReal::TFFTReal()
00066 {
00067
00068
00069 fIn = 0;
00070 fOut = 0;
00071 fPlan = 0;
00072 fN = 0;
00073 fKind = 0;
00074 fFlags = 0;
00075 fNdim = 0;
00076 fTotalSize = 0;
00077 }
00078
00079
00080 TFFTReal::TFFTReal(Int_t n, Bool_t inPlace)
00081 {
00082
00083
00084
00085 fIn = fftw_malloc(sizeof(Double_t)*n);
00086 if (inPlace) fOut = 0;
00087 else fOut = fftw_malloc(sizeof(Double_t)*n);
00088
00089 fPlan = 0;
00090 fNdim = 1;
00091 fN = new Int_t[1];
00092 fN[0] = n;
00093 fKind = 0;
00094 fTotalSize = n;
00095 fFlags = 0;
00096 }
00097
00098
00099 TFFTReal::TFFTReal(Int_t ndim, Int_t *n, Bool_t inPlace)
00100 {
00101
00102
00103
00104
00105 fTotalSize = 1;
00106 fNdim = ndim;
00107 fN = new Int_t[ndim];
00108 fKind = 0;
00109 fPlan = 0;
00110 fFlags = 0;
00111 for (Int_t i=0; i<ndim; i++){
00112 fTotalSize*=n[i];
00113 fN[i] = n[i];
00114 }
00115 fIn = fftw_malloc(sizeof(Double_t)*fTotalSize);
00116 if (!inPlace)
00117 fOut = fftw_malloc(sizeof(Double_t)*fTotalSize);
00118 else
00119 fOut = 0;
00120 }
00121
00122
00123 TFFTReal::~TFFTReal()
00124 {
00125
00126
00127 fftw_destroy_plan((fftw_plan)fPlan);
00128 fPlan = 0;
00129 fftw_free(fIn);
00130 fIn = 0;
00131 if (fOut){
00132 fftw_free(fOut);
00133 }
00134 fOut = 0;
00135 if (fN)
00136 delete [] fN;
00137 fN = 0;
00138 if (fKind)
00139 fftw_free((fftw_r2r_kind*)fKind);
00140 fKind = 0;
00141 }
00142
00143
00144 void TFFTReal::Init( Option_t* flags,Int_t , const Int_t *kind)
00145 {
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173 if (!fKind)
00174 fKind = (fftw_r2r_kind*)fftw_malloc(sizeof(fftw_r2r_kind)*fNdim);
00175
00176 if (MapOptions(kind)){
00177 if (fOut)
00178 fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fOut, (fftw_r2r_kind*)fKind, MapFlag(flags));
00179 else
00180 fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fIn, (fftw_r2r_kind*)fKind, MapFlag(flags));
00181 fFlags = flags;
00182 }
00183 }
00184
00185
00186 void TFFTReal::Transform()
00187 {
00188
00189
00190 if (fPlan)
00191 fftw_execute((fftw_plan)fPlan);
00192 else {
00193 Error("Transform", "transform hasn't been initialised");
00194 return;
00195 }
00196 }
00197
00198
00199 Option_t *TFFTReal::GetType() const
00200 {
00201
00202
00203 if (!fKind) {
00204 Error("GetType", "Type not defined yet (kind not set)");
00205 return "";
00206 }
00207 if (((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC) return "R2HC";
00208 if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R) return "HC2R";
00209 if (((fftw_r2r_kind*)fKind)[0]==FFTW_DHT) return "DHT";
00210 else return "R2R";
00211 }
00212
00213
00214 void TFFTReal::GetPoints(Double_t *data, Bool_t fromInput) const
00215 {
00216
00217
00218
00219 const Double_t * array = GetPointsReal(fromInput);
00220 if (!array) return;
00221 std::copy(array, array+fTotalSize, data);
00222 }
00223
00224
00225 Double_t TFFTReal::GetPointReal(Int_t ipoint, Bool_t fromInput) const
00226 {
00227
00228
00229 if (ipoint<0 || ipoint>fTotalSize){
00230 Error("GetPointReal", "No such point");
00231 return 0;
00232 }
00233 const Double_t * array = GetPointsReal(fromInput);
00234 return ( array ) ? array[ipoint] : 0;
00235 }
00236
00237
00238 Double_t TFFTReal::GetPointReal(const Int_t *ipoint, Bool_t fromInput) const
00239 {
00240
00241
00242 Int_t ireal = ipoint[0];
00243 for (Int_t i=0; i<fNdim-1; i++)
00244 ireal=fN[i+1]*ireal + ipoint[i+1];
00245
00246 const Double_t * array = GetPointsReal(fromInput);
00247 return ( array ) ? array[ireal] : 0;
00248 }
00249
00250
00251 void TFFTReal::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00252 {
00253
00254
00255 const Double_t * array = GetPointsReal(fromInput);
00256 if (!array) return;
00257 if ( ( ((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC && !fromInput ) ||
00258 ( ((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R && fromInput ) )
00259 {
00260 if (ipoint<fN[0]/2+1){
00261 re = array[ipoint];
00262 im = array[fN[0]-ipoint];
00263 } else {
00264 re = array[fN[0]-ipoint];
00265 im = -array[ipoint];
00266 }
00267 if ((fN[0]%2)==0 && ipoint==fN[0]/2) im = 0;
00268 }
00269 }
00270
00271 void TFFTReal::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00272 {
00273
00274
00275 GetPointComplex(ipoint[0], re, im, fromInput);
00276 }
00277
00278
00279 Double_t* TFFTReal::GetPointsReal(Bool_t fromInput) const
00280 {
00281
00282
00283
00284
00285
00286
00287
00288
00289 if (!fromInput && fOut)
00290 return (Double_t*)fOut;
00291 else if (fromInput && !fOut) {
00292 Error("GetPointsReal","Input array was destroyed");
00293 return 0;
00294 }
00295
00296 return (Double_t*)fIn;
00297 }
00298
00299
00300 void TFFTReal::SetPoint(Int_t ipoint, Double_t re, Double_t im)
00301 {
00302 if (ipoint<0 || ipoint>fTotalSize){
00303 Error("SetPoint", "illegal point index");
00304 return;
00305 }
00306 if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R){
00307 if ((fN[0]%2)==0 && ipoint==fN[0]/2)
00308 ((Double_t*)fIn)[ipoint] = re;
00309 else {
00310 ((Double_t*)fIn)[ipoint] = re;
00311 ((Double_t*)fIn)[fN[0]-ipoint]=im;
00312 }
00313 }
00314 else
00315 ((Double_t*)fIn)[ipoint]=re;
00316 }
00317
00318
00319 void TFFTReal::SetPoint(const Int_t *ipoint, Double_t re, Double_t )
00320 {
00321
00322
00323
00324 Int_t ireal = ipoint[0];
00325 for (Int_t i=0; i<fNdim-1; i++)
00326 ireal=fN[i+1]*ireal + ipoint[i+1];
00327 if (ireal < 0 || ireal >fTotalSize){
00328 Error("SetPoint", "illegal point index");
00329 return;
00330 }
00331 ((Double_t*)fIn)[ireal]=re;
00332 }
00333
00334
00335 void TFFTReal::SetPoints(const Double_t *data)
00336 {
00337
00338
00339 for (Int_t i=0; i<fTotalSize; i++)
00340 ((Double_t*)fIn)[i] = data[i];
00341 }
00342
00343
00344 Int_t TFFTReal::MapOptions(const Int_t *kind)
00345 {
00346
00347
00348 if (kind[0] == 10){
00349 if (fNdim>1){
00350 Error("Init", "Multidimensional R2HC transforms are not supported, use R2C interface instead");
00351 return 0;
00352 }
00353 ((fftw_r2r_kind*)fKind)[0] = FFTW_R2HC;
00354 }
00355 else if (kind[0] == 11) {
00356 if (fNdim>1){
00357 Error("Init", "Multidimensional HC2R transforms are not supported, use C2R interface instead");
00358 return 0;
00359 }
00360 ((fftw_r2r_kind*)fKind)[0] = FFTW_HC2R;
00361 }
00362 else if (kind[0] == 12) {
00363 for (Int_t i=0; i<fNdim; i++)
00364 ((fftw_r2r_kind*)fKind)[i] = FFTW_DHT;
00365 }
00366 else {
00367 for (Int_t i=0; i<fNdim; i++){
00368 switch (kind[i]) {
00369 case 0: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT00; break;
00370 case 1: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT01; break;
00371 case 2: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT10; break;
00372 case 3: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT11; break;
00373 case 4: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT00; break;
00374 case 5: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT01; break;
00375 case 6: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT10; break;
00376 case 7: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT11; break;
00377 default:
00378 ((fftw_r2r_kind*)fKind)[i] = FFTW_R2HC; break;
00379 }
00380 }
00381 }
00382 return 1;
00383 }
00384
00385
00386 UInt_t TFFTReal::MapFlag(Option_t *flag)
00387 {
00388
00389
00390
00391
00392
00393
00394 TString opt = flag;
00395 opt.ToUpper();
00396 if (opt.Contains("ES"))
00397 return FFTW_ESTIMATE;
00398 if (opt.Contains("M"))
00399 return FFTW_MEASURE;
00400 if (opt.Contains("P"))
00401 return FFTW_PATIENT;
00402 if (opt.Contains("EX"))
00403 return FFTW_EXHAUSTIVE;
00404 return FFTW_ESTIMATE;
00405 }