]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAODPidHF.h
Merge branch 'TPCdev' of https://git.cern.ch/reps/AliRoot into TPCdev
[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 SetTPC(Bool_t tpc){fTPC=tpc;return;}
92   void SetTOF(Bool_t tof){fTOF=tof;return;}
93   void SetITS(Bool_t its){fITS=its;return;}
94   void SetTRD(Bool_t trd){fTRD=trd;return;}
95   void SetMatch(Int_t match){fMatch=match;return;}
96   void SetCompat(Bool_t comp){fCompat=comp;return;}
97   void SetMC(Bool_t mc){fMC=mc;return;}
98   void SetMClowenpp2011(Bool_t mc){fMCLowEn2011=mc;return;}
99   void SetOnePad(Bool_t onepad){fOnePad=onepad;return;}
100   void SetppLowEn2011(Bool_t opt){fppLowEn2011=opt;return;}
101   void SetPbPb(Bool_t pbpb){fPbPb=pbpb;return;}
102   void SetPCompatTOF(Double_t pTOF){fPCompatTOF=pTOF;return;}
103   void SetTOFdecide(Bool_t tOFdecide){fTOFdecide=tOFdecide;return;}
104   void SetOldPid(Bool_t oldPid){fOldPid=oldPid;return;}
105   void SetPtThresholdTPC(Double_t ptThresholdTPC){fPtThresholdTPC=ptThresholdTPC;return;}
106   void SetMaxTrackMomForCombinedPID(Double_t mom){fMaxTrackMomForCombinedPID=mom;}
107   void SetPidResponse(AliPIDResponse *pidResp) {fPidResponse=pidResp;return;}
108   void SetCombDetectors(ECombDetectors pidComb) {
109     fCombDetectors=pidComb;
110   }
111   void SetPionPriorHisto(TH1F* histo){
112     if(fPriorsH[AliPID::kPion]) delete fPriorsH[AliPID::kPion];
113     fPriorsH[AliPID::kPion] = new TH1F(*histo);
114   }
115   void SetKaonPriorHisto(TH1F* histo){
116     if(fPriorsH[AliPID::kKaon]) delete fPriorsH[AliPID::kKaon];
117     fPriorsH[AliPID::kKaon] = new TH1F(*histo);
118   }
119   void SetProtonPriorHisto(TH1F* histo){
120     if(fPriorsH[AliPID::kProton]) delete fPriorsH[AliPID::kProton];
121     fPriorsH[AliPID::kProton] = new TH1F(*histo);
122   }
123
124  
125   //Getters
126  
127   Int_t GetnSigmaTPC(AliAODTrack *track, Int_t species, Double_t &sigma) const;
128   Int_t GetnSigmaTOF(AliAODTrack *track, Int_t species, Double_t &sigma) const;
129   Int_t GetnSigmaITS(AliAODTrack *track, Int_t species, Double_t &sigma) const;
130   Double_t GetSigma(Int_t idet) const{return fnSigma[idet];}
131   Double_t GetTofSigma() const{return fTOFSigma;}
132   //void GetPriors(Double_t *priors) const{priors=fPriors;return;}
133   //void GetPLimit(Double_t *plim) const{plim=fPLimit;}
134   void GetPriors(Double_t *priors) const{for(Int_t i=0;i<fnPriors;i++){priors[i]=fPriors[i];}return;}
135   void GetPLimit(Double_t *plim) const{for(Int_t i=0;i<fnPLimit;i++){plim[i]=fPLimit[i];}return;}
136   Bool_t GetAsym() const{return fAsym;}
137   Bool_t GetTPC() const{return fTPC;}
138   Bool_t GetTOF() const{return fTOF;}
139   Bool_t GetITS() const{return fITS;}
140   Bool_t GetTRD() const{return fTRD;}
141   Int_t GetMatch() const{return fMatch;}
142   Bool_t GetCompat() const{return fCompat;}
143   Bool_t GetMC() const{return fMC;}
144   Bool_t GetOnePad() const{return fOnePad;}
145   Bool_t GetppLowEn2011() const {return fppLowEn2011;}
146   Bool_t GetMCLowEn2011() const {return fMCLowEn2011;}
147   Bool_t GetPbPb() const{return fPbPb;}
148   Bool_t GetTOFdecide() const{return fTOFdecide;}
149   Double_t GetPCompatTOF() const{return fPCompatTOF;}
150   Double_t GetnSigmaCompatTPC() const{return fnSigmaCompat[0];}
151   Double_t GetnSigmaCompatTOF() const{return fnSigmaCompat[1];}
152   Bool_t GetOldPid(){return fOldPid;}
153   Double_t GetPtThresholdTPC(){return fPtThresholdTPC;}
154   Double_t GetMaxTrackMomForCombinedPID(){return fMaxTrackMomForCombinedPID;}
155   AliPIDResponse *GetPidResponse() const {return fPidResponse;}
156   AliPIDCombined *GetPidCombined() const {return fPidCombined;}
157   ECombDetectors GetCombDetectors() const {
158     return fCombDetectors;
159   }
160   Bool_t GetUseCombined() {return fUseCombined;}
161   Bool_t GetDefaultPriors() {return fDefaultPriors;}
162  
163   Int_t RawSignalPID (AliAODTrack *track, TString detector) const;
164   Bool_t IsKaonRaw (AliAODTrack *track, TString detector) const;
165   Bool_t IsPionRaw (AliAODTrack *track, TString detector) const;
166   Bool_t IsProtonRaw (AliAODTrack *track, TString detector) const;
167   Bool_t IsElectronRaw (AliAODTrack *track, TString detector) const;
168   void CombinedProbability(AliAODTrack *track,Bool_t *type) const; //0 = pion, 1 = kaon, 2 = proton
169   Bool_t CheckStatus(AliAODTrack *track,TString detectors) const;
170
171   Bool_t CheckITSPIDStatus(AliAODTrack *track) const;
172   Bool_t CheckTPCPIDStatus(AliAODTrack *track) const;
173   Bool_t CheckTOFPIDStatus(AliAODTrack *track) const;
174   Bool_t CheckTRDPIDStatus(AliAODTrack *track) const;
175
176   Bool_t TPCRawAsym(AliAODTrack* track,Int_t specie) const;
177   Int_t MatchTPCTOF(AliAODTrack *track,Int_t specie);
178
179   Int_t MakeRawPid(AliAODTrack *track,Int_t specie); //general method to perform PID using raw signals
180
181   Bool_t IsTOFPiKexcluded(AliAODTrack *track,Double_t nsigmaK);
182
183   Bool_t IsExcluded(AliAODTrack *track, Int_t labelTrack, Double_t nsigmaCut, TString detector);
184
185   void GetTPCBetheBlochParams(Double_t alephParameters[5]) const;
186   void SetBetheBloch();
187   // method for AliPIDCombined object
188   void SetSelectedSpecies(Int_t ispecies = AliPID::kSPECIES){GetPidCombined()->SetSelectedSpecies(ispecies);};
189   void SetPriorDistribution(AliPID::EParticleType type,TH1F *prior);
190   void DrawPrior(AliPID::EParticleType type);
191   void SetPriorsHistos(TString priorFileName);
192   void SetUpCombinedPID();
193   void SetUseCombined(Bool_t useCombined=kTRUE) {fUseCombined=useCombined;}
194   void SetUseDefaultPriors(Bool_t defaultP)         {fDefaultPriors=defaultP;}
195   Int_t ApplyPidTPCRaw(AliAODTrack *track,Int_t specie) const;
196   Int_t ApplyPidTOFRaw(AliAODTrack *track,Int_t specie) const;
197   Int_t ApplyPidITSRaw(AliAODTrack *track,Int_t specie) const;
198
199   void PrintAll() const;
200
201  protected:
202
203
204  private:
205   Int_t fnNSigma; // number of sigmas
206   Double_t *fnSigma; // [fnNSigma], sigma for the raw signal PID: 0-2 for TPC, 3 for TOF, 4 for ITS
207   Double_t fTOFSigma; // TOF precision 
208   Double_t fCutTOFmismatch; // Cut of TOF mismatch probability
209   UInt_t fMinNClustersTPCPID;       // Minimum TPC PID clusters cut
210   Int_t fnPriors; //number of priors
211   Double_t *fPriors; // [fnPriors], set of priors
212   Int_t fnPLimit; //number of Plimit
213   Double_t *fPLimit; // [fnPLimit], limit of p intervals for asimmetric PID: fPLimit<p[0], fPLimit[0]<p<fPLimit[1], p>fPLimit[1]
214   Bool_t fAsym; // asimmetric PID required
215   Bool_t fTPC; //switch to include or exclude TPC 
216   Bool_t fTOF; // switch to include or exclude TOF
217   Bool_t fITS; //switch to include or exclude ITS
218   Bool_t fTRD; // switch to include or exclude TRD
219   Int_t fMatch; //switch to combine the info from more detectors: 1 = || , 2 = &, 3 = p region
220   Bool_t fCompat; // compatibility region : useful only if fMatch=1
221   Double_t fPCompatTOF; //  compatibility p limit for TOF  
222   Int_t fnNSigmaCompat; // number of sigmas
223   Double_t *fnSigmaCompat; //[fnNSigmaCompat]  0: n sigma for TPC compatibility band, 1: for TOF  
224   Double_t fMaxnSigmaCombined[3];  // nSigma cut for pi,K,p (TPC^2+TOF^2)
225   Double_t fMinnSigmaTPC[3]; //min. of nSigma range for pi,K,p in TPC (match==5)
226   Double_t fMaxnSigmaTPC[3]; //max. of nSigma range for pi,K,p in TPC (match==5)
227   Double_t fMinnSigmaTOF[3]; //min. of nSigma range for pi,K,p in TOF (match==5)
228   Double_t fMaxnSigmaTOF[3]; //max. of nSigma range for pi,K,p in TOF (match==5)
229   Bool_t fMC; // MC(kTRUE) or real data (kFALSE, default option)
230   Bool_t fOnePad; //  real data with one pad clusters 
231   Bool_t fMCLowEn2011; //  MC for low energy MC
232   Bool_t fppLowEn2011; //  Data for low energy pp 2011
233   Bool_t fPbPb; //  real data PbPb 
234   Bool_t fTOFdecide; //  real data PbPb 
235   Bool_t fOldPid; //  old PID method implemented
236   Double_t fPtThresholdTPC; //  pT threshold to use TPC PID
237   Double_t fMaxTrackMomForCombinedPID; // momentum threshold to use PID
238   AliPIDResponse *fPidResponse; //! pid response
239   AliPIDCombined* fPidCombined; //! combined PID object 
240
241   AliTPCPIDResponse* fTPCResponse; //! TPC response 
242
243   TH1F* fPriorsH[AliPID::kSPECIES]; // priors histos
244   ECombDetectors fCombDetectors; // detectors to be involved for combined PID
245   Bool_t fUseCombined; // detectors to be involved for combined PID
246   Bool_t fDefaultPriors; // use default priors for combined PID
247
248   ClassDef(AliAODPidHF,21) // AliAODPid for heavy flavor PID
249
250     };
251
252 #endif
253