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