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 #include "TFFTComplex.h"
00040 #include "fftw3.h"
00041 #include "TComplex.h"
00042
00043
00044 ClassImp(TFFTComplex)
00045
00046
00047 TFFTComplex::TFFTComplex()
00048 {
00049
00050
00051 fIn = 0;
00052 fOut = 0;
00053 fPlan = 0;
00054 fN = 0;
00055 fFlags = 0;
00056 fNdim = 0;
00057 fTotalSize = 0;
00058 fSign = 1;
00059 }
00060
00061
00062 TFFTComplex::TFFTComplex(Int_t n, Bool_t inPlace)
00063 {
00064
00065
00066
00067 fIn = fftw_malloc(sizeof(fftw_complex) *n);
00068 if (!inPlace)
00069 fOut = fftw_malloc(sizeof(fftw_complex) * n);
00070 else
00071 fOut = 0;
00072 fN = new Int_t[1];
00073 fN[0] = n;
00074 fTotalSize = n;
00075 fNdim = 1;
00076 fSign = 1;
00077 fPlan = 0;
00078 fFlags = 0;
00079 }
00080
00081
00082 TFFTComplex::TFFTComplex(Int_t ndim, Int_t *n, Bool_t inPlace)
00083 {
00084
00085
00086
00087 fNdim = ndim;
00088 fTotalSize = 1;
00089 fN = new Int_t[fNdim];
00090 for (Int_t i=0; i<fNdim; i++){
00091 fN[i] = n[i];
00092 fTotalSize*=n[i];
00093 }
00094 fIn = fftw_malloc(sizeof(fftw_complex)*fTotalSize);
00095 if (!inPlace)
00096 fOut = fftw_malloc(sizeof(fftw_complex) * fTotalSize);
00097 else
00098 fOut = 0;
00099 fSign = 1;
00100 fPlan = 0;
00101 fFlags = 0;
00102 }
00103
00104
00105 TFFTComplex::~TFFTComplex()
00106 {
00107
00108
00109
00110
00111 fftw_destroy_plan((fftw_plan)fPlan);
00112 fPlan = 0;
00113 fftw_free((fftw_complex*)fIn);
00114 if (fOut)
00115 fftw_free((fftw_complex*)fOut);
00116 if (fN)
00117 delete [] fN;
00118 }
00119
00120
00121 void TFFTComplex::Init( Option_t *flags, Int_t sign,const Int_t* )
00122 {
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140 fSign = sign;
00141 fFlags = flags;
00142 if (fOut)
00143 fPlan = (void*)fftw_plan_dft(fNdim, fN, (fftw_complex*)fIn, (fftw_complex*)fOut, sign,MapFlag(flags));
00144 else
00145 fPlan = (void*)fftw_plan_dft(fNdim, fN, (fftw_complex*)fIn, (fftw_complex*)fIn, sign, MapFlag(flags));
00146 }
00147
00148
00149 void TFFTComplex::Transform()
00150 {
00151
00152
00153 if (fPlan)
00154 fftw_execute((fftw_plan)fPlan);
00155 else {
00156 Error("Transform", "transform not initialised");
00157 return;
00158 }
00159 }
00160
00161
00162 void TFFTComplex::GetPoints(Double_t *data, Bool_t fromInput) const
00163 {
00164
00165
00166 if (!fromInput){
00167 for (Int_t i=0; i<2*fTotalSize; i+=2){
00168 data[i] = ((fftw_complex*)fOut)[i/2][0];
00169 data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00170 }
00171 } else {
00172 for (Int_t i=0; i<2*fTotalSize; i+=2){
00173 data[i] = ((fftw_complex*)fIn)[i/2][0];
00174 data[i+1] = ((fftw_complex*)fIn)[i/2][1];
00175 }
00176 }
00177 }
00178
00179
00180 void TFFTComplex::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00181 {
00182
00183
00184 if (fOut && !fromInput){
00185 re = ((fftw_complex*)fOut)[ipoint][0];
00186 im = ((fftw_complex*)fOut)[ipoint][1];
00187 } else {
00188 re = ((fftw_complex*)fIn)[ipoint][0];
00189 im = ((fftw_complex*)fIn)[ipoint][1];
00190 }
00191 }
00192
00193
00194 void TFFTComplex::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00195 {
00196
00197
00198 Int_t ireal = ipoint[0];
00199 for (Int_t i=0; i<fNdim-1; i++)
00200 ireal=fN[i+1]*ireal + ipoint[i+1];
00201
00202 if (fOut && !fromInput){
00203 re = ((fftw_complex*)fOut)[ireal][0];
00204 im = ((fftw_complex*)fOut)[ireal][1];
00205 } else {
00206 re = ((fftw_complex*)fIn)[ireal][0];
00207 im = ((fftw_complex*)fIn)[ireal][1];
00208 }
00209 }
00210
00211
00212 void TFFTComplex::GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput) const
00213 {
00214
00215
00216 if (fOut && !fromInput){
00217 for (Int_t i=0; i<fTotalSize; i++){
00218 re[i] = ((fftw_complex*)fOut)[i][0];
00219 im[i] = ((fftw_complex*)fOut)[i][1];
00220 }
00221 } else {
00222 for (Int_t i=0; i<fTotalSize; i++){
00223 re[i] = ((fftw_complex*)fIn)[i][0];
00224 im[i] = ((fftw_complex*)fIn)[i][1];
00225 }
00226 }
00227 }
00228
00229
00230 void TFFTComplex::GetPointsComplex(Double_t *data, Bool_t fromInput) const
00231 {
00232
00233
00234 if (fOut && !fromInput){
00235 for (Int_t i=0; i<fTotalSize; i+=2){
00236 data[i] = ((fftw_complex*)fOut)[i/2][0];
00237 data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00238 }
00239 } else {
00240 for (Int_t i=0; i<fTotalSize; i+=2){
00241 data[i] = ((fftw_complex*)fIn)[i/2][0];
00242 data[i+1] = ((fftw_complex*)fIn)[i/2][1];
00243 }
00244 }
00245 }
00246
00247
00248 void TFFTComplex::SetPoint(Int_t ipoint, Double_t re, Double_t im)
00249 {
00250
00251
00252 ((fftw_complex*)fIn)[ipoint][0]=re;
00253 ((fftw_complex*)fIn)[ipoint][1]=im;
00254 }
00255
00256
00257 void TFFTComplex::SetPoint(const Int_t *ipoint, Double_t re, Double_t im)
00258 {
00259
00260
00261 Int_t ireal = ipoint[0];
00262 for (Int_t i=0; i<fNdim-1; i++)
00263 ireal=fN[i+1]*ireal + ipoint[i+1];
00264
00265 ((fftw_complex*)fIn)[ireal][0]=re;
00266 ((fftw_complex*)fIn)[ireal][1]=im;
00267 }
00268
00269
00270 void TFFTComplex::SetPointComplex(Int_t ipoint, TComplex &c)
00271 {
00272 ((fftw_complex*)fIn)[ipoint][0] = c.Re();
00273 ((fftw_complex*)fIn)[ipoint][1] = c.Im();
00274 }
00275
00276
00277 void TFFTComplex::SetPoints(const Double_t *data)
00278 {
00279
00280
00281
00282 for (Int_t i=0; i<2*fTotalSize-1; i+=2){
00283 ((fftw_complex*)fIn)[i/2][0]=data[i];
00284 ((fftw_complex*)fIn)[i/2][1]=data[i+1];
00285 }
00286 }
00287
00288
00289 void TFFTComplex::SetPointsComplex(const Double_t *re_data, const Double_t *im_data)
00290 {
00291
00292
00293 if (!fIn){
00294 Error("SetPointsComplex", "Size is not set yet");
00295 return;
00296 }
00297 for (Int_t i=0; i<fTotalSize; i++){
00298 ((fftw_complex*)fIn)[i][0]=re_data[i];
00299 ((fftw_complex*)fIn)[i][1]=im_data[i];
00300 }
00301 }
00302
00303
00304 UInt_t TFFTComplex::MapFlag(Option_t *flag)
00305 {
00306
00307
00308
00309
00310
00311
00312 TString opt = flag;
00313 opt.ToUpper();
00314 if (opt.Contains("ES"))
00315 return FFTW_ESTIMATE;
00316 if (opt.Contains("M"))
00317 return FFTW_MEASURE;
00318 if (opt.Contains("P"))
00319 return FFTW_PATIENT;
00320 if (opt.Contains("EX"))
00321 return FFTW_EXHAUSTIVE;
00322 return FFTW_ESTIMATE;
00323 }