AliPhysics  32b88a8 (32b88a8)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliHFMultiTrials.h
Go to the documentation of this file.
1 #ifndef ALIHFMULTITRIALS_H
2 #define ALIHFMULTITRIALS_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 
11 
12 class AliHFMultiTrials : public TNamed {
13 
14  public:
16  virtual ~AliHFMultiTrials();
17 
18  void ConfigureRebinSteps(Int_t nSteps, Int_t* values){
19  if(fRebinSteps) delete [] fRebinSteps;
20  fNumOfRebinSteps=nSteps;
22  for(Int_t ib=0; ib<fNumOfRebinSteps; ib++) fRebinSteps[ib]=values[ib];
23  }
24 
27  }
28 
29  void ConfigureLowLimFitSteps(Int_t nSteps, Double_t* values){
30  if(fLowLimFitSteps) delete [] fLowLimFitSteps;
31  fNumOfLowLimFitSteps=nSteps;
33  for(Int_t ib=0; ib<fNumOfLowLimFitSteps; ib++) fLowLimFitSteps[ib]=values[ib];
34  }
35 
36  void ConfigureUpLimFitSteps(Int_t nSteps, Double_t* values){
37  if(fUpLimFitSteps) delete [] fUpLimFitSteps;
38  fNumOfUpLimFitSteps=nSteps;
40  for(Int_t ib=0; ib<fNumOfUpLimFitSteps; ib++) fUpLimFitSteps[ib]=values[ib];
41  }
42 
43  void ConfigurenSigmaBinCSteps(Int_t nSteps, Double_t* values){
44  if(fnSigmaBinCSteps) delete [] fnSigmaBinCSteps;
45  fNumOfnSigmaBinCSteps=nSteps;
47  for(Int_t ib=0; ib<fNumOfnSigmaBinCSteps; ib++) fnSigmaBinCSteps[ib]=values[ib];
48  }
49 
51  min = fMinYieldGlob;
52  max = fMaxYieldGlob;
53  }
54 
58 
59  void SetSuffixForHistoNames(const Char_t* name){
60  fSuffix=name;
61  }
64  void SetUseExpoBackground(Bool_t opt=kTRUE){fUseExpoBkg=opt;}
65  void SetUseLinBackground(Bool_t opt=kTRUE){fUseLinBkg=opt;}
66  void SetUsePol2Background(Bool_t opt=kTRUE){fUsePol2Bkg=opt;}
67  void SetUsePol3Background(Bool_t opt=kTRUE){fUsePol3Bkg=opt;}
68  void SetUsePol4Background(Bool_t opt=kTRUE){fUsePol4Bkg=opt;}
69  void SetUsePol5Background(Bool_t opt=kTRUE){fUsePol5Bkg=opt;}
70 
73  void SetUseFreeS(Bool_t opt=kTRUE) {fUseFreeS=opt;}
77 
79 
80  Bool_t DoMultiTrials(TH1D* hInvMassHisto, TPad* thePad=0x0);
81  void SaveToRoot(TString fileName, TString option="recreate") const;
82  void DrawHistos(TCanvas* cry) const;
83 
86 
87  private:
88 
90  TH1F* RebinHisto(TH1D* hOrig, Int_t reb, Int_t firstUse) const;
91  void BinCount(TH1F* h, TF1* fB, Int_t rebin, Double_t minMass, Double_t maxMass, Double_t& count, Double_t& ecount) const;
92  Bool_t DoFitWithPol3Bkg(TH1F* histoToFit, Double_t hmin, Double_t hmax,
93  Int_t theCase);
94 
95  AliHFMultiTrials(const AliHFMultiTrials &source);
97 
99  Int_t* fRebinSteps; //[fNumOfRebinSteps] values of rebin
102  Double_t* fLowLimFitSteps; //[fNumOfLowLimFitSteps] values of low limits for fit
104  Double_t* fUpLimFitSteps; //[fNumOfUpLimFitSteps] values of up limits for fit
106  Double_t* fnSigmaBinCSteps; //[fNumOfnSigmaBinCSteps] values of nsigma for bin count
107 
125 
127 
134 
137 
144 
147 
149 
152 
154  ClassDef(AliHFMultiTrials,2);
155 };
157 
158 #endif
159 
void SetUseFixSigFixMean(Bool_t opt=kTRUE)
void GetGlobalMinMaxYield(Double_t &min, Double_t &max)
double Double_t
Definition: External.C:58
Double_t fMaxYieldGlob
minimum yield
void SetUseFixSigFreeMean(Bool_t opt=kTRUE)
Definition: External.C:236
void SetUsePol2Background(Bool_t opt=kTRUE)
TH1F * fHistoSigmaTrialAll
histo with yield from all trials
TNtuple * fNtupleMultiTrials
histo with bin counts from subsamples of trials
Double_t fMinYieldGlob
tree
TH1F * fHistoRawYieldDistAll
flag for drawing fits
Bool_t fUseLinBkg
switch for exponential background
void SetUsePol4Background(Bool_t opt=kTRUE)
TString fileName
void SetUseLogLikelihoodFit()
Double_t mass
Double_t * fLowLimFitSteps
number of steps on the min. mass for fit
Bool_t fUseFixSigFixMean
switch for FixSigFreeMean
char Char_t
Definition: External.C:18
void ConfigurenSigmaBinCSteps(Int_t nSteps, Double_t *values)
TString fSuffix
mass of D meson
TH1F * fHistoMeanTrialAll
histo with gauss sigma from all trials
AliHFMultiTrials & operator=(const AliHFMultiTrials &source)
TH1F ** fHistoRawYieldTrial
histo with yield from subsamples of trials
Double_t fMassD
relative variation of the sigma
Int_t * fRebinSteps
number of rebin steps
TH2F ** fHistoRawYieldTrialBinC
histo with bin counts from subsamples of trials
TH1F ** fHistoSigmaTrial
histo with yield from subsamples of trials
Int_t fNumOfLowLimFitSteps
number of steps in the first bin for rebin
void DrawHistos(TCanvas *cry) const
void SetUseFreeS(Bool_t opt=kTRUE)
TH1F ** fHistoSignifTrial
histo with chi2 from subsamples of trials
TH1F ** fHistoRawYieldDist
histo with bin counts from all trials
Double_t * fnSigmaBinCSteps
number of steps on the bin counting
Int_t fFitOption
name to characterize analysis case
int Int_t
Definition: External.C:63
TH1F * fHistoSignifTrialAll
histo with chi2 from all trials
void ConfigureLowLimFitSteps(Int_t nSteps, Double_t *values)
Bool_t fDrawIndividualFits
switch for FixSigFixMean
Bool_t fUsePol2Bkg
switch for linear background
Bool_t fUsePol5Bkg
switch for pol4 background
void SetUseExpoBackground(Bool_t opt=kTRUE)
void SetMass(Double_t mass)
Bool_t fUseFreeS
switch for FixSigDownFreeMean
TH1F ** fHistoRawYieldDistBinC
histo with chi2 from subsamples of trials
Bool_t fUseFixedMeanFreeS
switch for FreeSigma
void SetDrawIndividualFits(Bool_t opt=kTRUE)
Definition: External.C:212
void SetSigmaGaussMC(Double_t sig)
Bool_t fUseFixSigUpFreeMean
switch for pol5 background
TH1F * fHistoRawYieldDistBinCAll
histo with chi2 from all trials
Bool_t DoFitWithPol3Bkg(TH1F *histoToFit, Double_t hmin, Double_t hmax, Int_t theCase)
void SetSuffixForHistoNames(const Char_t *name)
TH1F ** fHistoMeanTrial
histo with gauss sigma from subsamples of trials
void SetSigmaMCVariation(Double_t var=0.15)
void SetUsePol5Background(Bool_t opt=kTRUE)
void SetUseFixedMeanFreeS(Bool_t opt=kTRUE)
TH1F * RebinHisto(TH1D *hOrig, Int_t reb, Int_t firstUse) const
void BinCount(TH1F *h, TF1 *fB, Int_t rebin, Double_t minMass, Double_t maxMass, Double_t &count, Double_t &ecount) const
void SetUseLinBackground(Bool_t opt=kTRUE)
TH2F * fHistoRawYieldTrialBinCAll
histo with bin counts from all trials
Double_t fSigmaMCVariation
sigma of D meson peak from MC
void SetUseFixSigUpFreeMean(Bool_t opt=kTRUE)
void SetUseFixSigDownFreeMean(Bool_t opt=kTRUE)
Bool_t fUsePol3Bkg
switch for pol2 background
Bool_t fUseFixSigDownFreeMean
switch for FixSigUpFreeMean
Bool_t fUsePol4Bkg
switch for pol3 background
Double_t minMass
Int_t rebin
TH1F ** fHistoChi2Trial
histo with gauss mean from subsamples of trials
void SetUsePol3Background(Bool_t opt=kTRUE)
TH1F * fHistoRawYieldTrialAll
histo with yield from all trials
void ConfigureRebinSteps(Int_t nSteps, Int_t *values)
void ConfigureUpLimFitSteps(Int_t nSteps, Double_t *values)
void SaveToRoot(TString fileName, TString option="recreate") const
Bool_t fUseFixSigFreeMean
switch for FixedMeanFreeS
Double_t maxMass
virtual ~AliHFMultiTrials()
bool Bool_t
Definition: External.C:53
Bool_t fUseExpoBkg
LL or chi2 fit.
Bool_t DoMultiTrials(TH1D *hInvMassHisto, TPad *thePad=0x0)
Double_t * fUpLimFitSteps
number of steps on the max. mass for fit
TH1F * fHistoChi2TrialAll
histo with gauss mean from all trials
void SetNumOfFirstBinSteps(Int_t nfst)