]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliESDpid.cxx
test
[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
38 #include <AliDetectorPID.h>
39
40 ClassImp(AliESDpid)
41
42 Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF) const {
43   //
44   //  Calculate probabilities for all detectors, except if TPConly==kTRUE
45   //  and combine PID
46   //  
47   //   Option TPConly==kTRUE is used during reconstruction, 
48   //  because ITS tracking uses TPC pid
49   //  HMPID and TRD pid are done in detector reconstructors
50   //
51
52   /*
53   Float_t timeZeroTOF = 0;
54   if (subtractT0) 
55     timeZeroTOF = event->GetT0();
56   */
57   Int_t nTrk=event->GetNumberOfTracks();
58   for (Int_t iTrk=0; iTrk<nTrk; iTrk++) {  
59     AliESDtrack *track=event->GetTrack(iTrk);
60     MakeTPCPID(track);
61     if (!TPConly) {
62       MakeITSPID(track);
63       MakeTOFPID(track, timeZeroTOF);
64       //MakeHMPIDPID(track);
65       //MakeTRDPID(track);
66     }
67     CombinePID(track);
68   }
69   return 0;
70 }
71 //_________________________________________________________________________
72 Float_t AliESDpid::GetTPCsignalTunedOnData(const AliVTrack *t) const {
73     AliESDtrack *track = (AliESDtrack *) t;
74     Float_t dedx = track->GetTPCsignalTunedOnData();
75
76     if(dedx > 0) return dedx;
77
78     dedx = t->GetTPCsignal();
79     track->SetTPCsignalTunedOnData(dedx);
80     if(dedx < 20) return dedx;
81
82     AliPID::EParticleType type = AliPID::kPion;
83
84     AliMCEventHandler* eventHandler=dynamic_cast<AliMCEventHandler*>(fEventHandler);
85     if (eventHandler) {
86         AliMCEvent* mcEvent = eventHandler->MCEvent();
87         if(mcEvent){
88             Bool_t kGood = kTRUE;
89             AliMCParticle *MCpart = (AliMCParticle *) mcEvent->GetTrack(TMath::Abs(t->GetLabel()));
90             TParticle *part = MCpart->Particle();
91             
92             Int_t iS = TMath::Abs(part->GetPdgCode());
93
94             if(iS==AliPID::ParticleCode(AliPID::kElectron)){
95                 type = AliPID::kElectron;
96             }
97             else if(iS==AliPID::ParticleCode(AliPID::kMuon)){
98                 type = AliPID::kMuon;
99             }
100             else if(iS==AliPID::ParticleCode(AliPID::kPion)){
101                 type = AliPID::kPion;
102             }
103             else if(iS==AliPID::ParticleCode(AliPID::kKaon)){
104                 type = AliPID::kKaon;
105             }
106             else if(iS==AliPID::ParticleCode(AliPID::kProton)){
107                 type = AliPID::kProton;
108             }
109             else if(iS==AliPID::ParticleCode(AliPID::kDeuteron)){ // d
110                 type = AliPID::kDeuteron;
111             }
112             else if(iS==AliPID::ParticleCode(AliPID::kTriton)){ // t
113                 type = AliPID::kTriton;
114             }
115             else if(iS==AliPID::ParticleCode(AliPID::kHe3)){ // 3He
116                 type = AliPID::kHe3;
117             }
118             else if(iS==AliPID::ParticleCode(AliPID::kAlpha)){ // 4He
119                 type = AliPID::kAlpha;
120             }
121             else
122                 kGood = kFALSE;
123
124             if(kGood){
125         //TODO maybe introduce different dEdxSources?
126         Double_t bethe = fTPCResponse.GetExpectedSignal(track, type, AliTPCPIDResponse::kdEdxDefault, this->UseTPCEtaCorrection(),
127                                                         this->UseTPCMultiplicityCorrection());
128         Double_t sigma = fTPCResponse.GetExpectedSigma(track, type, AliTPCPIDResponse::kdEdxDefault, this->UseTPCEtaCorrection(),
129                                                        this->UseTPCMultiplicityCorrection());
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 Float_t AliESDpid::GetTOFsignalTunedOnData(const AliVTrack *t) const {
141     AliESDtrack *track = (AliESDtrack *) t;
142     Double_t tofSignal = track->GetTOFsignalTunedOnData();
143
144     if(tofSignal <  99999) return (Float_t)tofSignal; // it has been already set
145
146     tofSignal = t->GetTOFsignal() + fTOFResponse.GetTailRandomValue();
147     track->SetTOFsignalTunedOnData(tofSignal);
148     return (Float_t)tofSignal;
149 }
150 //_________________________________________________________________________
151 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
152 {
153   //
154   //  TPC pid using bethe-bloch and gaussian response
155   //
156   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
157     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
158
159     Double_t mom = track->GetP();
160     const AliExternalTrackParam *in=track->GetInnerParam();
161     if (in) mom = in->GetP();
162
163     Double_t p[AliPID::kSPECIES];
164     Double_t dedx=track->GetTPCsignal(); 
165     Bool_t mismatch=kTRUE, heavy=kTRUE;
166
167     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
168       AliPID::EParticleType type=AliPID::EParticleType(j);
169       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
170       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
171       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
172         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
173       } else {
174         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
175         mismatch=kFALSE;
176       }
177
178       // Check for particles heavier than (AliPID::kSPECIES - 1)
179       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
180
181     }
182
183     if (mismatch)
184        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
185
186     track->SetTPCpid(p);
187
188     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
189
190 }
191 //_________________________________________________________________________
192 void AliESDpid::MakeITSPID(AliESDtrack *track) const
193 {
194   //
195   // ITS PID
196   // Two options, depending on fITSPIDmethod:
197   //  1) Truncated mean method
198   //  2) Likelihood, using charges measured in all 4 layers and 
199   //     Landau+gaus response functions
200   //
201
202   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
203       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
204
205   Double_t mom=track->GetP();  
206   if (fITSPIDmethod == kITSTruncMean) {
207     Double_t dedx=track->GetITSsignal();
208     Bool_t isSA=kTRUE;
209     Double_t momITS=mom;
210     ULong_t trStatus=track->GetStatus();
211     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
212     UChar_t clumap=track->GetITSClusterMap();
213     Int_t nPointsForPid=0;
214     for(Int_t i=2; i<6; i++){
215       if(clumap&(1<<i)) ++nPointsForPid;
216     }
217
218     if(nPointsForPid<3) { // track not to be used for combined PID purposes
219       track->ResetStatus(AliESDtrack::kITSpid);
220       return;
221     }
222
223     Double_t p[10];
224
225     Bool_t mismatch=kTRUE, heavy=kTRUE;
226     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
227       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
228       Double_t bethe=fITSResponse.Bethe(momITS,mass);
229       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
230       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
231         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
232       } else {
233         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
234         mismatch=kFALSE;
235       }
236
237       // Check for particles heavier than (AliPID::kSPECIES - 1)
238       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
239
240     }
241
242     if (mismatch)
243        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
244
245     track->SetITSpid(p);
246
247     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
248   }
249   else {  // Likelihood method
250     Double_t condprobfun[AliPID::kSPECIES];
251     Double_t qclu[4];
252     track->GetITSdEdxSamples(qclu);
253     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
254     track->SetITSpid(condprobfun);
255   }
256
257 }
258 //_________________________________________________________________________
259 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
260 {
261   //
262   //   TOF PID using gaussian response
263   //
264
265   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
266   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
267   if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
268
269   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
270   Float_t timezero = fTOFResponse.GetT0bin(ibin);
271
272   Double_t time[AliPID::kSPECIES];
273   track->GetIntegratedTimes(time);
274
275   Double_t sigma[AliPID::kSPECIES];
276   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
277     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
278   }
279
280   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
281            Form("Expected TOF signals [ps]: %f %f %f %f %f",
282                   time[AliPID::kElectron],
283                   time[AliPID::kMuon],
284                   time[AliPID::kPion],
285                   time[AliPID::kKaon],
286                   time[AliPID::kProton]));
287
288   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
289            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
290                   sigma[AliPID::kElectron],
291                   sigma[AliPID::kMuon],
292                   sigma[AliPID::kPion],
293                   sigma[AliPID::kKaon],
294                   sigma[AliPID::kProton]
295                   ));
296
297   Double_t tof = track->GetTOFsignal() - timezero;
298
299   Double_t p[AliPID::kSPECIES];
300 //   Bool_t mismatch = kTRUE;
301   Bool_t heavy = kTRUE;
302   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
303     Double_t sig = sigma[j];
304     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
305         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
306     } else
307       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
308
309     // Check the mismatching
310     
311 //     Double_t mass = AliPID::ParticleMass(j);
312 //     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
313 //     if (p[j]>pm) mismatch = kFALSE;
314
315     // Check for particles heavier than (AliPID::kSPECIES - 1)
316     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
317
318   }
319
320   /*
321     if (mismatch)
322     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
323   */
324
325   track->SetTOFpid(p);
326
327   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
328   
329   // kTOFmismatch flas is not set because deprecated from 18/02/2013
330   //  if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
331   //  else track->ResetStatus(AliESDtrack::kTOFmismatch);
332 }
333 //_________________________________________________________________________
334 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
335 {
336   //
337   // Method to recalculate the TRD PID probabilities
338   //
339   Double_t prob[AliPID::kSPECIES];
340   GetComputeTRDProbability(track, AliPID::kSPECIES, prob);
341   track->SetTRDpid(prob);
342 }
343 //_________________________________________________________________________
344 void AliESDpid::CombinePID(AliESDtrack *track) const
345 {
346   //
347   // Combine the information of various detectors
348   // to determine the Particle Identification
349   //
350   Int_t ns=AliPID::kSPECIES;
351   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
352
353   if (track->IsOn(AliESDtrack::kITSpid)) {
354     Double_t d[10];
355     track->GetITSpid(d);
356     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
357   }
358
359   if (track->IsOn(AliESDtrack::kTPCpid)) {
360     Double_t d[10];
361     track->GetTPCpid(d);
362     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
363   }
364
365   if (track->IsOn(AliESDtrack::kTRDpid)) {
366     Double_t d[10];
367     track->GetTRDpid(d);
368     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
369   }
370
371   if (track->IsOn(AliESDtrack::kTOFpid)) {
372     Double_t d[10];
373     track->GetTOFpid(d);
374     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
375   }
376
377   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
378     Double_t d[10];
379     track->GetHMPIDpid(d);
380     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
381   }
382
383   track->SetESDpid(p);
384 }
385 //_________________________________________________________________________
386 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
387   //
388   // Check pid matching of TOF with TPC as reference
389   //
390     Bool_t status = kFALSE;
391     
392     Double_t exptimes[5];
393     track->GetIntegratedTimes(exptimes);
394     
395     Float_t p = track->P();
396     
397     Float_t dedx = track->GetTPCsignal();
398     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
399     
400     Double_t ptpc[3];
401     track->GetInnerPxPyPz(ptpc);
402     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
403     
404     for(Int_t i=0;i < 5;i++){
405         AliPID::EParticleType type=AliPID::EParticleType(i);
406         
407         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
408         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
409             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
410             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
411             
412             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
413                 status = kTRUE;
414             }
415         }
416     }
417     
418     // for nuclei
419     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
420     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
421     
422     
423     return status;
424 }
425
426 //_________________________________________________________________________
427 Float_t AliESDpid::GetSignalDeltaTOFold(const AliVParticle *track, AliPID::EParticleType type, Bool_t ratio/*=kFALSE*/) const
428 {
429   //
430   // TOF signal - expected
431   //
432   AliVTrack *vtrack=(AliVTrack*)track;
433   
434   const Double_t expTime = fTOFResponse.GetExpectedSignal(vtrack,type);
435   Double_t tofTime = 99999;
436   if (fTuneMConData && ((fTuneMConDataMask & kDetTOF) == kDetTOF) ) tofTime = (Double_t)this->GetTOFsignalTunedOnData((AliVTrack*)vtrack);
437   else tofTime=vtrack->GetTOFsignal();
438   tofTime = tofTime  - fTOFResponse.GetStartTime(vtrack->P());
439   Double_t delta=-9999.;
440
441   if (!ratio) delta=tofTime-expTime;
442   else if (expTime>1.e-20) delta=tofTime/expTime;
443   
444   return delta;
445 }
446
447 //_________________________________________________________________________
448 Float_t AliESDpid::GetNumberOfSigmasTOFold(const AliVParticle *track, AliPID::EParticleType type) const
449 {
450   //
451   // Number of sigma implementation for the TOF
452   //
453
454   AliVTrack *vtrack=(AliVTrack*)track;
455   Double_t tofTime = 99999;
456   if (fTuneMConData && ((fTuneMConDataMask & kDetTOF) == kDetTOF) ) tofTime = (Double_t)this->GetTOFsignalTunedOnData((AliVTrack*)vtrack);
457   else tofTime=vtrack->GetTOFsignal();
458   
459   Double_t expTime = fTOFResponse.GetExpectedSignal(vtrack,type);
460   return (tofTime - fTOFResponse.GetStartTime(vtrack->P()) - expTime)/fTOFResponse.GetExpectedSigma(vtrack->P(),expTime,AliPID::ParticleMassZ(type));
461 }