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
105 Int_t realSize =
static_cast<Int_t
>(
fReal.size());
147 ROOT::Minuit2::MnUserParameters upar;
149 upar.Add(
"c0",
fPh_c0, 2.0);
150 upar.Add(
"c1",
fPh_c1, 2.0);
153 ROOT::Minuit2::MnMinimize mn_min(*
this, upar);
156 ROOT::Minuit2::FunctionMinimum min = mn_min();
159 fPh_c0 = min.UserState().Value(
"c0");
160 fPh_c1 = min.UserState().Value(
"c1");
165 std::cout << std::endl <<
">> **WARNING** minimize failed to find a minimum for the real FT phase correction ..." << std::endl;
204 std::cerr <<
">> **ERROR** requested phase correction parameter with index=" << idx <<
" does not exist!" << std::endl;
226 std::cerr <<
">> **ERROR** requested minimum is invalid!" << std::endl;
281 for (UInt_t i=0; i<
fRealPh.size(); i++) {
282 w =
static_cast<Double_t
>(i) /
static_cast<Double_t
>(
fReal.size());
314 for (UInt_t i=1; i<
fRealPh.size()-1; i++)
341 Double_t penalty = 0.0;
383 Double_t entropy = 0.0;
384 Double_t dval = 0.0, hh = 0.0;
388 if (dval > 1.0e-15) {
390 entropy -= hh * log(hh);
482PFourier::PFourier(TH1F *data, Int_t unitTag, Double_t startTime, Double_t endTime, Bool_t dcCorrected, UInt_t zeroPaddingPower) :
487 if (
fData ==
nullptr) {
488 std::cerr << std::endl <<
"**ERROR** PFourier::PFourier: no valid data" << std::endl << std::endl;
505 Int_t last =
fData->GetNbinsX()-1;
521 UInt_t noOfBins =
static_cast<UInt_t
>(pow(2.0,
static_cast<Double_t
>(
fZeroPaddingPower)));
551 fIn =
static_cast<fftw_complex *
>(fftw_malloc(
sizeof(fftw_complex)*
fNoOfBins));
552 fOut =
static_cast<fftw_complex *
>(fftw_malloc(
sizeof(fftw_complex)*
fNoOfBins));
555 if ((
fIn ==
nullptr) || (
fOut ==
nullptr)) {
643 Double_t shiftTime = 0.0;
644 for (Int_t i=1; i<
fData->GetXaxis()->GetNbins(); i++) {
646 shiftTime =
fData->GetXaxis()->GetBinCenter(i);
651 Double_t phase, re, im;
654 re =
fOut[i][0] * cos(phase) +
fOut[i][1] * sin(phase);
655 im = -
fOut[i][0] * sin(phase) +
fOut[i][1] * cos(phase);
687 UInt_t noOfFourierBins = 0;
740 snprintf(name,
sizeof(name),
"%s_Fourier_Re",
fData->GetName());
741 snprintf(title,
sizeof(title),
"%s_Fourier_Re",
fData->GetTitle());
743 UInt_t noOfFourierBins = 0;
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) {
752 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the real part of the Fourier transform." << std::endl;
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);
818 const Double_t scale,
const Double_t min,
const Double_t max)
820 if ((re ==
nullptr) || (im ==
nullptr))
825 const TAxis *axis = re->GetXaxis();
828 Int_t maxBin = axis->GetNbins();
829 Int_t noOfBins = axis->GetNbins();
830 Double_t res = axis->GetBinWidth(1);
834 minBin = axis->FindFixBin(min);
835 if ((minBin == 0) || (minBin > maxBin)) {
837 std::cerr <<
"**WARNING** minimum frequency/field out of range. Will adopted it." << std::endl;
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;
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));
859 if (phCorrectedReFT ==
nullptr) {
860 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't invoke PFTPhaseCorrection object." << std::endl;
865 if (!phCorrectedReFT->
IsValid()) {
866 std::cerr << std::endl <<
"**ERROR** could not find a valid phase correction minimum." << std::endl;
873 if (phCorrectedReFT) {
874 delete phCorrectedReFT;
875 phCorrectedReFT =
nullptr;
883 snprintf(name,
sizeof(name),
"%s_Fourier_PhOptRe", re->GetName());
884 snprintf(title,
sizeof(title),
"%s_Fourier_PhOptRe", re->GetTitle());
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;
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);
900 return realPhaseOptFourier;
940 snprintf(name,
sizeof(name),
"%s_Fourier_Im",
fData->GetName());
941 snprintf(title,
sizeof(title),
"%s_Fourier_Im",
fData->GetTitle());
943 UInt_t noOfFourierBins = 0;
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) {
952 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the imaginary part of the Fourier transform." << std::endl;
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);
962 return imaginaryFourier;
1010 snprintf(name,
sizeof(name),
"%s_Fourier_Pwr",
fData->GetName());
1011 snprintf(title,
sizeof(title),
"%s_Fourier_Pwr",
fData->GetTitle());
1013 UInt_t noOfFourierBins = 0;
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) {
1022 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the power part of the Fourier transform." << std::endl;
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);
1082 snprintf(name,
sizeof(name),
"%s_Fourier_Phase",
fData->GetName());
1083 snprintf(title,
sizeof(title),
"%s_Fourier_Phase",
fData->GetTitle());
1085 UInt_t noOfFourierBins = 0;
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) {
1094 std::cerr << std::endl <<
"**SEVERE ERROR** couldn't allocate memory for the phase part of the Fourier transform." << std::endl;
1099 Double_t value = 0.0;
1100 for (Int_t i=0; i<static_cast<Int_t>(noOfFourierBins); i++) {
1102 if (
fOut[i][0] == 0.0) {
1103 if (
fOut[i][1] >= 0.0)
1108 value = atan(
fOut[i][1]/
fOut[i][0]);
1110 if (
fOut[i][0] < 0.0) {
1111 if (
fOut[i][1] > 0.0)
1117 phaseFourier->SetBinContent(i+1, scale*value);
1118 phaseFourier->SetBinError(i+1, 0.0);
1121 return phaseFourier;
1159 for (Int_t i=1; i<
fData->GetNbinsX(); i++) {
1160 if (
fData->GetBinCenter(i) >= 0.0) {
1169 start =
static_cast<UInt_t
>(ival);
1172 Double_t mean = 0.0;
1175 mean +=
fData->GetBinContent(
static_cast<Int_t
>(i+start));
1177 mean /=
static_cast<Double_t
>(
fNoOfData);
1182 fIn[i][0] =
fData->GetBinContent(
static_cast<Int_t
>(i+start)) - mean;
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 };
1233 for (UInt_t i=0; i<5; i++) {
1237 switch (apodizationTag) {
1241 c[0] = cweak[0]+cweak[1]+cweak[2];
1242 c[1] = -(cweak[1]+2.0*cweak[2]);
1246 c[0] = cmedium[0]+cmedium[1]+cmedium[2];
1247 c[1] = -(cmedium[1]+2.0*cmedium[2]);
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];
1258 std::cerr << std::endl <<
">> **ERROR** User Apodization tag " << apodizationTag <<
" unknown, sorry ..." << std::endl;
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));
#define F_APODIZATION_WEAK
Weak apodization (gentle roll-off at edges)
#define F_APODIZATION_STRONG
Strong apodization (heavy roll-off for best frequency resolution)
#define F_APODIZATION_NONE
No apodization (rectangular window)
#define F_APODIZATION_MEDIUM
Medium apodization (moderate roll-off)
#define FOURIER_UNIT_FREQ
Frequency in MHz.
#define FOURIER_UNIT_GAUSS
Magnetic field in Gauss (G)
#define FOURIER_UNIT_CYCLES
Angular frequency in Mc/s (Mega-cycles per second)
#define FOURIER_UNIT_TESLA
Magnetic field in Tesla (T)
PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1)
virtual void CalcRealPhFTDerivative() const
Int_t fMinBin
1st derivative of fRealPh
std::vector< Double_t > fRealPh
original imag Fourier data set
virtual void Init()
keeps the minimum of the entropy/penalty minimization
std::vector< Double_t > fRealPhD
phased imag Fourier data set
virtual Double_t Entropy() const
std::vector< Double_t > fImagPh
phased real Fourier data set
std::vector< Double_t > fReal
virtual void CalcPhasedFT() const
virtual Double_t operator()(const std::vector< Double_t > &) const
Double_t fPh_c1
constant part of the phase dispersion used for the phase correction
Double_t fMin
gamma parameter to balance between entropy and penalty
Double_t fGamma
linear part of the phase dispersion used for the phase correction
Double_t fPh_c0
maximum bin from Fourier range to be used for the phase correction estimate
virtual Double_t Penalty() const
virtual Double_t GetMinimum()
Int_t fMaxBin
minimum bin from Fourier range to be used for the phase correction estimate
virtual Double_t GetPhaseCorrectionParam(UInt_t idx)
std::vector< Double_t > fImag
original real Fourier data set
PFourier(TH1F *data, Int_t unitTag, Double_t startTime=0.0, Double_t endTime=0.0, Bool_t dcCorrected=false, UInt_t zeroPaddingPower=0)
fftw_plan fFFTwPlan
fftw plan (see FFTW3 User Manual)
virtual void ApodizeData(Int_t apodizationTag)
Int_t fApodization
0=none, 1=weak, 2=medium, 3=strong
Bool_t fValid
true = all boundary conditions fullfilled and hence a Fourier transform can be performed.
UInt_t fZeroPaddingPower
power for zero padding, if set < 0 no zero padding will be done
UInt_t fNoOfBins
number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding
virtual TH1F * GetImaginaryFourier(const Double_t scale=1.0)
UInt_t fNoOfData
number of bins in the time interval between fStartTime and fStopTime
Double_t fResolution
Fourier resolution (field, frequency, or angular frequency)
virtual void Transform(UInt_t apodizationTag=0)
virtual TH1F * GetRealFourier(const Double_t scale=1.0)
virtual Double_t GetMaxFreq()
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)
TH1F * fData
data histogram to be Fourier transformed.
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)
Double_t fStartTime
start time of the data histogram
Double_t fTimeResolution
time resolution of the data histogram in (us)
fftw_complex * fIn
real part of the Fourier transform
fftw_complex * fOut
imaginary part of the Fourier transform
Double_t fEndTime
end time of the data histogram
Bool_t fDCCorrected
if true, removed DC offset from signal before Fourier transformation, otherwise not