]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliFlowBayesianPID.h
From Francesco: Added AliFlowVZEROResults, PID-VZERO class restructuring
[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 AliTOFGeometry;
23 class AliTOFT0maker;
24
25 /*
26 HOW TO
27
28 1) in your main program:
29
30   AliFlowBayesianPID *mypid = new AliFlowBayesianPID();
31
32 or (if you want a global AliESDpid object)
33
34   AliESDpid *PIDesd = new AliESDpid();
35   ....
36   AliFlowBayesianPID *mypid = new AliFlowBayesianPID(PIDesd);
37
38 for data starting from the PbPb pass2 reconstruction
39
40   mypid->SetNewTrackParam();
41
42
43 before to loop on the on ESD tracks 
44  mypid->SetDetResponse(esdEvent, centrality,AliESDpid::kTOF_T0,kFALSE); // centrality > 0 = PbPb or -1 for pp collisions
45
46 before to loop on the on AOD tracks 
47  mypid->SetDetResponse(aodEvent, centrality); // centrality > 0 = PbPb or -1 for pp collisions
48
49 if you want to use a dE/dx depdence on DeltaPhi set the EP with its resolution otherwise it will be skipped
50  mypid->SetPsiCorrectionDeDx(psiEP,resEP);
51
52
53      for(...){ // track loop
54        mypid->ComputeProb(track,centrality); // both for ESD and AOD tracks
55        Float_t *prob = mypid->GetProb(); // Bayesian Probability (from 0 to 4) (Combined TPC || TOF) including a tuning of priors and TOF mismatch parameterization
56      }
57
58
59 More details:
60      // for the single detector weights (no priors)
61      Float_t *tpcWeight = mypid->GetWeights(0); // TPC weights (equal weights in case of no kTPCpid)
62      Float_t *tofWeight = mypid->GetWeights(1); // TOF weights (equal weights in case of no kTOFpid)
63
64      // more tools
65      TSpline3 *mismatchShape = mypid->GetMismatch(); // time distribution of mismatched tracks - L/c (L is the distance of the pad from the PV)
66      TF1 *tpcprob = mypid->GetTPCprob(); // normalized TPC response function (in Nsigmas)  
67      TF1 *tofprob = mypid->GetTOFprob(); // normalized TOF response function (in Nsigmas)  
68
69      TH2D *hPr = mypid->GetHistoPriors(isp); // 2D (centrality - pT) histo for the priors of specie-isp (centrality < 0 means pp collisions)
70                                              // all the priors are normalized to the pion ones
71
72 */
73
74 class AliFlowBayesianPID : public AliPIDResponse{
75  public:
76   AliFlowBayesianPID(AliESDpid *esdpid=NULL); 
77   virtual ~AliFlowBayesianPID();
78   
79   // virtual method of AliPIDResponse
80   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
81
82   // setter
83   void SetDetResponse(AliESDEvent *esd,Float_t centrality=-1.0,EStartTimeType_t flagStart=AliESDpid::kTOF_T0,Bool_t recomputeT0TOF=kFALSE);
84   void SetDetResponse(AliAODEvent *aod,Float_t centrality=-1.0);
85   void SetNewTrackParam(Bool_t flag=kTRUE){fNewTrackParam=flag;};
86   void SetDetAND(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskAND[idet] = kTRUE;};
87   void SetDetOR(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskOR[idet] = kTRUE;};
88   void ResetDetAND(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskAND[idet] = kFALSE;};
89   void ResetDetOR(Int_t idet){if(idet < fgkNdetectors && idet >= 0) fMaskOR[idet] = kFALSE;};
90   void SetPsiCorrectionDeDx(Float_t psi,Float_t res);
91   void SetMC(Bool_t flag){fIsMC=flag;};
92
93   // getter
94   AliESDpid* GetESDpid(){return fPIDesd;};
95   const TH2D *GetHistoPriors(Int_t specie) const {if(specie >=0 && specie < fgkNspecies) return fghPriors[specie]; else return NULL;};
96   TSpline3 *GetMismatch();  
97   const TF1 *GetTOFprob() const {return fTOFResponseF;};
98   const TF1 *GetTPCprob() const {return fTPCResponseF;};
99   const Float_t *GetWeights(Int_t det) const {if(det < fgkNdetectors && det >= 0){return fWeights[det];} else{return NULL;}};
100   Float_t *GetProb() {return fProb;};
101   Float_t GetTOFMismWeight() const {return fWTofMism;};
102   Float_t GetTOFMismProb() const {return fProbTofMism;};
103   Float_t GetMassOverZ() const {return fMassTOF;};
104   Float_t GetZ() const {return fZ;};
105   Bool_t GetDetANDstatus(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskAND[idet];} else{return kFALSE;} };
106   Bool_t GetDetORstatus(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskOR[idet];} else{return kFALSE;} };
107   Bool_t GetCurrentMask(Int_t idet) const {if(idet < fgkNdetectors && idet >= 0){return fMaskCurrent[idet];} else{return kFALSE;} };
108   Float_t GetExpDeDx(const AliESDtrack *t,Int_t iS) const;
109   Float_t GetExpDeDx(const AliAODTrack *t,Int_t iS) 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 = 8;// 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   static AliTOFGeometry *fgTofGeo; // TOF geometry needed to reproduce mismatch shape
130
131   AliESDpid *fPIDesd;//ESDpid object
132   TDatabasePDG *fDB; // Database pdg
133   Double_t fMass[fgkNspecies]; // mass for el(0),mu(1),pi(2),K(3),p(4)
134
135   Bool_t fNewTrackParam; // switch for new tracking resolution TOF parameterization
136   Float_t fTOFresolution; // TOF res needed only if T0-TOF should be recomputed
137
138   AliFlowBayesianPID(const AliFlowBayesianPID&); // not implemented
139   AliFlowBayesianPID& operator=(const AliFlowBayesianPID&); // not implemented 
140
141   TF1 *fTOFResponseF; // TOF Gaussian+tail response function (tail at 1.1 sigma)
142   TF1 *fTPCResponseF; // TPC Gaussian+tail response function (tail at 1.8 sigma)
143
144   AliTOFT0maker *fTOFmaker; //TOF-T0 maker object
145
146   Float_t fWeights[fgkNdetectors][fgkNspecies]; // weights: 0=tpc,1=tof
147   Float_t fProb[fgkNspecies],fWTofMism,fProbTofMism; // Bayesian Combined PID + mismatch weights and probability 
148
149   Float_t fZ,fMassTOF; //measured charge(Z) and mass/Z
150   TF1 *fBBdata; // Bethe Bloch function (needed to compute the charge of the particle)
151
152   Bool_t fMaskAND[fgkNdetectors],fMaskOR[fgkNdetectors],fMaskCurrent[fgkNdetectors]; // mask detector should be used
153
154   Float_t fCurrCentrality; // Centrality in current event
155
156   Float_t fPsi,fPsiRes;    // RP and its resolution for the event (999 if not available) to correct dEdx for p > 1
157
158   Bool_t fIsMC; // switch for MC analysis
159
160   Float_t fDedx; // dE/dx tuned for MC
161
162   ClassDef(AliFlowBayesianPID, 6); // example of analysis
163 };
164
165 #endif
166
167
168