musrfit  1.9.2
PFourier.h
Go to the documentation of this file.
1 /***************************************************************************
2 
3  PFourier.h
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 #ifndef _PFOURIER_H_
31 #define _PFOURIER_H_
32 
33 #include <vector>
34 
35 #include "fftw3.h"
36 
37 #include <TH1F.h>
38 
39 #include "Minuit2/FCNBase.h"
40 
41 #include "PMusr.h"
42 
43 #define F_APODIZATION_NONE 1
44 #define F_APODIZATION_WEAK 2
45 #define F_APODIZATION_MEDIUM 3
46 #define F_APODIZATION_STRONG 4
47 
51 class PFTPhaseCorrection : public ROOT::Minuit2::FCNBase
52 {
53  public:
54  PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1);
55  PFTPhaseCorrection(std::vector<Double_t> &reFT, std::vector<Double_t> &imFT, const Int_t minBin=-1, const Int_t maxBin=-1);
56  virtual ~PFTPhaseCorrection() {}
57 
58  virtual Bool_t IsValid() { return fValid; }
59  virtual void Minimize();
60 
61  virtual void SetGamma(const Double_t gamma) { fGamma = gamma; }
62  virtual void SetPh(const Double_t c0, const Double_t c1) { fPh_c0 = c0; fPh_c1 = c1; CalcPhasedFT(); CalcRealPhFTDerivative(); }
63 
64  virtual Double_t GetGamma() { return fGamma; }
65  virtual Double_t GetPhaseCorrectionParam(UInt_t idx);
66  virtual Double_t GetMinimum();
67 
68  private:
69  Bool_t fValid;
70 
71  std::vector<Double_t> fReal;
72  std::vector<Double_t> fImag;
73  mutable std::vector<Double_t> fRealPh;
74  mutable std::vector<Double_t> fImagPh;
75  mutable std::vector<Double_t> fRealPhD;
76 
77  Int_t fMinBin;
78  Int_t fMaxBin;
79  mutable Double_t fPh_c0;
80  mutable Double_t fPh_c1;
81  Double_t fGamma;
82  Double_t fMin;
83 
84  virtual void Init();
85  virtual void CalcPhasedFT() const;
86  virtual void CalcRealPhFTDerivative() const;
87  virtual Double_t Penalty() const;
88  virtual Double_t Entropy() const;
89 
90  virtual Double_t Up() const { return 1.0; }
91  virtual Double_t operator()(const std::vector<Double_t>&) const;
92 };
93 
97 class PFourier
98 {
99  public:
100  PFourier(TH1F *data, Int_t unitTag,
101  Double_t startTime = 0.0, Double_t endTime = 0.0,
102  Bool_t dcCorrected = false, UInt_t zeroPaddingPower = 0);
103  virtual ~PFourier();
104 
105  virtual void Transform(UInt_t apodizationTag = 0);
106 
107  virtual const char* GetDataTitle() { return fData->GetTitle(); }
108  virtual const Int_t GetUnitTag() { return fUnitTag; }
109  virtual Double_t GetResolution() { return fResolution; }
110  virtual Double_t GetMaxFreq();
111  virtual TH1F* GetRealFourier(const Double_t scale = 1.0);
112 //as virtual TH1F* GetPhaseOptRealFourier(std::vector<Double_t> &phase, const Double_t scale = 1.0, const Double_t min = -1.0, const Double_t max = -1.0);
113  virtual TH1F* GetImaginaryFourier(const Double_t scale = 1.0);
114  virtual TH1F* GetPowerFourier(const Double_t scale = 1.0);
115  virtual TH1F* GetPhaseFourier(const Double_t scale = 1.0);
116 
117  static TH1F* GetPhaseOptRealFourier(const TH1F *re, const TH1F *im, std::vector<Double_t> &phase,
118  const Double_t scale = 1.0, const Double_t min = -1.0, const Double_t max = -1.0);
119 
120  virtual Bool_t IsValid() { return fValid; }
121 
122  private:
123  TH1F *fData;
124 
125  Bool_t fValid;
126  Int_t fUnitTag;
127 
128  Int_t fApodization;
129 
130  Double_t fTimeResolution;
131  Double_t fStartTime;
132  Double_t fEndTime;
133  Bool_t fDCCorrected;
135  Double_t fResolution;
136 
137  UInt_t fNoOfData;
138  UInt_t fNoOfBins;
139  fftw_plan fFFTwPlan;
140  fftw_complex *fIn;
141  fftw_complex *fOut;
142 
143 //as PFTPhaseCorrection *fPhCorrectedReFT;
144 
145  virtual void PrepareFFTwInputData(UInt_t apodizationTag);
146  virtual void ApodizeData(Int_t apodizationTag);
147 };
148 
149 #endif // _PFOURIER_H_
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
virtual const Int_t GetUnitTag()
Definition: PFourier.h:108
UInt_t fNoOfData
number of bins in the time interval between fStartTime and fStopTime
Definition: PFourier.h:137
virtual Double_t GetGamma()
Definition: PFourier.h:64
virtual Double_t Entropy() const
Definition: PFourier.cpp:228
virtual Bool_t IsValid()
Definition: PFourier.h:120
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
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
virtual Double_t Up() const
Definition: PFourier.h:90
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
virtual Double_t GetMaxFreq()
Definition: PFourier.cpp:433
virtual void Transform(UInt_t apodizationTag=0)
Definition: PFourier.cpp:398
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
virtual Double_t GetResolution()
Definition: PFourier.h:109
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
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
virtual void SetPh(const Double_t c0, const Double_t c1)
Definition: PFourier.h:62
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
virtual const char * GetDataTitle()
Definition: PFourier.h:107
virtual Double_t GetPhaseCorrectionParam(UInt_t idx)
Definition: PFourier.cpp:123
virtual ~PFTPhaseCorrection()
Definition: PFourier.h:56
Bool_t fValid
true = all boundary conditions fullfilled and hence a Fourier transform can be performed.
Definition: PFourier.h:125
Int_t fMinBin
1st derivative of fRealPh
Definition: PFourier.h:77
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
virtual void SetGamma(const Double_t gamma)
Definition: PFourier.h:61
Double_t fMin
gamma parameter to balance between entropy and penalty
Definition: PFourier.h:82
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