AliPhysics  85eb05e (85eb05e)
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 #include <set>
10 #include <vector>
11 
12 class TNtuple;
13 class AliHFMassFitterVAR;
14 
16 
17 class AliHFMultiTrials : public TNamed {
18 
19  public:
21  virtual ~AliHFMultiTrials();
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  }
78 
85 
87 
89 
90  Bool_t DoMultiTrials(TH1D* hInvMassHisto, TPad* thePad=0x0);
91  void SaveToRoot(TString fileName, TString option="recreate") const;
92  void DrawHistos(TCanvas* cry) const;
93  TH1F* SetTemplateRefl(const TH1F *hTemplRefl);
94 
95  void SetFixRefoS(Float_t refloS){fFixRefloS=refloS;}
96 
97  void AddInvMassFitSaveAsFormat(std::string format) { fInvMassFitSaveAsFormats.insert(format); }
99 
100 
103 
104  private:
105 
107  TH1F* RebinHisto(TH1D* hOrig, Int_t reb, Int_t firstUse) const;
108  void BinCount(TH1F* h, TF1* fB, Int_t rebin, Double_t minMass, Double_t maxMass, Double_t& count, Double_t& ecount) const;
109  Bool_t DoFitWithPol3Bkg(TH1F* histoToFit, Double_t hmin, Double_t hmax,
110  Int_t theCase);
111 
112  AliHFMultiTrials(const AliHFMultiTrials &source);
114 
115  std::set<std::string> fInvMassFitSaveAsFormats;
117  Int_t* fRebinSteps; //[fNumOfRebinSteps] values of rebin
120  Double_t* fLowLimFitSteps; //[fNumOfLowLimFitSteps] values of low limits for fit
122  Double_t* fUpLimFitSteps; //[fNumOfUpLimFitSteps] values of up limits for fit
124  Double_t* fnSigmaBinCSteps; //[fNumOfnSigmaBinCSteps] values of nsigma for bin count
125  Double_t fnSigmaForBkgEval; //value of sigma in which to extract bkg value
126 
146 
148 
150 
159 
162 
169  TH1F** fHistoBkgTrial;
171 
174  TH1F *fhTemplRefl;
177 
180 
181  std::vector<AliHFMassFitterVAR*> fMassFitters;
182 
184  ClassDef(AliHFMultiTrials,5);
185 };
187 
188 #endif
189 
void SetUseFixSigFixMean(Bool_t opt=kTRUE)
void AddInvMassFitSaveAsFormat(std::string format)
void GetGlobalMinMaxYield(Double_t &min, Double_t &max)
TH1F ** fHistoBkgTrial
histo with chi2 from subsamples of trials
double Double_t
Definition: External.C:58
Double_t fMaxYieldGlob
minimum yield
Bool_t fUsePowLawBkg
switch for pol5 background
void SetUseFixSigFreeMean(Bool_t opt=kTRUE)
Definition: External.C:236
void SetUsePol2Background(Bool_t opt=kTRUE)
TH1F * fhTemplRefl
histo with bin counts from subsamples of trials
TH1F * fHistoSigmaTrialAll
histo with yield from all trials
TNtuple * fNtupleMultiTrials
Bool_t fSaveBkgVal
switch for FixSigFixMean
Double_t fMinYieldGlob
tree
TH1F * SetTemplateRefl(const TH1F *hTemplRefl)
TH1F * fHistoRawYieldDistAll
flag for drawing fits
Bool_t fUseLinBkg
switch for exponential background
void SetUsePol4Background(Bool_t opt=kTRUE)
TString fileName
void SetUseLogLikelihoodFit()
TString format
file names tag, basically the trigger and calorimeter combination
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 SetSaveBkgValue(Bool_t opt=kTRUE, Double_t nsigma=3)
void ConfigurenSigmaBinCSteps(Int_t nSteps, Double_t *values)
TString fSuffix
mass of D meson
void DisableInvMassFitSaveAs()
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
void SetUsePowerLawTimesExpoBackground(Bool_t opt=kTRUE)
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
AliHFMassFitterVAR for the fit of invariant mass distribution of charmed mesons.
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
void SetFixRefoS(Float_t refloS)
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 saving bkg values in nsigma
Bool_t fUsePol2Bkg
switch for linear background
Bool_t fUsePol5Bkg
switch for pol4 background
float Float_t
Definition: External.C:68
void SetUseExpoBackground(Bool_t opt=kTRUE)
void SetMass(Double_t mass)
Bool_t fUseFreeS
switch for FixSigDownFreeMean
TH1F ** fHistoRawYieldDistBinC
histo with bkg in mass bin edges from subsamples of trials
Double_t fnSigmaForBkgEval
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 power law background
TH1F * fHistoRawYieldDistBinCAll
histo with bkg in mass bin edges from all trials
Bool_t DoFitWithPol3Bkg(TH1F *histoToFit, Double_t hmin, Double_t hmax, Int_t theCase)
Double_t nsigma
void SetSuffixForHistoNames(const Char_t *name)
TH1F ** fHistoMeanTrial
histo with gauss sigma from subsamples of trials
TH1F * fHistoBkgInBinEdgesTrialAll
histo with bkg from all trials
void SetSigmaMCVariation(Double_t var=0.15)
std::set< std::string > fInvMassFitSaveAsFormats
void SetUsePol5Background(Bool_t opt=kTRUE)
void SetUseFixedMeanFreeS(Bool_t opt=kTRUE)
TH1F * RebinHisto(TH1D *hOrig, Int_t reb, Int_t firstUse) const
Int_t fNumOfRebinSteps
saves the invariant mass fit canvases in the file formats listed in this vector (if empty...
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 SetUsePowerLawBackground(Bool_t opt=kTRUE)
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
TH1F * fHistoBkgTrialAll
histo with chi2 from all trials
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
TH1F ** fHistoBkgInBinEdgesTrial
histo with bkg from subsamples of trials
Double_t maxMass
bool Bool_t
Definition: External.C:53
void SetUseLikelihoodWithWeightsFit()
Bool_t fUseExpoBkg
LL or chi2 fit.
Bool_t DoMultiTrials(TH1D *hInvMassHisto, TPad *thePad=0x0)
std::vector< AliHFMassFitterVAR * > fMassFitters
maximum yield
Double_t * fUpLimFitSteps
number of steps on the max. mass for fit
Float_t fFixRefloS
template of reflection contribution
TH1F * fHistoChi2TrialAll
histo with gauss mean from all trials
void SetNumOfFirstBinSteps(Int_t nfst)
Bool_t fUsePowLawTimesExpoBkg
switch for power law background