]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAODPidHF.h
PWGHF converted to native cmake
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliAODPidHF.h
1 #ifndef ALIAODPIDHF_H
2 #define ALIAODPIDHF_H
3
4 /* Copyright(c) 1998-2006, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 /* $Id$ */
8
9 //***********************************************************
10 //// Class AliAODPidHF
11 //// class for PID with AliAODRecoDecayHF
12 //// Authors: D. Caffarri caffarri@pd.infn.it, A.Dainese andrea.dainese@pd.infn.it, S. Dash dash@to.infn.it, F. Prino prino@to.infn.it, R. Romita r.romita@gsi.de, Y. Wang yifei@pi0.physi.uni-heidelberg.de,  J. van der Maarel j.vandermaarel@cern.ch
13 ////***********************************************************
14
15 #include <TString.h>
16 #include <TH1F.h>
17 #include <TObject.h>
18 #include "AliAODTrack.h"
19 #include "AliPIDResponse.h"
20 #include "AliPIDCombined.h"
21 #include "AliPID.h"
22
23 class AliAODPidHF : public TObject{
24   
25 public:
26   
27   enum ECombDetectors {
28     kTPC,
29     kTOF,
30     kTPCTOF,
31     kTPCITS
32   };
33   
34   AliAODPidHF();
35   AliAODPidHF(const AliAODPidHF& pid);
36   virtual ~AliAODPidHF();
37   
38   //Setters
39   void SetSigma(Double_t *sigma){
40     for(Int_t i=0; i<fnNSigma; i++) fnSigma[i]=sigma[i];
41   }
42   void SetSigma(Int_t idet,Double_t sigma){fnSigma[idet]=sigma;return;}
43   void SetSigmaForTPC(Double_t *sigma){for(Int_t i=0;i<3;i++) fnSigma[i]=sigma[i];return;}
44   void SetSigmaForTPCCompat(Double_t sigma){fnSigmaCompat[0]=sigma;return;}
45   void SetSigmaForTOFCompat(Double_t sigma){fnSigmaCompat[1]=sigma;return;}
46   void SetSigmaForTOF(Double_t sigma){fnSigma[3]=sigma;return;}
47   void SetSigmaForITS(Double_t sigma){fnSigma[4]=sigma;return;}
48   void SetTofSigma(Double_t sigma){fTOFSigma=sigma;return;}
49   
50   void SetCutOnTOFmismatchProb(Double_t cut=0.01){fCutTOFmismatch=cut;}
51   void DisableCutOnTOFmismatchProb(){fCutTOFmismatch=999.;}
52   
53   void SetMinNClustersTPCPID(Int_t minc) {fMinNClustersTPCPID=minc;}
54   
55   void SetCombinednSigmaCutForPiKP(Float_t sigpi, Float_t sigk, Float_t sigp){
56     fMaxnSigmaCombined[0]=sigpi;
57     fMaxnSigmaCombined[1]=sigk;
58     fMaxnSigmaCombined[2]=sigp;
59   }
60   void SetTPCnSigmaRangeForPions(Float_t smin, Float_t smax){
61     fMinnSigmaTPC[0]=smin;
62     fMaxnSigmaTPC[0]=smax;
63   }
64   void SetTOFnSigmaRangeForPions(Float_t smin, Float_t smax){
65     fMinnSigmaTOF[0]=smin;
66     fMaxnSigmaTOF[0]=smax;
67   }
68   void SetTPCnSigmaRangeForKaons(Float_t smin, Float_t smax){
69     fMinnSigmaTPC[1]=smin;
70     fMaxnSigmaTPC[1]=smax;
71   }
72   void SetTOFnSigmaRangeForKaons(Float_t smin, Float_t smax){
73     fMinnSigmaTOF[1]=smin;
74     fMaxnSigmaTOF[1]=smax;
75   }
76   void SetTPCnSigmaRangeForProtons(Float_t smin, Float_t smax){
77     fMinnSigmaTPC[2]=smin;
78     fMaxnSigmaTPC[2]=smax;
79   }
80   void SetTOFnSigmaRangeForProtons(Float_t smin, Float_t smax){
81     fMinnSigmaTOF[2]=smin;
82     fMaxnSigmaTOF[2]=smax;
83   }
84   
85   void SetPriors(Double_t *priors){fPriors=priors;return;}
86   void SetPLimit(Double_t *plim){for(Int_t i=0;i<fnPLimit;i++) fPLimit[i]=plim[i];return;}
87   void SetPLimit(Double_t *plim,Int_t npLim){fnPLimit=npLim;for(Int_t i=0;i<fnPLimit;i++) fPLimit[i]=plim[i];return;}
88   void SetAsym(Bool_t asym){fAsym=asym;return;}
89   void SetUseAsymmnSigmaTOF(Double_t nsmin, Double_t nsmax, Double_t nscompmin, Double_t nscompmax){
90     fUseAsymTOF=kTRUE;
91     fLownSigmaTOF=nsmin; fUpnSigmaTOF=nsmax;
92     fLownSigmaCompatTOF=nscompmin; fUpnSigmaCompatTOF=nscompmax;
93   }
94   void SetTPC(Bool_t tpc){fTPC=tpc;return;}
95   void SetTOF(Bool_t tof){fTOF=tof;return;}
96   void SetITS(Bool_t its){fITS=its;return;}
97   void SetTRD(Bool_t trd){fTRD=trd;return;}
98   void SetMatch(Int_t match){fMatch=match;return;}
99   void SetForceTOFforKaons(Bool_t forceTOF){fForceTOFforKaons=forceTOF;return;}
100   void SetCompat(Bool_t comp){fCompat=comp;return;}
101   void SetMC(Bool_t mc){fMC=mc;return;}
102   void SetMClowenpp2011(Bool_t mc){fMCLowEn2011=mc;return;}
103   void SetOnePad(Bool_t onepad){fOnePad=onepad;return;}
104   void SetppLowEn2011(Bool_t opt){fppLowEn2011=opt;return;}
105   void SetPbPb(Bool_t pbpb){fPbPb=pbpb;return;}
106   void SetPCompatTOF(Double_t pTOF){fPCompatTOF=pTOF;return;}
107   void SetTOFdecide(Bool_t tOFdecide){fTOFdecide=tOFdecide;return;}
108   void SetOldPid(Bool_t oldPid){fOldPid=oldPid;return;}
109   void SetPtThresholdTPC(Double_t ptThresholdTPC){fPtThresholdTPC=ptThresholdTPC;return;}
110   void SetMaxTrackMomForCombinedPID(Double_t mom){fMaxTrackMomForCombinedPID=mom;}
111   void SetPidResponse(AliPIDResponse *pidResp) {fPidResponse=pidResp;return;}
112   void SetCombDetectors(ECombDetectors pidComb) {
113     fCombDetectors=pidComb;
114   }
115   void SetPionPriorHisto(TH1F* histo){
116     if(fPriorsH[AliPID::kPion]) delete fPriorsH[AliPID::kPion];
117     fPriorsH[AliPID::kPion] = new TH1F(*histo);
118   }
119   void SetKaonPriorHisto(TH1F* histo){
120     if(fPriorsH[AliPID::kKaon]) delete fPriorsH[AliPID::kKaon];
121     fPriorsH[AliPID::kKaon] = new TH1F(*histo);
122   }
123   void SetProtonPriorHisto(TH1F* histo){
124     if(fPriorsH[AliPID::kProton]) delete fPriorsH[AliPID::kProton];
125     fPriorsH[AliPID::kProton] = new TH1F(*histo);
126   }
127   
128   
129   //Getters
130   
131   Int_t GetnSigmaTPC(AliAODTrack *track, Int_t species, Double_t &sigma) const;
132   Int_t GetnSigmaTOF(AliAODTrack *track, Int_t species, Double_t &sigma) const;
133   Int_t GetnSigmaITS(AliAODTrack *track, Int_t species, Double_t &sigma) const;
134   Double_t GetSigma(Int_t idet) const{return fnSigma[idet];}
135   Double_t GetTofSigma() const{return fTOFSigma;}
136   //void GetPriors(Double_t *priors) const{priors=fPriors;return;}
137   //void GetPLimit(Double_t *plim) const{plim=fPLimit;}
138   void GetPriors(Double_t *priors) const{for(Int_t i=0;i<fnPriors;i++){priors[i]=fPriors[i];}return;}
139   void GetPLimit(Double_t *plim) const{for(Int_t i=0;i<fnPLimit;i++){plim[i]=fPLimit[i];}return;}
140   Bool_t GetAsym() const{return fAsym;}
141   Bool_t GetTPC() const{return fTPC;}
142   Bool_t GetTOF() const{return fTOF;}
143   Bool_t GetITS() const{return fITS;}
144   Bool_t GetTRD() const{return fTRD;}
145   Int_t GetMatch() const{return fMatch;}
146   Bool_t GetForceTOFforKaons() const{return fForceTOFforKaons;}
147   Bool_t GetCompat() const{return fCompat;}
148   Bool_t GetMC() const{return fMC;}
149   Bool_t GetOnePad() const{return fOnePad;}
150   Bool_t GetppLowEn2011() const {return fppLowEn2011;}
151   Bool_t GetMCLowEn2011() const {return fMCLowEn2011;}
152   Bool_t GetPbPb() const{return fPbPb;}
153   Bool_t GetTOFdecide() const{return fTOFdecide;}
154   Double_t GetPCompatTOF() const{return fPCompatTOF;}
155   Double_t GetnSigmaCompatTPC() const{return fnSigmaCompat[0];}
156   Double_t GetnSigmaCompatTOF() const{return fnSigmaCompat[1];}
157   Bool_t GetOldPid(){return fOldPid;}
158   Double_t GetPtThresholdTPC(){return fPtThresholdTPC;}
159   Double_t GetMaxTrackMomForCombinedPID(){return fMaxTrackMomForCombinedPID;}
160   AliPIDResponse *GetPidResponse() const {return fPidResponse;}
161   AliPIDCombined *GetPidCombined() const {return fPidCombined;}
162   ECombDetectors GetCombDetectors() const {
163     return fCombDetectors;
164   }
165   Bool_t GetUseCombined() {return fUseCombined;}
166   Bool_t GetDefaultPriors() {return fDefaultPriors;}
167   
168   Int_t RawSignalPID (AliAODTrack *track, TString detector) const;
169   Bool_t IsKaonRaw (AliAODTrack *track, TString detector) const;
170   Bool_t IsPionRaw (AliAODTrack *track, TString detector) const;
171   Bool_t IsProtonRaw (AliAODTrack *track, TString detector) const;
172   Bool_t IsElectronRaw (AliAODTrack *track, TString detector) const;
173   void CombinedProbability(AliAODTrack *track,Bool_t *type) const; //0 = pion, 1 = kaon, 2 = proton
174   Bool_t CheckStatus(AliAODTrack *track,TString detectors) const;
175   
176   Bool_t CheckITSPIDStatus(AliAODTrack *track) const;
177   Bool_t CheckTPCPIDStatus(AliAODTrack *track) const;
178   Bool_t CheckTOFPIDStatus(AliAODTrack *track) const;
179   Bool_t CheckTRDPIDStatus(AliAODTrack *track) const;
180   
181   Bool_t TPCRawAsym(AliAODTrack* track,Int_t specie) const;
182   Int_t MatchTPCTOF(AliAODTrack *track,Int_t specie);
183   
184   Int_t MakeRawPid(AliAODTrack *track,Int_t specie); //general method to perform PID using raw signals
185   
186   Bool_t IsTOFPiKexcluded(AliAODTrack *track,Double_t nsigmaK);
187   
188   Bool_t IsExcluded(AliAODTrack *track, Int_t labelTrack, Double_t nsigmaCut, TString detector);
189   
190   void GetTPCBetheBlochParams(Double_t alephParameters[5]) const;
191   void SetBetheBloch();
192   // method for AliPIDCombined object
193   void SetSelectedSpecies(Int_t ispecies = AliPID::kSPECIES){GetPidCombined()->SetSelectedSpecies(ispecies);};
194   void SetPriorDistribution(AliPID::EParticleType type,TH1F *prior);
195   void DrawPrior(AliPID::EParticleType type);
196   void SetPriorsHistos(TString priorFileName);
197   void SetUpCombinedPID();
198   void SetUseCombined(Bool_t useCombined=kTRUE) {fUseCombined=useCombined;}
199   void SetUseDefaultPriors(Bool_t defaultP)         {fDefaultPriors=defaultP;}
200   Int_t ApplyPidTPCRaw(AliAODTrack *track,Int_t specie) const;
201   Int_t ApplyPidTOFRaw(AliAODTrack *track,Int_t specie) const;
202   Int_t ApplyPidITSRaw(AliAODTrack *track,Int_t specie) const;
203   Int_t ApplyTOFCompatibilityBand(AliAODTrack *track,Int_t specie) const;
204   
205   void PrintAll() const;
206
207   //Assymetric PID using histograms
208   void SetIdBand(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, TH1F *min, TH1F *max);
209   void SetIdBand(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, TF1 *min, TF1 *max);
210   void SetCompBand(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, TH1F *min, TH1F *max);
211   void SetCompBand(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, TF1 *min, TF1 *max);
212   Bool_t CheckDetectorPIDStatus(AliPIDResponse::EDetector detector, AliAODTrack *track);
213   Float_t NumberOfSigmas(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, AliAODTrack *track);
214   Int_t CheckBands(AliPID::EParticleType specie, AliPIDResponse::EDetector detector, AliAODTrack *track);
215   TF1 *GetIdBandMin(AliPID::EParticleType specie, AliPIDResponse::EDetector detector) { return fIdBandMin[((int) specie)][((int) detector)]; }
216   TF1 *GetIdBandMax(AliPID::EParticleType specie, AliPIDResponse::EDetector detector) { return fIdBandMax[((int) specie)][((int) detector)]; }
217   TF1 *GetCompBandMin(AliPID::EParticleType specie, AliPIDResponse::EDetector detector) { return fCompBandMin[((int) specie)][((int) detector)]; }
218   TF1 *GetCompBandMax(AliPID::EParticleType specie, AliPIDResponse::EDetector detector) { return fCompBandMax[((int) specie)][((int) detector)]; }
219
220   //Some suggested asymmetric PID
221   void SetShiftedAsymmetricPID();
222   void SetIdAsymmetricPID();
223   void SetIdCompAsymmetricPID();
224   
225 protected:
226   
227   
228 private:
229
230   AliAODPidHF& operator=(const AliAODPidHF& pid);
231
232   Int_t fnNSigma; // number of sigmas
233   Double_t *fnSigma; // [fnNSigma], sigma for the raw signal PID: 0-2 for TPC, 3 for TOF, 4 for ITS
234   Double_t fTOFSigma; // TOF precision
235   Double_t fCutTOFmismatch; // Cut of TOF mismatch probability
236   UInt_t fMinNClustersTPCPID;       // Minimum TPC PID clusters cut
237   Int_t fnPriors; //number of priors
238   Double_t *fPriors; // [fnPriors], set of priors
239   Int_t fnPLimit; //number of Plimit
240   Double_t *fPLimit; // [fnPLimit], limit of p intervals for asimmetric PID: fPLimit<p[0], fPLimit[0]<p<fPLimit[1], p>fPLimit[1]
241   Bool_t fAsym; // asimmetric PID required
242   Bool_t fTPC; //switch to include or exclude TPC
243   Bool_t fTOF; // switch to include or exclude TOF
244   Bool_t fITS; //switch to include or exclude ITS
245   Bool_t fTRD; // switch to include or exclude TRD
246   Int_t fMatch; //switch to combine the info from more detectors: 1 = || , 2 = &, 3 = p region
247   Bool_t fForceTOFforKaons; // force TOF for kaons in mode fMatch==5
248   Bool_t fCompat; // compatibility region : useful only if fMatch=1
249   Double_t fPCompatTOF; //  compatibility p limit for TOF
250   Bool_t fUseAsymTOF; // flag for using asymmetrig nSigmaCut in TOF for fMatch==1
251   Double_t fLownSigmaTOF;  // lower nsigma TOF (for fUseAsymTOF)
252   Double_t fUpnSigmaTOF;  // upper nsigma TOF (for fUseAsymTOF)
253   Double_t fLownSigmaCompatTOF;  // lower nsigma TOF (for fUseAsymTOF)
254   Double_t fUpnSigmaCompatTOF;  // upper nsigma TOF (for fUseAsymTOF)
255   Int_t fnNSigmaCompat; // number of sigmas
256   Double_t *fnSigmaCompat; //[fnNSigmaCompat]  0: n sigma for TPC compatibility band, 1: for TOF
257   Double_t fMaxnSigmaCombined[3];  // nSigma cut for pi,K,p (TPC^2+TOF^2)
258   Double_t fMinnSigmaTPC[3]; //min. of nSigma range for pi,K,p in TPC (match==5)
259   Double_t fMaxnSigmaTPC[3]; //max. of nSigma range for pi,K,p in TPC (match==5)
260   Double_t fMinnSigmaTOF[3]; //min. of nSigma range for pi,K,p in TOF (match==5)
261   Double_t fMaxnSigmaTOF[3]; //max. of nSigma range for pi,K,p in TOF (match==5)
262   Bool_t fMC; // MC(kTRUE) or real data (kFALSE, default option)
263   Bool_t fOnePad; //  real data with one pad clusters
264   Bool_t fMCLowEn2011; //  MC for low energy MC
265   Bool_t fppLowEn2011; //  Data for low energy pp 2011
266   Bool_t fPbPb; //  real data PbPb
267   Bool_t fTOFdecide; //  real data PbPb
268   Bool_t fOldPid; //  old PID method implemented
269   Double_t fPtThresholdTPC; //  pT threshold to use TPC PID
270   Double_t fMaxTrackMomForCombinedPID; // momentum threshold to use PID
271   AliPIDResponse *fPidResponse; //! pid response
272   AliPIDCombined* fPidCombined; //! combined PID object
273   
274   AliTPCPIDResponse* fTPCResponse; //! TPC response
275   
276   TH1F* fPriorsH[AliPID::kSPECIES]; // priors histos
277   ECombDetectors fCombDetectors; // detectors to be involved for combined PID
278   Bool_t fUseCombined; // detectors to be involved for combined PID
279   Bool_t fDefaultPriors; // use default priors for combined PID
280
281   //Storage of identification/compatibility band for different species and detectors:
282   TF1 *fIdBandMin[AliPID::kSPECIES][4];
283   TF1 *fIdBandMax[AliPID::kSPECIES][4];
284   TF1 *fCompBandMin[AliPID::kSPECIES][4];
285   TF1 *fCompBandMax[AliPID::kSPECIES][4];
286
287   ClassDef(AliAODPidHF,24) // AliAODPid for heavy flavor PID
288
289 };
290
291 struct HistFunc {
292    HistFunc(TH1F *f): fHist(f) {}
293    double operator() (double *x, double * ) const {
294       TAxis *axis = fHist->GetXaxis();
295       Int_t bin = axis->FindBin(x[0]);
296       if (x[0] == axis->GetXmax()) {
297         bin = axis->GetNbins();
298       }
299       return fHist->GetBinContent(bin);
300
301    }
302    TH1F *fHist;
303 };
304
305 #endif
306