]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAODPidHF.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[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
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     kTPCAndTOF
33   };
34   
35   AliAODPidHF();
36   AliAODPidHF(const AliAODPidHF& pid);
37   AliAODPidHF& operator=(const AliAODPidHF& pid);
38   virtual ~AliAODPidHF();
39   
40   //Setters
41   void SetSigma(Double_t *sigma){
42     for(Int_t i=0; i<fnNSigma; i++) fnSigma[i]=sigma[i];
43   }
44   void SetSigma(Int_t idet,Double_t sigma){fnSigma[idet]=sigma;return;}
45   void SetSigmaForTPC(Double_t *sigma){for(Int_t i=0;i<3;i++) fnSigma[i]=sigma[i];return;}
46   void SetSigmaForTPCCompat(Double_t sigma){fnSigmaCompat[0]=sigma;return;}
47   void SetSigmaForTOFCompat(Double_t sigma){fnSigmaCompat[1]=sigma;return;}
48   void SetSigmaForTOF(Double_t sigma){fnSigma[3]=sigma;return;}
49   void SetSigmaForITS(Double_t sigma){fnSigma[4]=sigma;return;}
50   void SetTofSigma(Double_t sigma){fTOFSigma=sigma;return;}
51   
52   void SetCutOnTOFmismatchProb(Double_t cut=0.01){fCutTOFmismatch=cut;}
53   void DisableCutOnTOFmismatchProb(){fCutTOFmismatch=999.;}
54   
55   void SetMinNClustersTPCPID(Int_t minc) {fMinNClustersTPCPID=minc;}
56   
57   void SetCombinednSigmaCutForPiKP(Float_t sigpi, Float_t sigk, Float_t sigp){
58     fMaxnSigmaCombined[0]=sigpi;
59     fMaxnSigmaCombined[1]=sigk;
60     fMaxnSigmaCombined[2]=sigp;
61   }
62   void SetTPCnSigmaRangeForPions(Float_t smin, Float_t smax){
63     fMinnSigmaTPC[0]=smin;
64     fMaxnSigmaTPC[0]=smax;
65   }
66   void SetTOFnSigmaRangeForPions(Float_t smin, Float_t smax){
67     fMinnSigmaTOF[0]=smin;
68     fMaxnSigmaTOF[0]=smax;
69   }
70   void SetTPCnSigmaRangeForKaons(Float_t smin, Float_t smax){
71     fMinnSigmaTPC[1]=smin;
72     fMaxnSigmaTPC[1]=smax;
73   }
74   void SetTOFnSigmaRangeForKaons(Float_t smin, Float_t smax){
75     fMinnSigmaTOF[1]=smin;
76     fMaxnSigmaTOF[1]=smax;
77   }
78   void SetTPCnSigmaRangeForProtons(Float_t smin, Float_t smax){
79     fMinnSigmaTPC[2]=smin;
80     fMaxnSigmaTPC[2]=smax;
81   }
82   void SetTOFnSigmaRangeForProtons(Float_t smin, Float_t smax){
83     fMinnSigmaTOF[2]=smin;
84     fMaxnSigmaTOF[2]=smax;
85   }
86   
87   void SetPriors(Double_t *priors){fPriors=priors;return;}
88   void SetPLimit(Double_t *plim){for(Int_t i=0;i<fnPLimit;i++) fPLimit[i]=plim[i];return;}
89   void SetPLimit(Double_t *plim,Int_t npLim){fnPLimit=npLim;for(Int_t i=0;i<fnPLimit;i++) fPLimit[i]=plim[i];return;}
90   void SetAsym(Bool_t asym){fAsym=asym;return;}
91   void SetUseAsymmnSigmaTOF(Double_t nsmin, Double_t nsmax, Double_t nscompmin, Double_t nscompmax){
92     fUseAsymTOF=kTRUE;
93     fLownSigmaTOF=nsmin; fUpnSigmaTOF=nsmax;
94     fLownSigmaCompatTOF=nscompmin; fUpnSigmaCompatTOF=nscompmax;
95   }
96   void SetTPC(Bool_t tpc){fTPC=tpc;return;}
97   void SetTOF(Bool_t tof){fTOF=tof;return;}
98   void SetITS(Bool_t its){fITS=its;return;}
99   void SetTRD(Bool_t trd){fTRD=trd;return;}
100   void SetMatch(Int_t match){fMatch=match;return;}
101   void SetForceTOFforKaons(Bool_t forceTOF){fForceTOFforKaons=forceTOF;return;}
102   void SetCompat(Bool_t comp){fCompat=comp;return;}
103   void SetMC(Bool_t mc){fMC=mc;return;}
104   void SetMClowenpp2011(Bool_t mc){fMCLowEn2011=mc;return;}
105   void SetOnePad(Bool_t onepad){fOnePad=onepad;return;}
106   void SetppLowEn2011(Bool_t opt){fppLowEn2011=opt;return;}
107   void SetPbPb(Bool_t pbpb){fPbPb=pbpb;return;}
108   void SetPCompatTOF(Double_t pTOF){fPCompatTOF=pTOF;return;}
109   void SetTOFdecide(Bool_t tOFdecide){fTOFdecide=tOFdecide;return;}
110   void SetOldPid(Bool_t oldPid){fOldPid=oldPid;return;}
111   void SetPtThresholdTPC(Double_t ptThresholdTPC){fPtThresholdTPC=ptThresholdTPC;return;}
112   void SetMaxTrackMomForCombinedPID(Double_t mom){fMaxTrackMomForCombinedPID=mom;}
113   void SetPidResponse(AliPIDResponse *pidResp) {fPidResponse=pidResp;return;}
114   void SetCombDetectors(ECombDetectors pidComb) {
115     fCombDetectors=pidComb;
116   }
117   void SetPionPriorHisto(TH1F* histo){
118     if(fPriorsH[AliPID::kPion]) delete fPriorsH[AliPID::kPion];
119     fPriorsH[AliPID::kPion] = new TH1F(*histo);
120   }
121   void SetKaonPriorHisto(TH1F* histo){
122     if(fPriorsH[AliPID::kKaon]) delete fPriorsH[AliPID::kKaon];
123     fPriorsH[AliPID::kKaon] = new TH1F(*histo);
124   }
125   void SetProtonPriorHisto(TH1F* histo){
126     if(fPriorsH[AliPID::kProton]) delete fPriorsH[AliPID::kProton];
127     fPriorsH[AliPID::kProton] = new TH1F(*histo);
128   }
129   
130   
131   //Getters
132   
133   Int_t GetnSigmaTPC(AliAODTrack *track, Int_t species, Double_t &sigma) const;
134   Int_t GetnSigmaTOF(AliAODTrack *track, Int_t species, Double_t &sigma) const;
135   Int_t GetnSigmaITS(AliAODTrack *track, Int_t species, Double_t &sigma) const;
136   Double_t GetSigma(Int_t idet) const{return fnSigma[idet];}
137   Double_t GetTofSigma() const{return fTOFSigma;}
138   //void GetPriors(Double_t *priors) const{priors=fPriors;return;}
139   //void GetPLimit(Double_t *plim) const{plim=fPLimit;}
140   void GetPriors(Double_t *priors) const{for(Int_t i=0;i<fnPriors;i++){priors[i]=fPriors[i];}return;}
141   void GetPLimit(Double_t *plim) const{for(Int_t i=0;i<fnPLimit;i++){plim[i]=fPLimit[i];}return;}
142   Bool_t GetAsym() const{return fAsym;}
143   Bool_t GetTPC() const{return fTPC;}
144   Bool_t GetTOF() const{return fTOF;}
145   Bool_t GetITS() const{return fITS;}
146   Bool_t GetTRD() const{return fTRD;}
147   Int_t GetMatch() const{return fMatch;}
148   Bool_t GetForceTOFforKaons() const{return fForceTOFforKaons;}
149   Bool_t GetCompat() const{return fCompat;}
150   Bool_t GetMC() const{return fMC;}
151   Bool_t GetOnePad() const{return fOnePad;}
152   Bool_t GetppLowEn2011() const {return fppLowEn2011;}
153   Bool_t GetMCLowEn2011() const {return fMCLowEn2011;}
154   Bool_t GetPbPb() const{return fPbPb;}
155   Bool_t GetTOFdecide() const{return fTOFdecide;}
156   Double_t GetPCompatTOF() const{return fPCompatTOF;}
157   Double_t GetnSigmaCompatTPC() const{return fnSigmaCompat[0];}
158   Double_t GetnSigmaCompatTOF() const{return fnSigmaCompat[1];}
159   Bool_t GetOldPid(){return fOldPid;}
160   Double_t GetPtThresholdTPC(){return fPtThresholdTPC;}
161   Double_t GetMaxTrackMomForCombinedPID(){return fMaxTrackMomForCombinedPID;}
162   AliPIDResponse *GetPidResponse() const {return fPidResponse;}
163   AliPIDCombined *GetPidCombined() const {return fPidCombined;}
164   ECombDetectors GetCombDetectors() const {
165     return fCombDetectors;
166   }
167   Bool_t GetUseCombined() {return fUseCombined;}
168   Bool_t GetDefaultPriors() {return fDefaultPriors;}
169   
170   Int_t RawSignalPID (AliAODTrack *track, TString detector) const;
171   Bool_t IsKaonRaw (AliAODTrack *track, TString detector) const;
172   Bool_t IsPionRaw (AliAODTrack *track, TString detector) const;
173   Bool_t IsProtonRaw (AliAODTrack *track, TString detector) const;
174   Bool_t IsElectronRaw (AliAODTrack *track, TString detector) const;
175   void CombinedProbability(AliAODTrack *track,Bool_t *type) const; //0 = pion, 1 = kaon, 2 = proton
176   Bool_t CheckStatus(AliAODTrack *track,TString detectors) const;
177   
178   Bool_t CheckITSPIDStatus(AliAODTrack *track) const;
179   Bool_t CheckTPCPIDStatus(AliAODTrack *track) const;
180   Bool_t CheckTOFPIDStatus(AliAODTrack *track) const;
181   Bool_t CheckTRDPIDStatus(AliAODTrack *track) const;
182   
183   Bool_t TPCRawAsym(AliAODTrack* track,Int_t specie) const;
184   Int_t MatchTPCTOF(AliAODTrack *track,Int_t specie);
185   
186   Int_t MakeRawPid(AliAODTrack *track,Int_t specie); //general method to perform PID using raw signals
187   
188   Bool_t IsTOFPiKexcluded(AliAODTrack *track,Double_t nsigmaK);
189   
190   Bool_t IsExcluded(AliAODTrack *track, Int_t labelTrack, Double_t nsigmaCut, TString detector);
191   
192   void GetTPCBetheBlochParams(Double_t alephParameters[5]) const;
193   void SetBetheBloch();
194   // method for AliPIDCombined object
195   void SetSelectedSpecies(Int_t ispecies = AliPID::kSPECIES){GetPidCombined()->SetSelectedSpecies(ispecies);};
196   void SetPriorDistribution(AliPID::EParticleType type,TH1F *prior);
197   void DrawPrior(AliPID::EParticleType type);
198   void SetPriorsHistos(TString priorFileName);
199   void SetUpCombinedPID();
200   void SetUseCombined(Bool_t useCombined=kTRUE) {fUseCombined=useCombined;}
201   void SetUseDefaultPriors(Bool_t defaultP)         {fDefaultPriors=defaultP;}
202   Int_t ApplyPidTPCRaw(AliAODTrack *track,Int_t specie) const;
203   Int_t ApplyPidTOFRaw(AliAODTrack *track,Int_t specie) const;
204   Int_t ApplyPidITSRaw(AliAODTrack *track,Int_t specie) const;
205   Int_t ApplyTOFCompatibilityBand(AliAODTrack *track,Int_t specie) const;
206   
207   void PrintAll() const;
208   
209 protected:
210   
211   
212 private:
213   Int_t fnNSigma; // number of sigmas
214   Double_t *fnSigma; // [fnNSigma], sigma for the raw signal PID: 0-2 for TPC, 3 for TOF, 4 for ITS
215   Double_t fTOFSigma; // TOF precision
216   Double_t fCutTOFmismatch; // Cut of TOF mismatch probability
217   UInt_t fMinNClustersTPCPID;       // Minimum TPC PID clusters cut
218   Int_t fnPriors; //number of priors
219   Double_t *fPriors; // [fnPriors], set of priors
220   Int_t fnPLimit; //number of Plimit
221   Double_t *fPLimit; // [fnPLimit], limit of p intervals for asimmetric PID: fPLimit<p[0], fPLimit[0]<p<fPLimit[1], p>fPLimit[1]
222   Bool_t fAsym; // asimmetric PID required
223   Bool_t fTPC; //switch to include or exclude TPC
224   Bool_t fTOF; // switch to include or exclude TOF
225   Bool_t fITS; //switch to include or exclude ITS
226   Bool_t fTRD; // switch to include or exclude TRD
227   Int_t fMatch; //switch to combine the info from more detectors: 1 = || , 2 = &, 3 = p region
228   Bool_t fForceTOFforKaons; // force TOF for kaons in mode fMatch==5
229   Bool_t fCompat; // compatibility region : useful only if fMatch=1
230   Double_t fPCompatTOF; //  compatibility p limit for TOF
231   Bool_t fUseAsymTOF; // flag for using asymmetrig nSigmaCut in TOF for fMatch==1
232   Double_t fLownSigmaTOF;  // lower nsigma TOF (for fUseAsymTOF)
233   Double_t fUpnSigmaTOF;  // upper nsigma TOF (for fUseAsymTOF)
234   Double_t fLownSigmaCompatTOF;  // lower nsigma TOF (for fUseAsymTOF)
235   Double_t fUpnSigmaCompatTOF;  // upper nsigma TOF (for fUseAsymTOF)
236   Int_t fnNSigmaCompat; // number of sigmas
237   Double_t *fnSigmaCompat; //[fnNSigmaCompat]  0: n sigma for TPC compatibility band, 1: for TOF
238   Double_t fMaxnSigmaCombined[3];  // nSigma cut for pi,K,p (TPC^2+TOF^2)
239   Double_t fMinnSigmaTPC[3]; //min. of nSigma range for pi,K,p in TPC (match==5)
240   Double_t fMaxnSigmaTPC[3]; //max. of nSigma range for pi,K,p in TPC (match==5)
241   Double_t fMinnSigmaTOF[3]; //min. of nSigma range for pi,K,p in TOF (match==5)
242   Double_t fMaxnSigmaTOF[3]; //max. of nSigma range for pi,K,p in TOF (match==5)
243   Bool_t fMC; // MC(kTRUE) or real data (kFALSE, default option)
244   Bool_t fOnePad; //  real data with one pad clusters
245   Bool_t fMCLowEn2011; //  MC for low energy MC
246   Bool_t fppLowEn2011; //  Data for low energy pp 2011
247   Bool_t fPbPb; //  real data PbPb
248   Bool_t fTOFdecide; //  real data PbPb
249   Bool_t fOldPid; //  old PID method implemented
250   Double_t fPtThresholdTPC; //  pT threshold to use TPC PID
251   Double_t fMaxTrackMomForCombinedPID; // momentum threshold to use PID
252   AliPIDResponse *fPidResponse; //! pid response
253   AliPIDCombined* fPidCombined; //! combined PID object
254   
255   AliTPCPIDResponse* fTPCResponse; //! TPC response
256   
257   TH1F* fPriorsH[AliPID::kSPECIES]; // priors histos
258   ECombDetectors fCombDetectors; // detectors to be involved for combined PID
259   Bool_t fUseCombined; // detectors to be involved for combined PID
260   Bool_t fDefaultPriors; // use default priors for combined PID
261   
262   ClassDef(AliAODPidHF,23) // AliAODPid for heavy flavor PID
263   
264 };
265
266 #endif
267