musrfit 1.10.0
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-2026 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//--------------------------------------------------------------------------
67PFTPhaseCorrection::PFTPhaseCorrection(const Int_t minBin, const Int_t maxBin) :
68 fMinBin(minBin), fMaxBin(maxBin)
69{
70 Init();
71}
72
73//--------------------------------------------------------------------------
74// Constructor
75//--------------------------------------------------------------------------
100PFTPhaseCorrection::PFTPhaseCorrection(std::vector<Double_t> &reFT, std::vector<Double_t> &imFT, const Int_t minBin, const Int_t maxBin) :
101 fReal(reFT), fImag(imFT), fMinBin(minBin), fMaxBin(maxBin)
102{
103 Init();
104
105 Int_t realSize = static_cast<Int_t>(fReal.size());
106 if (fMinBin == -1)
107 fMinBin = 0;
108 if (fMaxBin == -1)
109 fMaxBin = realSize;
110 if (fMaxBin > realSize)
111 fMaxBin = realSize;
112
113 fRealPh.resize(fReal.size());
114 fImagPh.resize(fReal.size());
115}
116
117//--------------------------------------------------------------------------
118// Minimize (public)
119//--------------------------------------------------------------------------
145{
146 // create Minuit2 parameters
147 ROOT::Minuit2::MnUserParameters upar;
148
149 upar.Add("c0", fPh_c0, 2.0);
150 upar.Add("c1", fPh_c1, 2.0);
151
152 // create minimizer
153 ROOT::Minuit2::MnMinimize mn_min(*this, upar);
154
155 // minimize
156 ROOT::Minuit2::FunctionMinimum min = mn_min();
157
158 if (min.IsValid()) {
159 fPh_c0 = min.UserState().Value("c0");
160 fPh_c1 = min.UserState().Value("c1");
161 fMin = min.Fval();
162 } else {
163 fMin = -1.0;
164 fValid = false;
165 std::cout << std::endl << ">> **WARNING** minimize failed to find a minimum for the real FT phase correction ..." << std::endl;
166 return;
167 }
168}
169
170//--------------------------------------------------------------------------
171// GetPhaseCorrectionParam (public)
172//--------------------------------------------------------------------------
196{
197 Double_t result=0.0;
198
199 if (idx == 0)
200 result = fPh_c0;
201 else if (idx == 1)
202 result = fPh_c1;
203 else
204 std::cerr << ">> **ERROR** requested phase correction parameter with index=" << idx << " does not exist!" << std::endl;
205
206 return result;
207}
208
209//--------------------------------------------------------------------------
210// GetMinimum (public)
211//--------------------------------------------------------------------------
224{
225 if (!fValid) {
226 std::cerr << ">> **ERROR** requested minimum is invalid!" << std::endl;
227 return -1.0;
228 }
229
230 return fMin;
231}
232
233//--------------------------------------------------------------------------
234// Init (private)
235//--------------------------------------------------------------------------
249{
250 fValid = true;
251 fPh_c0 = 0.0;
252 fPh_c1 = 0.0;
253 fGamma = 1.0;
254 fMin = -1.0;
255}
256
257//--------------------------------------------------------------------------
258// CalcPhasedFT (private)
259//--------------------------------------------------------------------------
277{
278 Double_t phi=0.0;
279 Double_t w=0.0;
280
281 for (UInt_t i=0; i<fRealPh.size(); i++) {
282 w = static_cast<Double_t>(i) / static_cast<Double_t>(fReal.size());
283 phi = fPh_c0 + fPh_c1 * w;
284 fRealPh[i] = fReal[i]*cos(phi) - fImag[i]*sin(phi);
285 fImagPh[i] = fReal[i]*sin(phi) + fImag[i]*cos(phi);
286 }
287}
288
289//--------------------------------------------------------------------------
290// CalcRealPhFTDerivative (private)
291//--------------------------------------------------------------------------
308{
309 fRealPhD.resize(fRealPh.size());
310
311 fRealPhD[0] = 1.0;
312 fRealPhD[fRealPh.size()-1] = 1.0;
313
314 for (UInt_t i=1; i<fRealPh.size()-1; i++)
315 fRealPhD[i] = fRealPh[i+1]-fRealPh[i];
316}
317
318//--------------------------------------------------------------------------
319// Penalty (private)
320//--------------------------------------------------------------------------
340{
341 Double_t penalty = 0.0;
342
343 for (UInt_t i=fMinBin; i<fMaxBin; i++) {
344 if (fRealPh[i] < 0.0)
345 penalty += fRealPh[i]*fRealPh[i];
346 }
347
348 return fGamma*penalty;
349}
350
351//--------------------------------------------------------------------------
352// Entropy (private)
353//--------------------------------------------------------------------------
377{
378 Double_t norm = 0.0;
379
380 for (UInt_t i=fMinBin; i<fMaxBin; i++)
381 norm += fabs(fRealPhD[i]);
382
383 Double_t entropy = 0.0;
384 Double_t dval = 0.0, hh = 0.0;
385
386 for (UInt_t i=fMinBin; i<fMaxBin; i++) {
387 dval = fabs(fRealPhD[i]);
388 if (dval > 1.0e-15) {
389 hh = dval / norm;
390 entropy -= hh * log(hh);
391 }
392 }
393
394 return entropy;
395}
396
397//--------------------------------------------------------------------------
398// operator() (private)
399//--------------------------------------------------------------------------
425double PFTPhaseCorrection::operator()(const std::vector<double> &par) const
426{
427 // par : [0]: c0, [1]: c1
428
429 fPh_c0 = par[0];
430 fPh_c1 = par[1];
431
432 CalcPhasedFT();
434
435 return Entropy()+Penalty();
436}
437
438//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
439// PFourier
440//--------------------------------------------------------------------------
441// Constructor
442//--------------------------------------------------------------------------
482PFourier::PFourier(TH1F *data, Int_t unitTag, Double_t startTime, Double_t endTime, Bool_t dcCorrected, UInt_t zeroPaddingPower) :
483 fData(data), fUnitTag(unitTag), fStartTime(startTime), fEndTime(endTime),
484 fDCCorrected(dcCorrected), fZeroPaddingPower(zeroPaddingPower)
485{
486 // some necessary checks and initialization
487 if (fData == nullptr) {
488 std::cerr << std::endl << "**ERROR** PFourier::PFourier: no valid data" << std::endl << std::endl;
489 fValid = false;
490 return;
491 }
492
493 fValid = true;
494 fIn = nullptr;
495 fOut = nullptr;
496//as fPhCorrectedReFT = 0;
497
499
500 // calculate time resolution in (us)
501 fTimeResolution = fData->GetBinWidth(1);
502
503 // if endTime == 0 set it to the last time slot
504 if (fEndTime == 0.0) {
505 Int_t last = fData->GetNbinsX()-1;
506 fEndTime = fData->GetBinCenter(last);
507 }
508
509 // swap start and end time if necessary
510 if (fStartTime > fEndTime) {
511 Double_t keep = fStartTime;
513 fEndTime = keep;
514 }
515
516 // calculate start and end bin
517 fNoOfData = static_cast<UInt_t>((fEndTime-fStartTime)/fTimeResolution);
518
519 // check if zero padding is whished
520 if (fZeroPaddingPower > 0) {
521 UInt_t noOfBins = static_cast<UInt_t>(pow(2.0, static_cast<Double_t>(fZeroPaddingPower)));
522 if (noOfBins > fNoOfData)
523 fNoOfBins = noOfBins;
524 else
526 } else {
528 }
529
530 // calculate fourier resolution, depending on the units
531 Double_t resolution = 1.0/(fTimeResolution*fNoOfBins); // in MHz
532 switch (fUnitTag) {
534 fResolution = resolution/GAMMA_BAR_MUON;
535 break;
537 fResolution = 1e-4*resolution/GAMMA_BAR_MUON;
538 break;
540 fResolution = resolution;
541 break;
543 fResolution = 2.0*PI*resolution;
544 break;
545 default:
546 fValid = false;
547 return;
548 }
549
550 // allocate necessary memory
551 fIn = static_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex)*fNoOfBins));
552 fOut = static_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex)*fNoOfBins));
553
554 // check if memory allocation has been successful
555 if ((fIn == nullptr) || (fOut == nullptr)) {
556 fValid = false;
557 return;
558 }
559
560 // get the FFTW3 plan (see FFTW3 manual)
561 fFFTwPlan = fftw_plan_dft_1d(static_cast<Int_t>(fNoOfBins), fIn, fOut, FFTW_FORWARD, FFTW_ESTIMATE);
562
563 // check if a valid plan has been generated
564 if (!fFFTwPlan) {
565 fValid = false;
566 }
567}
568
569//--------------------------------------------------------------------------
570// Destructor
571//--------------------------------------------------------------------------
585{
586 if (fFFTwPlan)
587 fftw_destroy_plan(fFFTwPlan);
588 if (fIn)
589 fftw_free(fIn);
590 if (fOut)
591 fftw_free(fOut);
592//as if (fPhCorrectedReFT)
593//as delete fPhCorrectedReFT;
594}
595
596//--------------------------------------------------------------------------
597// Transform (public)
598//--------------------------------------------------------------------------
632void PFourier::Transform(UInt_t apodizationTag)
633{
634 if (!fValid)
635 return;
636
637 PrepareFFTwInputData(apodizationTag);
638
639 fftw_execute(fFFTwPlan);
640
641 // correct the phase for tstart != 0.0
642 // find the first bin >= fStartTime
643 Double_t shiftTime = 0.0;
644 for (Int_t i=1; i<fData->GetXaxis()->GetNbins(); i++) {
645 if (fData->GetXaxis()->GetBinCenter(i) >= fStartTime) {
646 shiftTime = fData->GetXaxis()->GetBinCenter(i);
647 break;
648 }
649 }
650
651 Double_t phase, re, im;
652 for (UInt_t i=0; i<fNoOfBins; i++) {
653 phase = 2.0*PI/(fTimeResolution*fNoOfBins) * i * shiftTime;
654 re = fOut[i][0] * cos(phase) + fOut[i][1] * sin(phase);
655 im = -fOut[i][0] * sin(phase) + fOut[i][1] * cos(phase);
656 fOut[i][0] = re;
657 fOut[i][1] = im;
658 }
659}
660
661//--------------------------------------------------------------------------
662// GetMaxFreq (public)
663//--------------------------------------------------------------------------
686{
687 UInt_t noOfFourierBins = 0;
688 if (fNoOfBins % 2 == 0)
689 noOfFourierBins = fNoOfBins/2;
690 else
691 noOfFourierBins = (fNoOfBins+1)/2;
692
693 return fResolution*noOfFourierBins;
694}
695
696//--------------------------------------------------------------------------
697// GetRealFourier (public)
698//--------------------------------------------------------------------------
731TH1F* PFourier::GetRealFourier(const Double_t scale)
732{
733 // check if valid flag is set
734 if (!fValid)
735 return nullptr;
736
737 // invoke realFourier
738 Char_t name[256];
739 Char_t title[256];
740 snprintf(name, sizeof(name), "%s_Fourier_Re", fData->GetName());
741 snprintf(title, sizeof(title), "%s_Fourier_Re", fData->GetTitle());
742
743 UInt_t noOfFourierBins = 0;
744 if (fNoOfBins % 2 == 0)
745 noOfFourierBins = fNoOfBins/2;
746 else
747 noOfFourierBins = (fNoOfBins+1)/2;
748
749 TH1F *realFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
750 if (realFourier == nullptr) {
751 fValid = false;
752 std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
753 return nullptr;
754 }
755
756 // fill realFourier vector
757 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
758 realFourier->SetBinContent(i+1, scale*fOut[i][0]);
759 realFourier->SetBinError(i+1, 0.0);
760 }
761 return realFourier;
762}
763
764//--------------------------------------------------------------------------
765// GetPhaseOptRealFourier (public, static)
766//--------------------------------------------------------------------------
817TH1F* PFourier::GetPhaseOptRealFourier(const TH1F *re, const TH1F *im, std::vector<Double_t> &phase,
818 const Double_t scale, const Double_t min, const Double_t max)
819{
820 if ((re == nullptr) || (im == nullptr))
821 return nullptr;
822
823 phase.resize(2); // c0, c1
824
825 const TAxis *axis = re->GetXaxis();
826
827 Int_t minBin = 1;
828 Int_t maxBin = axis->GetNbins();
829 Int_t noOfBins = axis->GetNbins();
830 Double_t res = axis->GetBinWidth(1);
831
832 // check if minimum frequency is given. If yes, get the proper minBin
833 if (min != -1.0) {
834 minBin = axis->FindFixBin(min);
835 if ((minBin == 0) || (minBin > maxBin)) {
836 minBin = 1;
837 std::cerr << "**WARNING** minimum frequency/field out of range. Will adopted it." << std::endl;
838 }
839 }
840
841 // check if maximum frequency is given. If yes, get the proper maxBin
842 if (max != -1.0) {
843 maxBin = axis->FindFixBin(max);
844 if ((maxBin == 0) || (maxBin > axis->GetNbins())) {
845 maxBin = axis->GetNbins();
846 std::cerr << "**WARNING** maximum frequency/field out of range. Will adopted it." << std::endl;
847 }
848 }
849
850 // copy the real/imag Fourier from min to max
851 std::vector<Double_t> realF, imagF;
852 for (Int_t i=minBin; i<=maxBin; i++) {
853 realF.push_back(re->GetBinContent(i));
854 imagF.push_back(im->GetBinContent(i));
855 }
856
857 // optimize real FT phase
858 PFTPhaseCorrection *phCorrectedReFT = new PFTPhaseCorrection(realF, imagF);
859 if (phCorrectedReFT == nullptr) {
860 std::cerr << std::endl << "**SEVERE ERROR** couldn't invoke PFTPhaseCorrection object." << std::endl;
861 return nullptr;
862 }
863
864 phCorrectedReFT->Minimize();
865 if (!phCorrectedReFT->IsValid()) {
866 std::cerr << std::endl << "**ERROR** could not find a valid phase correction minimum." << std::endl;
867 return nullptr;
868 }
869 phase[0] = phCorrectedReFT->GetPhaseCorrectionParam(0);
870 phase[1] = phCorrectedReFT->GetPhaseCorrectionParam(1);
871
872 // clean up
873 if (phCorrectedReFT) {
874 delete phCorrectedReFT;
875 phCorrectedReFT = nullptr;
876 }
877 realF.clear();
878 imagF.clear();
879
880 // invoke the real phase optimised histo to be filled. Caller is the owner!
881 Char_t name[256];
882 Char_t title[256];
883 snprintf(name, sizeof(name), "%s_Fourier_PhOptRe", re->GetName());
884 snprintf(title, sizeof(title), "%s_Fourier_PhOptRe", re->GetTitle());
885
886 TH1F *realPhaseOptFourier = new TH1F(name, title, noOfBins, -res/2.0, static_cast<Double_t>(noOfBins-1)*res+res/2.0);
887 if (realPhaseOptFourier == nullptr) {
888 std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
889 return nullptr;
890 }
891
892 // fill realFourier vector
893 Double_t ph;
894 for (Int_t i=0; i<noOfBins; i++) {
895 ph = phase[0] + phase[1] * static_cast<Double_t>(i-static_cast<Int_t>(minBin)) / static_cast<Double_t>(maxBin-minBin);
896 realPhaseOptFourier->SetBinContent(i+1, scale*(re->GetBinContent(i+1)*cos(ph) - im->GetBinContent(i+1)*sin(ph)));
897 realPhaseOptFourier->SetBinError(i+1, 0.0);
898 }
899
900 return realPhaseOptFourier;
901}
902
903//--------------------------------------------------------------------------
904// GetImaginaryFourier (public)
905//--------------------------------------------------------------------------
931TH1F* PFourier::GetImaginaryFourier(const Double_t scale)
932{
933 // check if valid flag is set
934 if (!fValid)
935 return nullptr;
936
937 // invoke imaginaryFourier
938 Char_t name[256];
939 Char_t title[256];
940 snprintf(name, sizeof(name), "%s_Fourier_Im", fData->GetName());
941 snprintf(title, sizeof(title), "%s_Fourier_Im", fData->GetTitle());
942
943 UInt_t noOfFourierBins = 0;
944 if (fNoOfBins % 2 == 0)
945 noOfFourierBins = fNoOfBins/2;
946 else
947 noOfFourierBins = (fNoOfBins+1)/2;
948
949 TH1F* imaginaryFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
950 if (imaginaryFourier == nullptr) {
951 fValid = false;
952 std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the imaginary part of the Fourier transform." << std::endl;
953 return nullptr;
954 }
955
956 // fill imaginaryFourier vector
957 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
958 imaginaryFourier->SetBinContent(i+1, scale*fOut[i][1]);
959 imaginaryFourier->SetBinError(i+1, 0.0);
960 }
961
962 return imaginaryFourier;
963}
964
965//--------------------------------------------------------------------------
966// GetPowerFourier (public)
967//--------------------------------------------------------------------------
1001TH1F* PFourier::GetPowerFourier(const Double_t scale)
1002{
1003 // check if valid flag is set
1004 if (!fValid)
1005 return nullptr;
1006
1007 // invoke powerFourier
1008 Char_t name[256];
1009 Char_t title[256];
1010 snprintf(name, sizeof(name), "%s_Fourier_Pwr", fData->GetName());
1011 snprintf(title, sizeof(title), "%s_Fourier_Pwr", fData->GetTitle());
1012
1013 UInt_t noOfFourierBins = 0;
1014 if (fNoOfBins % 2 == 0)
1015 noOfFourierBins = fNoOfBins/2;
1016 else
1017 noOfFourierBins = (fNoOfBins+1)/2;
1018
1019 TH1F* pwrFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
1020 if (pwrFourier == nullptr) {
1021 fValid = false;
1022 std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the power part of the Fourier transform." << std::endl;
1023 return nullptr;
1024 }
1025
1026 // fill powerFourier vector
1027 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
1028 pwrFourier->SetBinContent(i+1, scale*sqrt(fOut[i][0]*fOut[i][0]+fOut[i][1]*fOut[i][1]));
1029 pwrFourier->SetBinError(i+1, 0.0);
1030 }
1031
1032 return pwrFourier;
1033}
1034
1035//--------------------------------------------------------------------------
1036// GetPhaseFourier (public)
1037//--------------------------------------------------------------------------
1073TH1F* PFourier::GetPhaseFourier(const Double_t scale)
1074{
1075 // check if valid flag is set
1076 if (!fValid)
1077 return nullptr;
1078
1079 // invoke phaseFourier
1080 Char_t name[256];
1081 Char_t title[256];
1082 snprintf(name, sizeof(name), "%s_Fourier_Phase", fData->GetName());
1083 snprintf(title, sizeof(title), "%s_Fourier_Phase", fData->GetTitle());
1084
1085 UInt_t noOfFourierBins = 0;
1086 if (fNoOfBins % 2 == 0)
1087 noOfFourierBins = fNoOfBins/2;
1088 else
1089 noOfFourierBins = (fNoOfBins+1)/2;
1090
1091 TH1F* phaseFourier = new TH1F(name, title, static_cast<Int_t>(noOfFourierBins), -fResolution/2.0, static_cast<Double_t>(noOfFourierBins-1)*fResolution+fResolution/2.0);
1092 if (phaseFourier == nullptr) {
1093 fValid = false;
1094 std::cerr << std::endl << "**SEVERE ERROR** couldn't allocate memory for the phase part of the Fourier transform." << std::endl;
1095 return nullptr;
1096 }
1097
1098 // fill phaseFourier vector
1099 Double_t value = 0.0;
1100 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
1101 // calculate the phase
1102 if (fOut[i][0] == 0.0) {
1103 if (fOut[i][1] >= 0.0)
1104 value = PI_HALF;
1105 else
1106 value = -PI_HALF;
1107 } else {
1108 value = atan(fOut[i][1]/fOut[i][0]);
1109 // check sector
1110 if (fOut[i][0] < 0.0) {
1111 if (fOut[i][1] > 0.0)
1112 value = PI + value;
1113 else
1114 value = PI - value;
1115 }
1116 }
1117 phaseFourier->SetBinContent(i+1, scale*value);
1118 phaseFourier->SetBinError(i+1, 0.0);
1119 }
1120
1121 return phaseFourier;
1122}
1123
1124//--------------------------------------------------------------------------
1125// PrepareFFTwInputData (private)
1126//--------------------------------------------------------------------------
1155void PFourier::PrepareFFTwInputData(UInt_t apodizationTag)
1156{
1157 // 1st find t==0. fData start at times t<0!!
1158 Int_t t0bin = -1;
1159 for (Int_t i=1; i<fData->GetNbinsX(); i++) {
1160 if (fData->GetBinCenter(i) >= 0.0) {
1161 t0bin = i;
1162 break;
1163 }
1164 }
1165
1166 Int_t ival = static_cast<Int_t>(fStartTime/fTimeResolution) + t0bin;
1167 UInt_t start = 0;
1168 if (ival >= 0) {
1169 start = static_cast<UInt_t>(ival);
1170 }
1171
1172 Double_t mean = 0.0;
1173 if (fDCCorrected) {
1174 for (UInt_t i=0; i<fNoOfData; i++) {
1175 mean += fData->GetBinContent(static_cast<Int_t>(i+start));
1176 }
1177 mean /= static_cast<Double_t>(fNoOfData);
1178 }
1179
1180 // 2nd fill fIn
1181 for (UInt_t i=0; i<fNoOfData; i++) {
1182 fIn[i][0] = fData->GetBinContent(static_cast<Int_t>(i+start)) - mean;
1183 fIn[i][1] = 0.0;
1184 }
1185 for (UInt_t i=fNoOfData; i<fNoOfBins; i++) {
1186 fIn[i][0] = 0.0;
1187 fIn[i][1] = 0.0;
1188 }
1189
1190 // 3rd apodize data (if wished)
1191 ApodizeData(static_cast<Int_t>(apodizationTag));
1192}
1193
1194//--------------------------------------------------------------------------
1195// ApodizeData (private)
1196//--------------------------------------------------------------------------
1226void PFourier::ApodizeData(Int_t apodizationTag) {
1227
1228 const Double_t cweak[3] = { 0.384093, -0.087577, 0.703484 };
1229 const Double_t cmedium[3] = { 0.152442, -0.136176, 0.983734 };
1230 const Double_t cstrong[3] = { 0.045335, 0.554883, 0.399782 };
1231
1232 Double_t c[5];
1233 for (UInt_t i=0; i<5; i++) {
1234 c[i] = 0.0;
1235 }
1236
1237 switch (apodizationTag) {
1238 case F_APODIZATION_NONE:
1239 return;
1240 case F_APODIZATION_WEAK:
1241 c[0] = cweak[0]+cweak[1]+cweak[2];
1242 c[1] = -(cweak[1]+2.0*cweak[2]);
1243 c[2] = cweak[2];
1244 break;
1246 c[0] = cmedium[0]+cmedium[1]+cmedium[2];
1247 c[1] = -(cmedium[1]+2.0*cmedium[2]);
1248 c[2] = cmedium[2];
1249 break;
1251 c[0] = cstrong[0]+cstrong[1]+cstrong[2];
1252 c[1] = -2.0*(cstrong[1]+2.0*cstrong[2]);
1253 c[2] = cstrong[1]+6.0*cstrong[2];
1254 c[3] = -4.0*cstrong[2];
1255 c[4] = cstrong[2];
1256 break;
1257 default:
1258 std::cerr << std::endl << ">> **ERROR** User Apodization tag " << apodizationTag << " unknown, sorry ..." << std::endl;
1259 break;
1260 }
1261
1262 Double_t q;
1263 for (UInt_t i=0; i<fNoOfData; i++) {
1264 q = c[0];
1265 for (UInt_t j=1; j<5; j++) {
1266 q += c[j] * pow(static_cast<Double_t>(i)/static_cast<Double_t>(fNoOfData), 2.0*static_cast<Double_t>(j));
1267 }
1268 fIn[i][0] *= q;
1269 }
1270}
#define PI
Definition PFourier.cpp:45
#define PI_HALF
Definition PFourier.cpp:46
#define F_APODIZATION_WEAK
Weak apodization (gentle roll-off at edges)
Definition PFourier.h:55
#define F_APODIZATION_STRONG
Strong apodization (heavy roll-off for best frequency resolution)
Definition PFourier.h:59
#define F_APODIZATION_NONE
No apodization (rectangular window)
Definition PFourier.h:53
#define F_APODIZATION_MEDIUM
Medium apodization (moderate roll-off)
Definition PFourier.h:57
#define FOURIER_UNIT_FREQ
Frequency in MHz.
Definition PMusr.h:276
#define FOURIER_UNIT_GAUSS
Magnetic field in Gauss (G)
Definition PMusr.h:272
#define GAMMA_BAR_MUON
Definition PMusr.h:138
#define FOURIER_UNIT_CYCLES
Angular frequency in Mc/s (Mega-cycles per second)
Definition PMusr.h:278
#define FOURIER_UNIT_TESLA
Magnetic field in Tesla (T)
Definition PMusr.h:274
virtual void Minimize()
Definition PFourier.cpp:144
PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1)
Definition PFourier.cpp:67
virtual void CalcRealPhFTDerivative() const
Definition PFourier.cpp:307
Int_t fMinBin
1st derivative of fRealPh
Definition PFourier.h:153
std::vector< Double_t > fRealPh
original imag Fourier data set
Definition PFourier.h:149
virtual void Init()
keeps the minimum of the entropy/penalty minimization
Definition PFourier.cpp:248
std::vector< Double_t > fRealPhD
phased imag Fourier data set
Definition PFourier.h:151
virtual Bool_t IsValid()
Definition PFourier.h:109
virtual Double_t Entropy() const
Definition PFourier.cpp:376
std::vector< Double_t > fImagPh
phased real Fourier data set
Definition PFourier.h:150
std::vector< Double_t > fReal
Definition PFourier.h:147
virtual void CalcPhasedFT() const
Definition PFourier.cpp:276
virtual Double_t operator()(const std::vector< Double_t > &) const
Definition PFourier.cpp:425
Double_t fPh_c1
constant part of the phase dispersion used for the phase correction
Definition PFourier.h:156
Double_t fMin
gamma parameter to balance between entropy and penalty
Definition PFourier.h:158
Double_t fGamma
linear part of the phase dispersion used for the phase correction
Definition PFourier.h:157
Double_t fPh_c0
maximum bin from Fourier range to be used for the phase correction estimate
Definition PFourier.h:155
virtual Double_t Penalty() const
Definition PFourier.cpp:339
virtual Double_t GetMinimum()
Definition PFourier.cpp:223
Int_t fMaxBin
minimum bin from Fourier range to be used for the phase correction estimate
Definition PFourier.h:154
virtual Double_t GetPhaseCorrectionParam(UInt_t idx)
Definition PFourier.cpp:195
std::vector< Double_t > fImag
original real Fourier data set
Definition PFourier.h:148
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:482
fftw_plan fFFTwPlan
fftw plan (see FFTW3 User Manual)
Definition PFourier.h:322
virtual void ApodizeData(Int_t apodizationTag)
Int_t fApodization
0=none, 1=weak, 2=medium, 3=strong
Definition PFourier.h:311
Bool_t fValid
true = all boundary conditions fullfilled and hence a Fourier transform can be performed.
Definition PFourier.h:308
UInt_t fZeroPaddingPower
power for zero padding, if set < 0 no zero padding will be done
Definition PFourier.h:317
UInt_t fNoOfBins
number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding
Definition PFourier.h:321
virtual TH1F * GetImaginaryFourier(const Double_t scale=1.0)
Definition PFourier.cpp:931
UInt_t fNoOfData
number of bins in the time interval between fStartTime and fStopTime
Definition PFourier.h:320
Double_t fResolution
Fourier resolution (field, frequency, or angular frequency)
Definition PFourier.h:318
virtual void Transform(UInt_t apodizationTag=0)
Definition PFourier.cpp:632
virtual TH1F * GetRealFourier(const Double_t scale=1.0)
Definition PFourier.cpp:731
virtual Double_t GetMaxFreq()
Definition PFourier.cpp:685
virtual TH1F * GetPowerFourier(const Double_t scale=1.0)
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:817
TH1F * fData
data histogram to be Fourier transformed.
Definition PFourier.h:306
virtual void PrepareFFTwInputData(UInt_t apodizationTag)
virtual TH1F * GetPhaseFourier(const Double_t scale=1.0)
Int_t fUnitTag
1=Field Units (G), 2=Field Units (T), 3=Frequency Units (MHz), 4=Angular Frequency Units (Mc/s)
Definition PFourier.h:309
Double_t fStartTime
start time of the data histogram
Definition PFourier.h:314
Double_t fTimeResolution
time resolution of the data histogram in (us)
Definition PFourier.h:313
fftw_complex * fIn
real part of the Fourier transform
Definition PFourier.h:323
fftw_complex * fOut
imaginary part of the Fourier transform
Definition PFourier.h:324
Double_t fEndTime
end time of the data histogram
Definition PFourier.h:315
Bool_t fDCCorrected
if true, removed DC offset from signal before Fourier transformation, otherwise not
Definition PFourier.h:316
virtual ~PFourier()
Definition PFourier.cpp:584