]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliFlowBayesianPID.h
From Francesco: Updated Bayesian PID code
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliFlowBayesianPID.h
1 /*
2 -----------------------------------------------------------------------------------
3  AliFlowBayesianPID class implemented by F. Noferini (noferini@bo.infn.it)
4  needed to use Bayesian probability in the flow package (used in AliFlowTrackCuts)
5 -----------------------------------------------------------------------------------
6 */
7 #ifndef ALIFLOWBAYESIANPID_H
8 #define ALIFLOWBAYESIANPID_H
9
10 #include "AliESDpid.h"
11 #include "AliPIDResponse.h"
12
13 class TObject;
14 class TDatabasePDG;
15 class AliESDEvent;
16 class AliESDtrack;
17 class AliAODEvent;
18 class AliAODTrack;
19 class TH2D;
20 class TSpline3;
21 class TF1;
22 class TH1D;
23
24 /*
25 HOW TO
26
27 1) in your main program:
28
29   AliFlowBayesianPID *mypid = new AliFlowBayesianPID();
30
31 or (if you want a global AliESDpid object)
32
33   AliESDpid *PIDesd = new AliESDpid();
34   ....
35   AliFlowBayesianPID *mypid = new AliFlowBayesianPID(PIDesd);
36
37 for data starting from the PbPb pass2 reconstruction
38
39   mypid->SetNewTrackParam();
40
41
42 before to loop on the on ESD tracks 
43  mypid->SetDetResponse(esdEvent, centrality,AliESDpid::kTOF_T0,kFALSE); // centrality > 0 = PbPb or -1 for pp collisions
44
45 before to loop on the on AOD tracks 
46  mypid->SetDetResponse(aodEvent, centrality); // centrality > 0 = PbPb or -1 for pp collisions
47
48 if you want to use a dE/dx depdence on DeltaPhi set the EP with its resolution otherwise it will be skipped
49  mypid->SetPsiCorrectionDeDx(psiEP,resEP);
50
51
52      for(...){ // track loop
53        mypid->ComputeProb(track,centrality); // both for ESD and AOD tracks
54        Float_t *prob = mypid->GetProb(); // Bayesian Probability (from 0 to 4) (Combined TPC || TOF) including a tuning of priors and TOF mismatch parameterization
55      }
56
57
58 More details:
59      // for the single detector weights (no priors)
60      Float_t *tpcWeight = mypid->GetWeights(0); // TPC weights (equal weights in case of no kTPCpid)
61      Float_t *tofWeight = mypid->GetWeights(1); // TOF weights (equal weights in case of no kTOFpid)
62
63      // more tools
64      TSpline3 *mismatchShape = mypid->GetMismatch(); // time distribution of mismatched tracks - L/c (L is the distance of the pad from the PV)
65      TF1 *tpcprob = mypid->GetTPCprob(); // normalized TPC response function (in Nsigmas)  
66      TF1 *tofprob = mypid->GetTOFprob(); // normalized TOF response function (in Nsigmas)  
67
68      TH2D *hPr = mypid->GetHistoPriors(isp); // 2D (centrality - pT) histo for the priors of specie-isp (centrality < 0 means pp collisions)
69                                              // all the priors are normalized to the pion ones
70
71 */
72
73 class AliFlowBayesianPID : public AliPIDResponse{
74  public:
75   AliFlowBayesianPID(AliESDpid *esdpid=NULL); 
76   virtual ~AliFlowBayesianPID();
77   
78   // virtual method of AliPIDResponse
79   virtual Float_t NumberOfSigmasTOF(const AliVParticle *vtrack, AliPID::EParticleType type) const {if(vtrack) printf("Don't call AliFlowBayesianPID::NumberOfSigmasTOF method (%i)\n",type); return 0.0;} // do not use it
80
81   // setter
82   void SetDetResponse(AliESDEvent *esd,Float_t centrality=-1.0,EStartTimeType_t flagStart=AliESDpid::kTOF_T0,Bool_t /*recomputeT0TOF*/=kFALSE);
83   void SetDetResponse(AliAODEvent *aod,Float_t centrality=-1.0);
84   void SetNewTrackParam(Bool_t flag=kTRUE){fNewTrackParam=flag;};
85   void SetDetAND(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskAND[idet] = kTRUE;};
86   void SetDetOR(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskOR[idet] = kTRUE;};
87   void ResetDetAND(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskAND[idet] = kFALSE;};
88   void ResetDetOR(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskOR[idet] = kFALSE;};
89   void SetPsiCorrectionDeDx(Float_t psi,Float_t res);
90   void SetMC(Bool_t flag){fIsMC=flag;};
91
92   // getter
93   AliESDpid* GetESDpid(){return fPIDesd;};
94   const TH2D *GetHistoPriors(Int_t specie) const {if(specie >=0 && specie < fgkNspecies) return fghPriors[specie]; else return NULL;};
95   TSpline3 *GetMismatch();  
96   const TF1 *GetTOFprob() const {return fTOFResponseF;};
97   const TF1 *GetTPCprob() const {return fTPCResponseF;};
98   const Float_t *GetWeights(Int_t det) const {if(det < fgkNdetectors && det >= 0){return fWeights[det];} else{return NULL;}};
99   Float_t *GetProb() {return fProb;};
100   Float_t GetTOFMismWeight() const {return fWTofMism;};
101   Float_t GetTOFMismProb() const {return fProbTofMism;};
102   Float_t GetMassOverZ() const {return fMassTOF;};
103   Float_t GetZ() const {return fZ;};
104   Bool_t GetDetANDstatus(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskAND[idet];} else{return kFALSE;} };
105   Bool_t GetDetORstatus(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskOR[idet];} else{return kFALSE;} };
106   Bool_t GetCurrentMask(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskCurrent[idet];} else{return kFALSE;} };
107
108   Float_t GetExpDeDx(const AliVTrack *t,Int_t iS) const;
109   Float_t GetExpDeDx(const AliVTrack *t,Float_t m) const;
110
111   // methods for Bayesina Combined PID
112   void ComputeWeights(const AliESDtrack *t);
113   void ComputeProb(const AliESDtrack *t,Float_t); // obsolete method
114   void ComputeProb(const AliESDtrack *t){ComputeProb(t,0.0);}; 
115   void ComputeWeights(const AliAODTrack *t,AliAODEvent *aod=NULL);
116   void ComputeProb(const AliAODTrack *t,AliAODEvent *aod=NULL); // obsolete method
117
118   void SetTOFres(Float_t res){fTOFresolution=res;};
119
120   Float_t GetDeDx() const {return fDedx;};
121
122  private: 
123   void SetPriors();
124
125   static const Int_t fgkNdetectors = 2; // Number of detector used for PID
126   static const Int_t fgkNspecies = 9;// 0=el, 1=mu, 2=pi, 3=ka, 4=pr, 5=deuteron, 6=triton, 7=He3 
127   static TH2D* fghPriors[fgkNspecies]; // histo with priors (hardcoded)
128   static TSpline3 *fgMism; // function for mismatch
129
130   AliESDpid *fPIDesd;//ESDpid object
131   TDatabasePDG *fDB; // Database pdg
132   Double_t fMass[fgkNspecies]; // mass for el(0),mu(1),pi(2),K(3),p(4)
133
134   Bool_t fNewTrackParam; // switch for new tracking resolution TOF parameterization
135   Float_t fTOFresolution; // TOF res needed only if T0-TOF should be recomputed
136
137   AliFlowBayesianPID(const AliFlowBayesianPID&); // not implemented
138   AliFlowBayesianPID& operator=(const AliFlowBayesianPID&); // not implemented 
139
140   TF1 *fTOFResponseF; // TOF Gaussian+tail response function (tail at 1.1 sigma)
141   TF1 *fTPCResponseF; // TPC Gaussian+tail response function (tail at 1.8 sigma)
142
143   Float_t fWeights[fgkNdetectors][fgkNspecies]; // weights: 0=tpc,1=tof
144   Float_t fProb[fgkNspecies],fWTofMism,fProbTofMism; // Bayesian Combined PID + mismatch weights and probability 
145
146   Float_t fZ,fMassTOF; //measured charge(Z) and mass/Z
147   TF1 *fBBdata; // Bethe Bloch function (needed to compute the charge of the particle)
148
149   Bool_t fMaskAND[fgkNdetectors],fMaskOR[fgkNdetectors],fMaskCurrent[fgkNdetectors]; // mask detector should be used
150
151   Float_t fCurrCentrality; // Centrality in current event
152
153   Float_t fPsi,fPsiRes;    // RP and its resolution for the event (999 if not available) to correct dEdx for p > 1
154
155   Bool_t fIsMC; // switch for MC analysis
156
157   Float_t fDedx; // dE/dx tuned for MC
158
159   static TH1D *fgHtofChannelDist; // channel distance from IP
160
161   ClassDef(AliFlowBayesianPID, 8); // example of analysis
162 };
163
164 #endif
165
166
167