8bf4d56f4e4c680ca778896e378f7a3be4f53e23
[u/mrichter/AliRoot.git] / STEER / ESD / AliESDpid.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------
19 //           Implementation of the combined PID class
20 //           For the Event Summary Data Class
21 //           produced by the reconstruction process
22 //           and containing information on the particle identification
23 //      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
24 //-----------------------------------------------------------------
25
26 #include "TArrayI.h"
27 #include "TArrayF.h"
28
29 #include "TRandom.h"
30 #include "AliLog.h"
31 #include "AliPID.h"
32 #include "AliTOFHeader.h"
33 #include "AliESDpid.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliMCEvent.h"
37 #include "AliMCEventHandler.h"
38 #include "AliAnalysisManager.h"
39
40
41 ClassImp(AliESDpid)
42
43 Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF) const {
44   //
45   //  Calculate probabilities for all detectors, except if TPConly==kTRUE
46   //  and combine PID
47   //  
48   //   Option TPConly==kTRUE is used during reconstruction, 
49   //  because ITS tracking uses TPC pid
50   //  HMPID and TRD pid are done in detector reconstructors
51   //
52
53   /*
54   Float_t timeZeroTOF = 0;
55   if (subtractT0) 
56     timeZeroTOF = event->GetT0();
57   */
58   Int_t nTrk=event->GetNumberOfTracks();
59   for (Int_t iTrk=0; iTrk<nTrk; iTrk++) {  
60     AliESDtrack *track=event->GetTrack(iTrk);
61     MakeTPCPID(track);
62     if (!TPConly) {
63       MakeITSPID(track);
64       MakeTOFPID(track, timeZeroTOF);
65       //MakeHMPIDPID(track);
66       //MakeTRDPID(track);
67     }
68     CombinePID(track);
69   }
70   return 0;
71 }
72 //_________________________________________________________________________
73 Float_t AliESDpid::GetTPCsignalTunedOnData(const AliVTrack *t) const {
74     AliESDtrack *track = (AliESDtrack *) t;
75     Float_t dedx = track->GetTPCsignalTunedOnData();
76
77     if(dedx > 0) return dedx;
78
79     Double_t mom = t->GetTPCmomentum();
80     dedx = t->GetTPCsignal();
81     track->SetTPCsignalTunedOnData(dedx);
82     if(dedx < 20) return dedx;
83
84     AliPID::EParticleType type = AliPID::kPion;
85
86     AliMCEventHandler* eventHandler=NULL;
87     eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
88     if (eventHandler) {
89         AliMCEvent* mcEvent = eventHandler->MCEvent();
90         if(mcEvent){
91             Bool_t kGood = kTRUE;
92             AliMCParticle *MCpart = (AliMCParticle *) mcEvent->GetTrack(TMath::Abs(t->GetLabel()));
93             TParticle *part = MCpart->Particle();
94             
95             Int_t iS = TMath::Abs(part->GetPdgCode());
96
97             if(iS==AliPID::ParticleCode(AliPID::kElectron)){
98                 type = AliPID::kElectron;
99             }
100             else if(iS==AliPID::ParticleCode(AliPID::kMuon)){
101                 type = AliPID::kMuon;
102             }
103             else if(iS==AliPID::ParticleCode(AliPID::kPion)){
104                 type = AliPID::kPion;
105             }
106             else if(iS==AliPID::ParticleCode(AliPID::kKaon)){
107                 type = AliPID::kKaon;
108             }
109             else if(iS==AliPID::ParticleCode(AliPID::kProton)){
110                 type = AliPID::kProton;
111             }
112             else if(iS==AliPID::ParticleCode(AliPID::kDeuteron)){ // d
113                 type = AliPID::kDeuteron;
114             }
115             else if(iS==AliPID::ParticleCode(AliPID::kTriton)){ // t
116                 type = AliPID::kTriton;
117             }
118             else if(iS==AliPID::ParticleCode(AliPID::kHe3)){ // 3He
119                 type = AliPID::kHe3;
120             }
121             else if(iS==AliPID::ParticleCode(AliPID::kAlpha)){ // 4He
122                 type = AliPID::kAlpha;
123             }
124             else
125                 kGood = kFALSE;
126
127             if(kGood){
128                 Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
129                 Double_t sigma=fTPCResponse.GetExpectedSigma(mom,t->GetTPCsignalN(),type);
130                 dedx = gRandom->Gaus(bethe,sigma);
131                 if(iS == AliPID::ParticleCode(AliPID::kHe3) || iS == AliPID::ParticleCode(AliPID::kAlpha)) dedx *= 5;
132             }
133         }
134     }
135
136     track->SetTPCsignalTunedOnData(dedx);
137     return dedx;
138 }
139 //_________________________________________________________________________
140 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
141 {
142   //
143   //  TPC pid using bethe-bloch and gaussian response
144   //
145   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
146     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
147
148     Double_t mom = track->GetP();
149     const AliExternalTrackParam *in=track->GetInnerParam();
150     if (in) mom = in->GetP();
151
152     Double_t p[AliPID::kSPECIES];
153     Double_t dedx=track->GetTPCsignal(); 
154     Bool_t mismatch=kTRUE, heavy=kTRUE;
155
156     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
157       AliPID::EParticleType type=AliPID::EParticleType(j);
158       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
159       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
160       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
161         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
162       } else {
163         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
164         mismatch=kFALSE;
165       }
166
167       // Check for particles heavier than (AliPID::kSPECIES - 1)
168       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
169
170     }
171
172     if (mismatch)
173        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
174
175     track->SetTPCpid(p);
176
177     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
178
179 }
180 //_________________________________________________________________________
181 void AliESDpid::MakeITSPID(AliESDtrack *track) const
182 {
183   //
184   // ITS PID
185   // Two options, depending on fITSPIDmethod:
186   //  1) Truncated mean method
187   //  2) Likelihood, using charges measured in all 4 layers and 
188   //     Landau+gaus response functions
189   //
190
191   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
192       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
193
194   Double_t mom=track->GetP();  
195   if (fITSPIDmethod == kITSTruncMean) {
196     Double_t dedx=track->GetITSsignal();
197     Bool_t isSA=kTRUE;
198     Double_t momITS=mom;
199     ULong_t trStatus=track->GetStatus();
200     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
201     UChar_t clumap=track->GetITSClusterMap();
202     Int_t nPointsForPid=0;
203     for(Int_t i=2; i<6; i++){
204       if(clumap&(1<<i)) ++nPointsForPid;
205     }
206
207     if(nPointsForPid<3) { // track not to be used for combined PID purposes
208       track->ResetStatus(AliESDtrack::kITSpid);
209       return;
210     }
211
212     Double_t p[10];
213
214     Bool_t mismatch=kTRUE, heavy=kTRUE;
215     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
216       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
217       Double_t bethe=fITSResponse.Bethe(momITS,mass);
218       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
219       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
220         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
221       } else {
222         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
223         mismatch=kFALSE;
224       }
225
226       // Check for particles heavier than (AliPID::kSPECIES - 1)
227       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
228
229     }
230
231     if (mismatch)
232        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
233
234     track->SetITSpid(p);
235
236     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
237   }
238   else {  // Likelihood method
239     Double_t condprobfun[AliPID::kSPECIES];
240     Double_t qclu[4];
241     track->GetITSdEdxSamples(qclu);
242     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
243     track->SetITSpid(condprobfun);
244   }
245
246 }
247 //_________________________________________________________________________
248 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
249 {
250   //
251   //   TOF PID using gaussian response
252   //
253
254   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
255   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
256   if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
257
258   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
259   Float_t timezero = fTOFResponse.GetT0bin(ibin);
260
261   Double_t time[AliPID::kSPECIESN];
262   track->GetIntegratedTimes(time);
263
264   Double_t sigma[AliPID::kSPECIES];
265   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
266     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
267   }
268
269   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
270            Form("Expected TOF signals [ps]: %f %f %f %f %f",
271                   time[AliPID::kElectron],
272                   time[AliPID::kMuon],
273                   time[AliPID::kPion],
274                   time[AliPID::kKaon],
275                   time[AliPID::kProton]));
276
277   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
278            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
279                   sigma[AliPID::kElectron],
280                   sigma[AliPID::kMuon],
281                   sigma[AliPID::kPion],
282                   sigma[AliPID::kKaon],
283                   sigma[AliPID::kProton]
284                   ));
285
286   Double_t tof = track->GetTOFsignal() - timezero;
287
288   Double_t p[AliPID::kSPECIES];
289   Bool_t mismatch = kTRUE, heavy = kTRUE;
290   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
291     Double_t sig = sigma[j];
292     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
293         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
294     } else
295       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
296
297     // Check the mismatching
298     
299     Double_t mass = AliPID::ParticleMass(j);
300     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
301     if (p[j]>pm) mismatch = kFALSE;
302
303     // Check for particles heavier than (AliPID::kSPECIES - 1)
304     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
305
306   }
307
308   /*
309     if (mismatch)
310     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
311   */
312
313   track->SetTOFpid(p);
314
315   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
316   
317   if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
318   else track->ResetStatus(AliESDtrack::kTOFmismatch);
319 }
320 //_________________________________________________________________________
321 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
322 {
323   //
324   // Method to recalculate the TRD PID probabilities
325   //
326   Double_t prob[AliPID::kSPECIES];
327   ComputeTRDProbability(track, AliPID::kSPECIES, prob);
328   track->SetTRDpid(prob);
329 }
330 //_________________________________________________________________________
331 void AliESDpid::CombinePID(AliESDtrack *track) const
332 {
333   //
334   // Combine the information of various detectors
335   // to determine the Particle Identification
336   //
337   Int_t ns=AliPID::kSPECIES;
338   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
339
340   if (track->IsOn(AliESDtrack::kITSpid)) {
341     Double_t d[10];
342     track->GetITSpid(d);
343     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
344   }
345
346   if (track->IsOn(AliESDtrack::kTPCpid)) {
347     Double_t d[10];
348     track->GetTPCpid(d);
349     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
350   }
351
352   if (track->IsOn(AliESDtrack::kTRDpid)) {
353     Double_t d[10];
354     track->GetTRDpid(d);
355     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
356   }
357
358   if (track->IsOn(AliESDtrack::kTOFpid)) {
359     Double_t d[10];
360     track->GetTOFpid(d);
361     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
362   }
363
364   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
365     Double_t d[10];
366     track->GetHMPIDpid(d);
367     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
368   }
369
370   track->SetESDpid(p);
371 }
372 //_________________________________________________________________________
373 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
374   //
375   // Check pid matching of TOF with TPC as reference
376   //
377     Bool_t status = kFALSE;
378     
379     Double_t exptimes[5];
380     track->GetIntegratedTimes(exptimes);
381     
382     Float_t p = track->P();
383     
384     Float_t dedx = track->GetTPCsignal();
385     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
386     
387     Double_t ptpc[3];
388     track->GetInnerPxPyPz(ptpc);
389     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
390     
391     for(Int_t i=0;i < 5;i++){
392         AliPID::EParticleType type=AliPID::EParticleType(i);
393         
394         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
395         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
396             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
397             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
398             
399             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
400                 status = kTRUE;
401             }
402         }
403     }
404     
405     // for nuclei
406     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
407     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
408     
409     
410     return status;
411 }