AliPhysics  4446124 (4446124)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliHFInvMassMultiTrialFit.h
Go to the documentation of this file.
1 #ifndef ALIHFINVMASSMULTITRIALFIT_H
2 #define ALIHFINVMASSMULTITRIALFIT_H
3 /* Copyright(c) 2008-2019, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice */
5 
6 #include <TNamed.h>
7 #include <TString.h>
8 #include <TPad.h>
9 #include <set>
10 #include <vector>
11 
12 class TNtuple;
13 class AliHFInvMassFitter;
14 
16 
18 
19  public:
22 
23  void ConfigureRebinSteps(Int_t nSteps, Int_t* values){
24  if(fRebinSteps) delete [] fRebinSteps;
25  fNumOfRebinSteps=nSteps;
27  for(Int_t ib=0; ib<fNumOfRebinSteps; ib++) fRebinSteps[ib]=values[ib];
28  }
29 
32  }
33 
34  void ConfigureLowLimFitSteps(Int_t nSteps, Double_t* values){
35  if(fLowLimFitSteps) delete [] fLowLimFitSteps;
36  fNumOfLowLimFitSteps=nSteps;
38  for(Int_t ib=0; ib<fNumOfLowLimFitSteps; ib++) fLowLimFitSteps[ib]=values[ib];
39  }
40 
41  void ConfigureUpLimFitSteps(Int_t nSteps, Double_t* values){
42  if(fUpLimFitSteps) delete [] fUpLimFitSteps;
43  fNumOfUpLimFitSteps=nSteps;
45  for(Int_t ib=0; ib<fNumOfUpLimFitSteps; ib++) fUpLimFitSteps[ib]=values[ib];
46  }
47 
48  void ConfigurenSigmaBinCSteps(Int_t nSteps, Double_t* values){
49  if(fnSigmaBinCSteps) delete [] fnSigmaBinCSteps;
50  fNumOfnSigmaBinCSteps=nSteps;
52  for(Int_t ib=0; ib<fNumOfnSigmaBinCSteps; ib++) fnSigmaBinCSteps[ib]=values[ib];
53  }
54 
56  min = fMinYieldGlob;
57  max = fMaxYieldGlob;
58  }
59 
63 
64  void SetSuffixForHistoNames(const Char_t* name){
65  fSuffix=name;
66  }
69 
70 
71  void SetUseExpoBackground(Bool_t opt=kTRUE){fUseExpoBkg=opt;}
72  void SetUseLinBackground(Bool_t opt=kTRUE){fUseLinBkg=opt;}
73  void SetUsePol2Background(Bool_t opt=kTRUE){fUsePol2Bkg=opt;}
74  void SetUsePol3Background(Bool_t opt=kTRUE){fUsePol3Bkg=opt;}
75  void SetUsePol4Background(Bool_t opt=kTRUE){fUsePol4Bkg=opt;}
76  void SetUsePol5Background(Bool_t opt=kTRUE){fUsePol5Bkg=opt;}
79 
82  void SetUseFreeS(Bool_t opt=kTRUE) {fUseFreeS=opt;}
86 
88 
90 
91  Bool_t DoMultiTrials(TH1D* hInvMassHisto, TPad* thePad=0x0);
92  void SaveToRoot(TString fileName, TString option="recreate") const;
93  void DrawHistos(TCanvas* cry) const;
94  void SetTemplatesForReflections(const TH1F *hTemplRefl, const TH1F *hTemplSig);
95 
96  void SetFixRefoS(Float_t refloS){fFixRefloS=refloS;}
97 
99  fUseSecondPeak=kTRUE;
102  }
103 
104  void AddInvMassFitSaveAsFormat(std::string format) { fInvMassFitSaveAsFormats.insert(format); }
106 
107 
110 
111  private:
112 
114  Bool_t DoFitWithPol3Bkg(TH1F* histoToFit, Double_t hmin, Double_t hmax,
115  Int_t theCase);
116 
119 
120  std::set<std::string> fInvMassFitSaveAsFormats;
122  Int_t* fRebinSteps; //[fNumOfRebinSteps] values of rebin
125  Double_t* fLowLimFitSteps; //[fNumOfLowLimFitSteps] values of low limits for fit
127  Double_t* fUpLimFitSteps; //[fNumOfUpLimFitSteps] values of up limits for fit
129  Double_t* fnSigmaBinCSteps; //[fNumOfnSigmaBinCSteps] values of nsigma for bin count
130  Double_t fnSigmaForBkgEval; //value of sigma in which to extract bkg value
131 
151 
157 
159 
161 
170 
175 
182  TH1F** fHistoBkgTrial;
184 
189 
190  TH1F *fhTemplRefl;
191  TH1F *fhTemplSign;
194 
197 
198  std::vector<AliHFInvMassFitter*> fMassFitters;
199 
201  ClassDef(AliHFInvMassMultiTrialFit,2);
202 };
204 
205 #endif
206 
Bool_t fUsePol4Bkg
switch for pol3 background
void SetUseFixSigFixMean(Bool_t opt=kTRUE)
Int_t * fRebinSteps
number of rebin steps
Bool_t fDrawIndividualFits
switch for saving bkg values in nsigma
TH1F ** fHistoChi2Trial
histo with gauss mean from subsamples of trials
TH1F * fHistoBkgInBinEdgesTrialAll
histo with bkg from all trials
Bool_t fUsePowLawTimesExpoBkg
switch for power law background
double Double_t
Definition: External.C:58
Bool_t fUseFixSigFixMean
switch for FixSigFreeMean
Definition: External.C:236
void GetGlobalMinMaxYield(Double_t &min, Double_t &max)
Double_t fSigmaSecondPeak
position of the 2nd peak
void SetSigmaMCVariation(Double_t var=0.15)
TString fileName
Float_t fFixRefloS
template of signal contribution
TString format
file names tag, basically the trigger and calorimeter combination
Double_t fSigmaMCVariation
sigma of D meson peak from MC
Double_t mass
TH1F ** fHistoRawYieldDist
histo with bin counts from all trials
TH2F ** fHistoRawYieldTrialBinC1
histo with bin counts from subsamples of trials
TH1F ** fHistoBkgInBinEdgesTrial
histo with bkg from subsamples of trials
void SetUseLinBackground(Bool_t opt=kTRUE)
Double_t * fnSigmaBinCSteps
number of steps on the bin counting
char Char_t
Definition: External.C:18
TH1F * fHistoBkgTrialAll
histo with chi2 from all trials
AliHFInvMassFitter class for the fit of invariant mass distribution of charm hadrons.
void SetUseFreeS(Bool_t opt=kTRUE)
void ConfigureLowLimFitSteps(Int_t nSteps, Double_t *values)
std::vector< AliHFInvMassFitter * > fMassFitters
maximum yield
TH1F * fHistoChi2TrialAll
histo with gauss mean from all trials
void AddInvMassFitSaveAsFormat(std::string format)
void SetUseExpoBackground(Bool_t opt=kTRUE)
void ConfigureUpLimFitSteps(Int_t nSteps, Double_t *values)
TH1F ** fHistoSignifTrial
histo with chi2 from subsamples of trials
Double_t fMaxYieldGlob
minimum yield
Bool_t fUseFixSigFreeMean
switch for FixedMeanFreeS
TH1F ** fHistoRawYieldDistBinC0
histo with bkg in mass bin edges from subsamples of trials
Bool_t fUsePol3Bkg
switch for pol2 background
AliHFInvMassMultiTrialFit & operator=(const AliHFInvMassMultiTrialFit &source)
Bool_t fUseFixSigUpFreeMean
switch for power law background
TH1F * fHistoMeanTrialAll
histo with gauss sigma from all trials
Double_t fMassD
relative variation of the sigma
void SetUseFixedMeanFreeS(Bool_t opt=kTRUE)
void SetUsePol2Background(Bool_t opt=kTRUE)
Bool_t fFixSigmaSecondPeak
flag to fix the position of the 2nd peak
TH1F * fHistoRawYieldDistBinC0All
histo with bkg in mass bin edges from all trials
TH1F * fHistoSignifTrialAll
histo with chi2 from all trials
Bool_t fUseSecondPeak
switch for FixSigFixMean
Bool_t fFixMassSecondPeak
width of the 2nd peak
int Int_t
Definition: External.C:63
TH1F ** fHistoRawYieldTrial
histo with yield from subsamples of trials
void SetTemplatesForReflections(const TH1F *hTemplRefl, const TH1F *hTemplSig)
Int_t fNumOfLowLimFitSteps
number of steps in the first bin for rebin
float Float_t
Definition: External.C:68
TH1F * fHistoRawYieldDistBinC1All
histo with bin counts from all trials
Bool_t DoFitWithPol3Bkg(TH1F *histoToFit, Double_t hmin, Double_t hmax, Int_t theCase)
TH1F * fHistoRawYieldTrialAll
histo with yield from all trials
Double_t * fUpLimFitSteps
number of steps on the max. mass for fit
Double_t fMassSecondPeak
switch off/on second peak (for D+->KKpi in Ds)
Bool_t DoMultiTrials(TH1D *hInvMassHisto, TPad *thePad=0x0)
Definition: External.C:212
TH1F ** fHistoRawYieldDistBinC1
histo with bin counts from subsamples of trials
void ConfigureRebinSteps(Int_t nSteps, Int_t *values)
TH1F * fHistoSigmaTrialAll
histo with yield from all trials
void SaveToRoot(TString fileName, TString option="recreate") const
TH1F * fhTemplRefl
histo with bin counts from subsamples of trials
TH1F * fHistoRawYieldDistAll
flag for drawing fits
Double_t nsigma
TH2F ** fHistoRawYieldTrialBinC0
histo with bin counts from subsamples of trials
TH1F ** fHistoSigmaTrial
histo with yield from subsamples of trials
void SetSaveBkgValue(Bool_t opt=kTRUE, Double_t nsigma=3)
Bool_t fUsePowLawBkg
switch for pol5 background
Bool_t fUseFreeS
switch for FixSigDownFreeMean
void SetUseFixSigFreeMean(Bool_t opt=kTRUE)
void IncludeSecondGausPeak(Double_t mass, Bool_t fixm, Double_t width, Bool_t fixw)
Int_t fFitOption
name to characterize analysis case
void SetUseFixSigUpFreeMean(Bool_t opt=kTRUE)
Bool_t fUsePol5Bkg
switch for pol4 background
Bool_t fSaveBkgVal
flag to fix the width of the 2nd peak
void SetUsePol5Background(Bool_t opt=kTRUE)
Bool_t fUseFixedMeanFreeS
switch for FreeSigma
TH1F * fhTemplSign
template of reflection contribution
Bool_t fUseLinBkg
switch for exponential background
Bool_t fUseFixSigDownFreeMean
switch for FixSigUpFreeMean
void SetSuffixForHistoNames(const Char_t *name)
Bool_t fUseExpoBkg
LL or chi2 fit.
Int_t fNumOfRebinSteps
saves the invariant mass fit canvases in the file formats listed in this vector (if empty...
void SetUsePol3Background(Bool_t opt=kTRUE)
Bool_t fUsePol2Bkg
switch for linear background
Double_t * fLowLimFitSteps
number of steps on the min. mass for fit
void DrawHistos(TCanvas *cry) const
void SetUseFixSigDownFreeMean(Bool_t opt=kTRUE)
TH1F ** fHistoMeanTrial
histo with gauss sigma from subsamples of trials
void SetUsePowerLawTimesExpoBackground(Bool_t opt=kTRUE)
bool Bool_t
Definition: External.C:53
std::set< std::string > fInvMassFitSaveAsFormats
void ConfigurenSigmaBinCSteps(Int_t nSteps, Double_t *values)
TH1F ** fHistoBkgTrial
histo with chi2 from subsamples of trials
TH2F * fHistoRawYieldTrialBinC0All
histo with bin counts from all trials
void SetDrawIndividualFits(Bool_t opt=kTRUE)
void SetUsePol4Background(Bool_t opt=kTRUE)
TH2F * fHistoRawYieldTrialBinC1All
histo with bin counts from all trials
void SetUsePowerLawBackground(Bool_t opt=kTRUE)