38 #include "Minuit2/FunctionMinimum.h" 39 #include "Minuit2/MnUserParameters.h" 40 #include "Minuit2/MnMinimize.h" 45 #define PI 3.14159265358979312 46 #define PI_HALF 1.57079632679489656 57 fMinBin(minBin), fMaxBin(maxBin)
69 fReal(reFT), fImag(imFT), fMinBin(minBin), fMaxBin(maxBin)
73 Int_t realSize =
static_cast<Int_t
>(
fReal.size());
94 ROOT::Minuit2::MnUserParameters upar;
96 upar.Add(
"c0",
fPh_c0, 2.0);
97 upar.Add(
"c1",
fPh_c1, 2.0);
100 ROOT::Minuit2::MnMinimize mn_min(*
this, upar);
103 ROOT::Minuit2::FunctionMinimum min = mn_min();
106 fPh_c0 = min.UserState().Value(
"c0");
107 fPh_c1 = min.UserState().Value(
"c1");
112 std::cout << std::endl <<
">> **WARNING** minimize failed to find a minimum for the real FT phase correction ..." << std::endl;
132 std::cerr <<
">> **ERROR** requested phase correction parameter with index=" << idx <<
" does not exist!" << std::endl;
146 std::cerr <<
">> **ERROR** requested minimum is invalid!" << std::endl;
179 for (UInt_t i=0; i<
fRealPh.size(); i++) {
180 w =
static_cast<Double_t
>(i) / static_cast<Double_t>(
fReal.size());
200 for (UInt_t i=1; i<
fRealPh.size()-1; i++)
212 Double_t penalty = 0.0;
235 Double_t entropy = 0.0;
236 Double_t dval = 0.0, hh = 0.0;
240 if (dval > 1.0e-15) {
242 entropy -= hh * log(hh);
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)
289 if (
fData ==
nullptr) {
290 std::cerr << std::endl <<
"**ERROR** PFourier::PFourier: no valid data" << std::endl << std::endl;
307 Int_t last =
fData->GetNbinsX()-1;
323 UInt_t noOfBins =
static_cast<UInt_t
>(pow(2.0, static_cast<Double_t>(
fZeroPaddingPower)));
353 fIn =
static_cast<fftw_complex *
>(fftw_malloc(
sizeof(fftw_complex)*
fNoOfBins));
354 fOut =
static_cast<fftw_complex *
>(fftw_malloc(
sizeof(fftw_complex)*
fNoOfBins));
357 if ((
fIn ==
nullptr) || (
fOut ==
nullptr)) {
409 Double_t shiftTime = 0.0;
410 for (Int_t i=1; i<
fData->GetXaxis()->GetNbins(); i++) {
412 shiftTime =
fData->GetXaxis()->GetBinCenter(i);
417 Double_t phase, re, im;
420 re =
fOut[i][0] * cos(phase) +
fOut[i][1] * sin(phase);
421 im = -
fOut[i][0] * sin(phase) +
fOut[i][1] * cos(phase);
435 UInt_t noOfFourierBins = 0;
461 snprintf(name,
sizeof(name),
"%s_Fourier_Re",
fData->GetName());
462 snprintf(title,
sizeof(title),
"%s_Fourier_Re",
fData->GetTitle());
464 UInt_t noOfFourierBins = 0;
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) {
473 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
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);
501 const Double_t scale,
const Double_t min,
const Double_t max)
503 if ((re ==
nullptr) || (im ==
nullptr))
508 const TAxis *axis = re->GetXaxis();
511 Int_t maxBin = axis->GetNbins();
512 Int_t noOfBins = axis->GetNbins();
513 Double_t res = axis->GetBinWidth(1);
517 minBin = axis->FindFixBin(min);
518 if ((minBin == 0) || (minBin > maxBin)) {
520 std::cerr <<
"**WARNING** minimum frequency/field out of range. Will adopted it." << std::endl;
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;
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));
542 if (phCorrectedReFT ==
nullptr) {
543 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't invoke PFTPhaseCorrection object." << std::endl;
548 if (!phCorrectedReFT->
IsValid()) {
549 std::cerr << std::endl <<
"**ERROR** could not find a valid phase correction minimum." << std::endl;
556 if (phCorrectedReFT) {
557 delete phCorrectedReFT;
558 phCorrectedReFT =
nullptr;
566 snprintf(name,
sizeof(name),
"%s_Fourier_PhOptRe", re->GetName());
567 snprintf(title,
sizeof(title),
"%s_Fourier_PhOptRe", re->GetTitle());
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;
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);
583 return realPhaseOptFourier;
603 snprintf(name,
sizeof(name),
"%s_Fourier_Im",
fData->GetName());
604 snprintf(title,
sizeof(title),
"%s_Fourier_Im",
fData->GetTitle());
606 UInt_t noOfFourierBins = 0;
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) {
615 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the imaginary part of the Fourier transform." << std::endl;
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);
625 return imaginaryFourier;
645 snprintf(name,
sizeof(name),
"%s_Fourier_Pwr",
fData->GetName());
646 snprintf(title,
sizeof(title),
"%s_Fourier_Pwr",
fData->GetTitle());
648 UInt_t noOfFourierBins = 0;
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) {
657 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the power part of the Fourier transform." << std::endl;
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);
687 snprintf(name,
sizeof(name),
"%s_Fourier_Phase",
fData->GetName());
688 snprintf(title,
sizeof(title),
"%s_Fourier_Phase",
fData->GetTitle());
690 UInt_t noOfFourierBins = 0;
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) {
699 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the phase part of the Fourier transform." << std::endl;
704 Double_t value = 0.0;
705 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
707 if (
fOut[i][0] == 0.0) {
708 if (
fOut[i][1] >= 0.0)
713 value = atan(
fOut[i][1]/
fOut[i][0]);
715 if (
fOut[i][0] < 0.0) {
716 if (
fOut[i][1] > 0.0)
722 phaseFourier->SetBinContent(i+1, scale*value);
723 phaseFourier->SetBinError(i+1, 0.0);
743 for (Int_t i=1; i<
fData->GetNbinsX(); i++) {
744 if (
fData->GetBinCenter(i) >= 0.0) {
753 start =
static_cast<UInt_t
>(ival);
759 mean +=
fData->GetBinContent(static_cast<Int_t>(i+start));
761 mean /=
static_cast<Double_t
>(
fNoOfData);
766 fIn[i][0] =
fData->GetBinContent(static_cast<Int_t>(i+start)) - mean;
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 };
795 for (UInt_t i=0; i<5; i++) {
799 switch (apodizationTag) {
803 c[0] = cweak[0]+cweak[1]+cweak[2];
804 c[1] = -(cweak[1]+2.0*cweak[2]);
808 c[0] = cmedium[0]+cmedium[1]+cmedium[2];
809 c[1] = -(cmedium[1]+2.0*cmedium[2]);
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];
820 std::cerr << std::endl <<
">> **ERROR** User Apodization tag " << apodizationTag <<
" unknown, sorry ..." << std::endl;
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));
fftw_complex * fOut
imaginary part of the Fourier transform
Double_t fPh_c0
maximum bin from Fourier range to be used for the phase correction estimate
virtual void PrepareFFTwInputData(UInt_t apodizationTag)
virtual void CalcRealPhFTDerivative() const
Double_t fPh_c1
constant part of the phase dispersion used for the phase correction
UInt_t fNoOfData
number of bins in the time interval between fStartTime and fStopTime
virtual Double_t Entropy() const
Double_t fGamma
linear part of the phase dispersion used for the phase correction
fftw_plan fFFTwPlan
fftw plan (see FFTW3 User Manual)
#define F_APODIZATION_STRONG
Int_t fApodization
0=none, 1=weak, 2=medium, 3=strong
Int_t fMaxBin
minimum bin from Fourier range to be used for the phase correction estimate
TH1F * fData
data histogram to be Fourier transformed.
fftw_complex * fIn
real part of the Fourier transform
PFourier(TH1F *data, Int_t unitTag, Double_t startTime=0.0, Double_t endTime=0.0, Bool_t dcCorrected=false, UInt_t zeroPaddingPower=0)
Double_t fResolution
Fourier resolution (field, frequency, or angular frequency)
std::vector< Double_t > fReal
#define FOURIER_UNIT_CYCLES
virtual Double_t GetMaxFreq()
virtual void Transform(UInt_t apodizationTag=0)
#define FOURIER_UNIT_FREQ
PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1)
virtual void CalcPhasedFT() const
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)
virtual TH1F * GetRealFourier(const Double_t scale=1.0)
#define F_APODIZATION_NONE
virtual Double_t operator()(const std::vector< Double_t > &) const
Bool_t fDCCorrected
if true, removed DC offset from signal before Fourier transformation, otherwise not ...
#define F_APODIZATION_WEAK
std::vector< Double_t > fImag
original real Fourier data set
Double_t fStartTime
start time of the data histogram
std::vector< Double_t > fImagPh
phased real Fourier data set
virtual void ApodizeData(Int_t apodizationTag)
virtual TH1F * GetPowerFourier(const Double_t scale=1.0)
virtual TH1F * GetPhaseFourier(const Double_t scale=1.0)
virtual Double_t GetMinimum()
#define FOURIER_UNIT_GAUSS
virtual Double_t GetPhaseCorrectionParam(UInt_t idx)
#define FOURIER_UNIT_TESLA
Bool_t fValid
true = all boundary conditions fullfilled and hence a Fourier transform can be performed.
Int_t fMinBin
1st derivative of fRealPh
#define F_APODIZATION_MEDIUM
std::vector< Double_t > fRealPhD
phased imag Fourier data set
Int_t fUnitTag
1=Field Units (G), 2=Field Units (T), 3=Frequency Units (MHz), 4=Angular Frequency Units (Mc/s) ...
virtual void Init()
keeps the minimum of the entropy/penalty minimization
std::vector< Double_t > fRealPh
original imag Fourier data set
virtual Double_t Penalty() const
Double_t fMin
gamma parameter to balance between entropy and penalty
virtual TH1F * GetImaginaryFourier(const Double_t scale=1.0)
UInt_t fNoOfBins
number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding ...
Double_t fTimeResolution
time resolution of the data histogram in (us)
Double_t fEndTime
end time of the data histogram
UInt_t fZeroPaddingPower
power for zero padding, if set < 0 no zero padding will be done