musrfit  1.9.2
PFourier.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 
3  PFourier.cpp
4 
5  Author: Andreas Suter
6  e-mail: andreas.suter@psi.ch
7 
8 ***************************************************************************/
9 
10 /***************************************************************************
11  * Copyright (C) 2007-2023 by Andreas Suter *
12  * andreas.suter@psi.ch *
13  * *
14  * This program is free software; you can redistribute it and/or modify *
15  * it under the terms of the GNU General Public License as published by *
16  * the Free Software Foundation; either version 2 of the License, or *
17  * (at your option) any later version. *
18  * *
19  * This program is distributed in the hope that it will be useful, *
20  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
21  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
22  * GNU General Public License for more details. *
23  * *
24  * You should have received a copy of the GNU General Public License *
25  * along with this program; if not, write to the *
26  * Free Software Foundation, Inc., *
27  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
28  ***************************************************************************/
29 
30 #include <cmath>
31 
32 #include <iostream>
33 #include <iomanip>
34 
35 #include "TF1.h"
36 #include "TAxis.h"
37 
38 #include "Minuit2/FunctionMinimum.h"
39 #include "Minuit2/MnUserParameters.h"
40 #include "Minuit2/MnMinimize.h"
41 
42 #include "PMusr.h"
43 #include "PFourier.h"
44 
45 #define PI 3.14159265358979312
46 #define PI_HALF 1.57079632679489656
47 
48 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
49 // PFTPhaseCorrection
50 //--------------------------------------------------------------------------
51 // Constructor
52 //--------------------------------------------------------------------------
56 PFTPhaseCorrection::PFTPhaseCorrection(const Int_t minBin, const Int_t maxBin) :
57  fMinBin(minBin), fMaxBin(maxBin)
58 {
59  Init();
60 }
61 
62 //--------------------------------------------------------------------------
63 // Constructor
64 //--------------------------------------------------------------------------
68 PFTPhaseCorrection::PFTPhaseCorrection(std::vector<Double_t> &reFT, std::vector<Double_t> &imFT, const Int_t minBin, const Int_t maxBin) :
69  fReal(reFT), fImag(imFT), fMinBin(minBin), fMaxBin(maxBin)
70 {
71  Init();
72 
73  Int_t realSize = static_cast<Int_t>(fReal.size());
74  if (fMinBin == -1)
75  fMinBin = 0;
76  if (fMaxBin == -1)
77  fMaxBin = realSize;
78  if (fMaxBin > realSize)
79  fMaxBin = realSize;
80 
81  fRealPh.resize(fReal.size());
82  fImagPh.resize(fReal.size());
83 }
84 
85 //--------------------------------------------------------------------------
86 // Minimize (public)
87 //--------------------------------------------------------------------------
92 {
93  // create Minuit2 parameters
94  ROOT::Minuit2::MnUserParameters upar;
95 
96  upar.Add("c0", fPh_c0, 2.0);
97  upar.Add("c1", fPh_c1, 2.0);
98 
99  // create minimizer
100  ROOT::Minuit2::MnMinimize mn_min(*this, upar);
101 
102  // minimize
103  ROOT::Minuit2::FunctionMinimum min = mn_min();
104 
105  if (min.IsValid()) {
106  fPh_c0 = min.UserState().Value("c0");
107  fPh_c1 = min.UserState().Value("c1");
108  fMin = min.Fval();
109  } else {
110  fMin = -1.0;
111  fValid = false;
112  std::cout << std::endl << ">> **WARNING** minimize failed to find a minimum for the real FT phase correction ..." << std::endl;
113  return;
114  }
115 }
116 
117 //--------------------------------------------------------------------------
118 // GetPhaseCorrectionParam (public)
119 //--------------------------------------------------------------------------
124 {
125  Double_t result=0.0;
126 
127  if (idx == 0)
128  result = fPh_c0;
129  else if (idx == 1)
130  result = fPh_c1;
131  else
132  std::cerr << ">> **ERROR** requested phase correction parameter with index=" << idx << " does not exist!" << std::endl;
133 
134  return result;
135 }
136 
137 //--------------------------------------------------------------------------
138 // GetMinimum (public)
139 //--------------------------------------------------------------------------
144 {
145  if (!fValid) {
146  std::cerr << ">> **ERROR** requested minimum is invalid!" << std::endl;
147  return -1.0;
148  }
149 
150  return fMin;
151 }
152 
153 //--------------------------------------------------------------------------
154 // Init (private)
155 //--------------------------------------------------------------------------
160 {
161  fValid = true;
162  fPh_c0 = 0.0;
163  fPh_c1 = 0.0;
164  fGamma = 1.0;
165  fMin = -1.0;
166 }
167 
168 //--------------------------------------------------------------------------
169 // CalcPhasedFT (private)
170 //--------------------------------------------------------------------------
175 {
176  Double_t phi=0.0;
177  Double_t w=0.0;
178 
179  for (UInt_t i=0; i<fRealPh.size(); i++) {
180  w = static_cast<Double_t>(i) / static_cast<Double_t>(fReal.size());
181  phi = fPh_c0 + fPh_c1 * w;
182  fRealPh[i] = fReal[i]*cos(phi) - fImag[i]*sin(phi);
183  fImagPh[i] = fReal[i]*sin(phi) + fImag[i]*cos(phi);
184  }
185 }
186 
187 //--------------------------------------------------------------------------
188 // CalcRealPhFTDerivative (private)
189 //--------------------------------------------------------------------------
194 {
195  fRealPhD.resize(fRealPh.size());
196 
197  fRealPhD[0] = 1.0;
198  fRealPhD[fRealPh.size()-1] = 1.0;
199 
200  for (UInt_t i=1; i<fRealPh.size()-1; i++)
201  fRealPhD[i] = fRealPh[i+1]-fRealPh[i];
202 }
203 
204 //--------------------------------------------------------------------------
205 // Penalty (private)
206 //--------------------------------------------------------------------------
211 {
212  Double_t penalty = 0.0;
213 
214  for (UInt_t i=fMinBin; i<fMaxBin; i++) {
215  if (fRealPh[i] < 0.0)
216  penalty += fRealPh[i]*fRealPh[i];
217  }
218 
219  return fGamma*penalty;
220 }
221 
222 //--------------------------------------------------------------------------
223 // Entropy (private)
224 //--------------------------------------------------------------------------
229 {
230  Double_t norm = 0.0;
231 
232  for (UInt_t i=fMinBin; i<fMaxBin; i++)
233  norm += fabs(fRealPhD[i]);
234 
235  Double_t entropy = 0.0;
236  Double_t dval = 0.0, hh = 0.0;
237 
238  for (UInt_t i=fMinBin; i<fMaxBin; i++) {
239  dval = fabs(fRealPhD[i]);
240  if (dval > 1.0e-15) {
241  hh = dval / norm;
242  entropy -= hh * log(hh);
243  }
244  }
245 
246  return entropy;
247 }
248 
249 //--------------------------------------------------------------------------
250 // operator() (private)
251 //--------------------------------------------------------------------------
255 double PFTPhaseCorrection::operator()(const std::vector<double> &par) const
256 {
257  // par : [0]: c0, [1]: c1
258 
259  fPh_c0 = par[0];
260  fPh_c1 = par[1];
261 
262  CalcPhasedFT();
264 
265  return Entropy()+Penalty();
266 }
267 
268 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
269 // PFourier
270 //--------------------------------------------------------------------------
271 // Constructor
272 //--------------------------------------------------------------------------
284 PFourier::PFourier(TH1F *data, Int_t unitTag, Double_t startTime, Double_t endTime, Bool_t dcCorrected, UInt_t zeroPaddingPower) :
285  fData(data), fUnitTag(unitTag), fStartTime(startTime), fEndTime(endTime),
286  fDCCorrected(dcCorrected), fZeroPaddingPower(zeroPaddingPower)
287 {
288  // some necessary checks and initialization
289  if (fData == nullptr) {
290  std::cerr << std::endl << "**ERROR** PFourier::PFourier: no valid data" << std::endl << std::endl;
291  fValid = false;
292  return;
293  }
294 
295  fValid = true;
296  fIn = nullptr;
297  fOut = nullptr;
298 //as fPhCorrectedReFT = 0;
299 
301 
302  // calculate time resolution in (us)
303  fTimeResolution = fData->GetBinWidth(1);
304 
305  // if endTime == 0 set it to the last time slot
306  if (fEndTime == 0.0) {
307  Int_t last = fData->GetNbinsX()-1;
308  fEndTime = fData->GetBinCenter(last);
309  }
310 
311  // swap start and end time if necessary
312  if (fStartTime > fEndTime) {
313  Double_t keep = fStartTime;
315  fEndTime = keep;
316  }
317 
318  // calculate start and end bin
319  fNoOfData = static_cast<UInt_t>((fEndTime-fStartTime)/fTimeResolution);
320 
321  // check if zero padding is whished
322  if (fZeroPaddingPower > 0) {
323  UInt_t noOfBins = static_cast<UInt_t>(pow(2.0, static_cast<Double_t>(fZeroPaddingPower)));
324  if (noOfBins > fNoOfData)
325  fNoOfBins = noOfBins;
326  else
328  } else {
330  }
331 
332  // calculate fourier resolution, depending on the units
333  Double_t resolution = 1.0/(fTimeResolution*fNoOfBins); // in MHz
334  switch (fUnitTag) {
335  case FOURIER_UNIT_GAUSS:
336  fResolution = resolution/GAMMA_BAR_MUON;
337  break;
338  case FOURIER_UNIT_TESLA:
339  fResolution = 1e-4*resolution/GAMMA_BAR_MUON;
340  break;
341  case FOURIER_UNIT_FREQ:
342  fResolution = resolution;
343  break;
344  case FOURIER_UNIT_CYCLES:
345  fResolution = 2.0*PI*resolution;
346  break;
347  default:
348  fValid = false;
349  return;
350  }
351 
352  // allocate necessary memory
353  fIn = static_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex)*fNoOfBins));
354  fOut = static_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex)*fNoOfBins));
355 
356  // check if memory allocation has been successful
357  if ((fIn == nullptr) || (fOut == nullptr)) {
358  fValid = false;
359  return;
360  }
361 
362  // get the FFTW3 plan (see FFTW3 manual)
363  fFFTwPlan = fftw_plan_dft_1d(static_cast<Int_t>(fNoOfBins), fIn, fOut, FFTW_FORWARD, FFTW_ESTIMATE);
364 
365  // check if a valid plan has been generated
366  if (!fFFTwPlan) {
367  fValid = false;
368  }
369 }
370 
371 //--------------------------------------------------------------------------
372 // Destructor
373 //--------------------------------------------------------------------------
378 {
379  if (fFFTwPlan)
380  fftw_destroy_plan(fFFTwPlan);
381  if (fIn)
382  fftw_free(fIn);
383  if (fOut)
384  fftw_free(fOut);
385 //as if (fPhCorrectedReFT)
386 //as delete fPhCorrectedReFT;
387 }
388 
389 //--------------------------------------------------------------------------
390 // Transform (public)
391 //--------------------------------------------------------------------------
398 void PFourier::Transform(UInt_t apodizationTag)
399 {
400  if (!fValid)
401  return;
402 
403  PrepareFFTwInputData(apodizationTag);
404 
405  fftw_execute(fFFTwPlan);
406 
407  // correct the phase for tstart != 0.0
408  // find the first bin >= fStartTime
409  Double_t shiftTime = 0.0;
410  for (Int_t i=1; i<fData->GetXaxis()->GetNbins(); i++) {
411  if (fData->GetXaxis()->GetBinCenter(i) >= fStartTime) {
412  shiftTime = fData->GetXaxis()->GetBinCenter(i);
413  break;
414  }
415  }
416 
417  Double_t phase, re, im;
418  for (UInt_t i=0; i<fNoOfBins; i++) {
419  phase = 2.0*PI/(fTimeResolution*fNoOfBins) * i * shiftTime;
420  re = fOut[i][0] * cos(phase) + fOut[i][1] * sin(phase);
421  im = -fOut[i][0] * sin(phase) + fOut[i][1] * cos(phase);
422  fOut[i][0] = re;
423  fOut[i][1] = im;
424  }
425 }
426 
427 //--------------------------------------------------------------------------
428 // GetMaxFreq (public)
429 //--------------------------------------------------------------------------
434 {
435  UInt_t noOfFourierBins = 0;
436  if (fNoOfBins % 2 == 0)
437  noOfFourierBins = fNoOfBins/2;
438  else
439  noOfFourierBins = (fNoOfBins+1)/2;
440 
441  return fResolution*noOfFourierBins;
442 }
443 
444 //--------------------------------------------------------------------------
445 // GetRealFourier (public)
446 //--------------------------------------------------------------------------
452 TH1F* PFourier::GetRealFourier(const Double_t scale)
453 {
454  // check if valid flag is set
455  if (!fValid)
456  return nullptr;
457 
458  // invoke realFourier
459  Char_t name[256];
460  Char_t title[256];
461  snprintf(name, sizeof(name), "%s_Fourier_Re", fData->GetName());
462  snprintf(title, sizeof(title), "%s_Fourier_Re", fData->GetTitle());
463 
464  UInt_t noOfFourierBins = 0;
465  if (fNoOfBins % 2 == 0)
466  noOfFourierBins = fNoOfBins/2;
467  else
468  noOfFourierBins = (fNoOfBins+1)/2;
469 
470  TH1F *realFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
471  if (realFourier == nullptr) {
472  fValid = false;
473  std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
474  return nullptr;
475  }
476 
477  // fill realFourier vector
478  for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
479  realFourier->SetBinContent(i+1, scale*fOut[i][0]);
480  realFourier->SetBinError(i+1, 0.0);
481  }
482  return realFourier;
483 }
484 
485 //--------------------------------------------------------------------------
486 // GetPhaseOptRealFourier (public, static)
487 //--------------------------------------------------------------------------
500 TH1F* PFourier::GetPhaseOptRealFourier(const TH1F *re, const TH1F *im, std::vector<Double_t> &phase,
501  const Double_t scale, const Double_t min, const Double_t max)
502 {
503  if ((re == nullptr) || (im == nullptr))
504  return nullptr;
505 
506  phase.resize(2); // c0, c1
507 
508  const TAxis *axis = re->GetXaxis();
509 
510  Int_t minBin = 1;
511  Int_t maxBin = axis->GetNbins();
512  Int_t noOfBins = axis->GetNbins();
513  Double_t res = axis->GetBinWidth(1);
514 
515  // check if minimum frequency is given. If yes, get the proper minBin
516  if (min != -1.0) {
517  minBin = axis->FindFixBin(min);
518  if ((minBin == 0) || (minBin > maxBin)) {
519  minBin = 1;
520  std::cerr << "**WARNING** minimum frequency/field out of range. Will adopted it." << std::endl;
521  }
522  }
523 
524  // check if maximum frequency is given. If yes, get the proper maxBin
525  if (max != -1.0) {
526  maxBin = axis->FindFixBin(max);
527  if ((maxBin == 0) || (maxBin > axis->GetNbins())) {
528  maxBin = axis->GetNbins();
529  std::cerr << "**WARNING** maximum frequency/field out of range. Will adopted it." << std::endl;
530  }
531  }
532 
533  // copy the real/imag Fourier from min to max
534  std::vector<Double_t> realF, imagF;
535  for (Int_t i=minBin; i<=maxBin; i++) {
536  realF.push_back(re->GetBinContent(i));
537  imagF.push_back(im->GetBinContent(i));
538  }
539 
540  // optimize real FT phase
541  PFTPhaseCorrection *phCorrectedReFT = new PFTPhaseCorrection(realF, imagF);
542  if (phCorrectedReFT == nullptr) {
543  std::cerr << std::endl << "**SEVERE ERROR** couldn't invoke PFTPhaseCorrection object." << std::endl;
544  return nullptr;
545  }
546 
547  phCorrectedReFT->Minimize();
548  if (!phCorrectedReFT->IsValid()) {
549  std::cerr << std::endl << "**ERROR** could not find a valid phase correction minimum." << std::endl;
550  return nullptr;
551  }
552  phase[0] = phCorrectedReFT->GetPhaseCorrectionParam(0);
553  phase[1] = phCorrectedReFT->GetPhaseCorrectionParam(1);
554 
555  // clean up
556  if (phCorrectedReFT) {
557  delete phCorrectedReFT;
558  phCorrectedReFT = nullptr;
559  }
560  realF.clear();
561  imagF.clear();
562 
563  // invoke the real phase optimised histo to be filled. Caller is the owner!
564  Char_t name[256];
565  Char_t title[256];
566  snprintf(name, sizeof(name), "%s_Fourier_PhOptRe", re->GetName());
567  snprintf(title, sizeof(title), "%s_Fourier_PhOptRe", re->GetTitle());
568 
569  TH1F *realPhaseOptFourier = new TH1F(name, title, noOfBins, -res/2.0, static_cast<Double_t>(noOfBins-1)*res+res/2.0);
570  if (realPhaseOptFourier == nullptr) {
571  std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
572  return nullptr;
573  }
574 
575  // fill realFourier vector
576  Double_t ph;
577  for (Int_t i=0; i<noOfBins; i++) {
578  ph = phase[0] + phase[1] * static_cast<Double_t>(i-static_cast<Int_t>(minBin)) / static_cast<Double_t>(maxBin-minBin);
579  realPhaseOptFourier->SetBinContent(i+1, scale*(re->GetBinContent(i+1)*cos(ph) - im->GetBinContent(i+1)*sin(ph)));
580  realPhaseOptFourier->SetBinError(i+1, 0.0);
581  }
582 
583  return realPhaseOptFourier;
584 }
585 
586 //--------------------------------------------------------------------------
587 // GetImaginaryFourier (public)
588 //--------------------------------------------------------------------------
594 TH1F* PFourier::GetImaginaryFourier(const Double_t scale)
595 {
596  // check if valid flag is set
597  if (!fValid)
598  return nullptr;
599 
600  // invoke imaginaryFourier
601  Char_t name[256];
602  Char_t title[256];
603  snprintf(name, sizeof(name), "%s_Fourier_Im", fData->GetName());
604  snprintf(title, sizeof(title), "%s_Fourier_Im", fData->GetTitle());
605 
606  UInt_t noOfFourierBins = 0;
607  if (fNoOfBins % 2 == 0)
608  noOfFourierBins = fNoOfBins/2;
609  else
610  noOfFourierBins = (fNoOfBins+1)/2;
611 
612  TH1F* imaginaryFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
613  if (imaginaryFourier == nullptr) {
614  fValid = false;
615  std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the imaginary part of the Fourier transform." << std::endl;
616  return nullptr;
617  }
618 
619  // fill imaginaryFourier vector
620  for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
621  imaginaryFourier->SetBinContent(i+1, scale*fOut[i][1]);
622  imaginaryFourier->SetBinError(i+1, 0.0);
623  }
624 
625  return imaginaryFourier;
626 }
627 
628 //--------------------------------------------------------------------------
629 // GetPowerFourier (public)
630 //--------------------------------------------------------------------------
636 TH1F* PFourier::GetPowerFourier(const Double_t scale)
637 {
638  // check if valid flag is set
639  if (!fValid)
640  return nullptr;
641 
642  // invoke powerFourier
643  Char_t name[256];
644  Char_t title[256];
645  snprintf(name, sizeof(name), "%s_Fourier_Pwr", fData->GetName());
646  snprintf(title, sizeof(title), "%s_Fourier_Pwr", fData->GetTitle());
647 
648  UInt_t noOfFourierBins = 0;
649  if (fNoOfBins % 2 == 0)
650  noOfFourierBins = fNoOfBins/2;
651  else
652  noOfFourierBins = (fNoOfBins+1)/2;
653 
654  TH1F* pwrFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
655  if (pwrFourier == nullptr) {
656  fValid = false;
657  std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the power part of the Fourier transform." << std::endl;
658  return nullptr;
659  }
660 
661  // fill powerFourier vector
662  for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
663  pwrFourier->SetBinContent(i+1, scale*sqrt(fOut[i][0]*fOut[i][0]+fOut[i][1]*fOut[i][1]));
664  pwrFourier->SetBinError(i+1, 0.0);
665  }
666 
667  return pwrFourier;
668 }
669 
670 //--------------------------------------------------------------------------
671 // GetPhaseFourier (public)
672 //--------------------------------------------------------------------------
678 TH1F* PFourier::GetPhaseFourier(const Double_t scale)
679 {
680  // check if valid flag is set
681  if (!fValid)
682  return nullptr;
683 
684  // invoke phaseFourier
685  Char_t name[256];
686  Char_t title[256];
687  snprintf(name, sizeof(name), "%s_Fourier_Phase", fData->GetName());
688  snprintf(title, sizeof(title), "%s_Fourier_Phase", fData->GetTitle());
689 
690  UInt_t noOfFourierBins = 0;
691  if (fNoOfBins % 2 == 0)
692  noOfFourierBins = fNoOfBins/2;
693  else
694  noOfFourierBins = (fNoOfBins+1)/2;
695 
696  TH1F* phaseFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
697  if (phaseFourier == nullptr) {
698  fValid = false;
699  std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the phase part of the Fourier transform." << std::endl;
700  return nullptr;
701  }
702 
703  // fill phaseFourier vector
704  Double_t value = 0.0;
705  for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
706  // calculate the phase
707  if (fOut[i][0] == 0.0) {
708  if (fOut[i][1] >= 0.0)
709  value = PI_HALF;
710  else
711  value = -PI_HALF;
712  } else {
713  value = atan(fOut[i][1]/fOut[i][0]);
714  // check sector
715  if (fOut[i][0] < 0.0) {
716  if (fOut[i][1] > 0.0)
717  value = PI + value;
718  else
719  value = PI - value;
720  }
721  }
722  phaseFourier->SetBinContent(i+1, scale*value);
723  phaseFourier->SetBinError(i+1, 0.0);
724  }
725 
726  return phaseFourier;
727 }
728 
729 //--------------------------------------------------------------------------
730 // PrepareFFTwInputData (private)
731 //--------------------------------------------------------------------------
739 void PFourier::PrepareFFTwInputData(UInt_t apodizationTag)
740 {
741  // 1st find t==0. fData start at times t<0!!
742  Int_t t0bin = -1;
743  for (Int_t i=1; i<fData->GetNbinsX(); i++) {
744  if (fData->GetBinCenter(i) >= 0.0) {
745  t0bin = i;
746  break;
747  }
748  }
749 
750  Int_t ival = static_cast<Int_t>(fStartTime/fTimeResolution) + t0bin;
751  UInt_t start = 0;
752  if (ival >= 0) {
753  start = static_cast<UInt_t>(ival);
754  }
755 
756  Double_t mean = 0.0;
757  if (fDCCorrected) {
758  for (UInt_t i=0; i<fNoOfData; i++) {
759  mean += fData->GetBinContent(static_cast<Int_t>(i+start));
760  }
761  mean /= static_cast<Double_t>(fNoOfData);
762  }
763 
764  // 2nd fill fIn
765  for (UInt_t i=0; i<fNoOfData; i++) {
766  fIn[i][0] = fData->GetBinContent(static_cast<Int_t>(i+start)) - mean;
767  fIn[i][1] = 0.0;
768  }
769  for (UInt_t i=fNoOfData; i<fNoOfBins; i++) {
770  fIn[i][0] = 0.0;
771  fIn[i][1] = 0.0;
772  }
773 
774  // 3rd apodize data (if wished)
775  ApodizeData(static_cast<Int_t>(apodizationTag));
776 }
777 
778 //--------------------------------------------------------------------------
779 // ApodizeData (private)
780 //--------------------------------------------------------------------------
788 void PFourier::ApodizeData(Int_t apodizationTag) {
789 
790  const Double_t cweak[3] = { 0.384093, -0.087577, 0.703484 };
791  const Double_t cmedium[3] = { 0.152442, -0.136176, 0.983734 };
792  const Double_t cstrong[3] = { 0.045335, 0.554883, 0.399782 };
793 
794  Double_t c[5];
795  for (UInt_t i=0; i<5; i++) {
796  c[i] = 0.0;
797  }
798 
799  switch (apodizationTag) {
800  case F_APODIZATION_NONE:
801  return;
802  case F_APODIZATION_WEAK:
803  c[0] = cweak[0]+cweak[1]+cweak[2];
804  c[1] = -(cweak[1]+2.0*cweak[2]);
805  c[2] = cweak[2];
806  break;
808  c[0] = cmedium[0]+cmedium[1]+cmedium[2];
809  c[1] = -(cmedium[1]+2.0*cmedium[2]);
810  c[2] = cmedium[2];
811  break;
813  c[0] = cstrong[0]+cstrong[1]+cstrong[2];
814  c[1] = -2.0*(cstrong[1]+2.0*cstrong[2]);
815  c[2] = cstrong[1]+6.0*cstrong[2];
816  c[3] = -4.0*cstrong[2];
817  c[4] = cstrong[2];
818  break;
819  default:
820  std::cerr << std::endl << ">> **ERROR** User Apodization tag " << apodizationTag << " unknown, sorry ..." << std::endl;
821  break;
822  }
823 
824  Double_t q;
825  for (UInt_t i=0; i<fNoOfData; i++) {
826  q = c[0];
827  for (UInt_t j=1; j<5; j++) {
828  q += c[j] * pow(static_cast<Double_t>(i)/static_cast<Double_t>(fNoOfData), 2.0*static_cast<Double_t>(j));
829  }
830  fIn[i][0] *= q;
831  }
832 }
fftw_complex * fOut
imaginary part of the Fourier transform
Definition: PFourier.h:141
Double_t fPh_c0
maximum bin from Fourier range to be used for the phase correction estimate
Definition: PFourier.h:79
virtual void PrepareFFTwInputData(UInt_t apodizationTag)
Definition: PFourier.cpp:739
virtual void CalcRealPhFTDerivative() const
Definition: PFourier.cpp:193
Double_t fPh_c1
constant part of the phase dispersion used for the phase correction
Definition: PFourier.h:80
UInt_t fNoOfData
number of bins in the time interval between fStartTime and fStopTime
Definition: PFourier.h:137
virtual Double_t Entropy() const
Definition: PFourier.cpp:228
Double_t fGamma
linear part of the phase dispersion used for the phase correction
Definition: PFourier.h:81
virtual Bool_t IsValid()
Definition: PFourier.h:58
fftw_plan fFFTwPlan
fftw plan (see FFTW3 User Manual)
Definition: PFourier.h:139
#define F_APODIZATION_STRONG
Definition: PFourier.h:46
Int_t fApodization
0=none, 1=weak, 2=medium, 3=strong
Definition: PFourier.h:128
Int_t fMaxBin
minimum bin from Fourier range to be used for the phase correction estimate
Definition: PFourier.h:78
TH1F * fData
data histogram to be Fourier transformed.
Definition: PFourier.h:123
fftw_complex * fIn
real part of the Fourier transform
Definition: PFourier.h:140
PFourier(TH1F *data, Int_t unitTag, Double_t startTime=0.0, Double_t endTime=0.0, Bool_t dcCorrected=false, UInt_t zeroPaddingPower=0)
Definition: PFourier.cpp:284
Double_t fResolution
Fourier resolution (field, frequency, or angular frequency)
Definition: PFourier.h:135
std::vector< Double_t > fReal
Definition: PFourier.h:71
#define FOURIER_UNIT_CYCLES
Definition: PMusr.h:135
virtual Double_t GetMaxFreq()
Definition: PFourier.cpp:433
virtual void Transform(UInt_t apodizationTag=0)
Definition: PFourier.cpp:398
#define FOURIER_UNIT_FREQ
Definition: PMusr.h:134
PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1)
Definition: PFourier.cpp:56
virtual void CalcPhasedFT() const
Definition: PFourier.cpp:174
static TH1F * GetPhaseOptRealFourier(const TH1F *re, const TH1F *im, std::vector< Double_t > &phase, const Double_t scale=1.0, const Double_t min=-1.0, const Double_t max=-1.0)
Definition: PFourier.cpp:500
virtual TH1F * GetRealFourier(const Double_t scale=1.0)
Definition: PFourier.cpp:452
#define F_APODIZATION_NONE
Definition: PFourier.h:43
virtual Double_t operator()(const std::vector< Double_t > &) const
Definition: PFourier.cpp:255
Bool_t fDCCorrected
if true, removed DC offset from signal before Fourier transformation, otherwise not ...
Definition: PFourier.h:133
#define F_APODIZATION_WEAK
Definition: PFourier.h:44
std::vector< Double_t > fImag
original real Fourier data set
Definition: PFourier.h:72
Double_t fStartTime
start time of the data histogram
Definition: PFourier.h:131
std::vector< Double_t > fImagPh
phased real Fourier data set
Definition: PFourier.h:74
virtual void ApodizeData(Int_t apodizationTag)
Definition: PFourier.cpp:788
virtual TH1F * GetPowerFourier(const Double_t scale=1.0)
Definition: PFourier.cpp:636
virtual TH1F * GetPhaseFourier(const Double_t scale=1.0)
Definition: PFourier.cpp:678
virtual Double_t GetMinimum()
Definition: PFourier.cpp:143
#define FOURIER_UNIT_GAUSS
Definition: PMusr.h:132
virtual Double_t GetPhaseCorrectionParam(UInt_t idx)
Definition: PFourier.cpp:123
#define PI_HALF
Definition: PFourier.cpp:46
#define FOURIER_UNIT_TESLA
Definition: PMusr.h:133
Bool_t fValid
true = all boundary conditions fullfilled and hence a Fourier transform can be performed.
Definition: PFourier.h:125
#define PI
Definition: PFourier.cpp:45
Int_t fMinBin
1st derivative of fRealPh
Definition: PFourier.h:77
#define F_APODIZATION_MEDIUM
Definition: PFourier.h:45
virtual ~PFourier()
Definition: PFourier.cpp:377
std::vector< Double_t > fRealPhD
phased imag Fourier data set
Definition: PFourier.h:75
Int_t fUnitTag
1=Field Units (G), 2=Field Units (T), 3=Frequency Units (MHz), 4=Angular Frequency Units (Mc/s) ...
Definition: PFourier.h:126
virtual void Init()
keeps the minimum of the entropy/penalty minimization
Definition: PFourier.cpp:159
std::vector< Double_t > fRealPh
original imag Fourier data set
Definition: PFourier.h:73
virtual void Minimize()
Definition: PFourier.cpp:91
virtual Double_t Penalty() const
Definition: PFourier.cpp:210
Double_t fMin
gamma parameter to balance between entropy and penalty
Definition: PFourier.h:82
#define GAMMA_BAR_MUON
Definition: PMusr.h:78
virtual TH1F * GetImaginaryFourier(const Double_t scale=1.0)
Definition: PFourier.cpp:594
UInt_t fNoOfBins
number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding ...
Definition: PFourier.h:138
Double_t fTimeResolution
time resolution of the data histogram in (us)
Definition: PFourier.h:130
Double_t fEndTime
end time of the data histogram
Definition: PFourier.h:132
UInt_t fZeroPaddingPower
power for zero padding, if set < 0 no zero padding will be done
Definition: PFourier.h:134