]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliESDpid.cxx
* add tuned on data functionality for TOF (non-gaussian tails)
[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         Double_t sigma = fTPCResponse.GetExpectedSigma(track, type, AliTPCPIDResponse::kdEdxDefault, this->UseTPCEtaCorrection());
128                 dedx = gRandom->Gaus(bethe,sigma);
129 //              if(iS == AliPID::ParticleCode(AliPID::kHe3) || iS == AliPID::ParticleCode(AliPID::kAlpha)) dedx *= 5;
130             }
131         }
132     }
133
134     track->SetTPCsignalTunedOnData(dedx);
135     return dedx;
136 }
137 //_________________________________________________________________________
138 Float_t AliESDpid::GetTOFsignalTunedOnData(const AliVTrack *t) const {
139     AliESDtrack *track = (AliESDtrack *) t;
140     Double_t tofSignal = track->GetTOFsignalTunedOnData();
141
142     if(tofSignal <  99999) return (Float_t)tofSignal; // it has been already set
143
144     tofSignal = t->GetTOFsignal() + fTOFResponse.GetTailRandomValue();
145     track->SetTOFsignalTunedOnData(tofSignal);
146     return (Float_t)tofSignal;
147 }
148 //_________________________________________________________________________
149 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
150 {
151   //
152   //  TPC pid using bethe-bloch and gaussian response
153   //
154   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
155     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
156
157     Double_t mom = track->GetP();
158     const AliExternalTrackParam *in=track->GetInnerParam();
159     if (in) mom = in->GetP();
160
161     Double_t p[AliPID::kSPECIES];
162     Double_t dedx=track->GetTPCsignal(); 
163     Bool_t mismatch=kTRUE, heavy=kTRUE;
164
165     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
166       AliPID::EParticleType type=AliPID::EParticleType(j);
167       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
168       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
169       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
170         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
171       } else {
172         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
173         mismatch=kFALSE;
174       }
175
176       // Check for particles heavier than (AliPID::kSPECIES - 1)
177       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
178
179     }
180
181     if (mismatch)
182        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
183
184     track->SetTPCpid(p);
185
186     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
187
188 }
189 //_________________________________________________________________________
190 void AliESDpid::MakeITSPID(AliESDtrack *track) const
191 {
192   //
193   // ITS PID
194   // Two options, depending on fITSPIDmethod:
195   //  1) Truncated mean method
196   //  2) Likelihood, using charges measured in all 4 layers and 
197   //     Landau+gaus response functions
198   //
199
200   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
201       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
202
203   Double_t mom=track->GetP();  
204   if (fITSPIDmethod == kITSTruncMean) {
205     Double_t dedx=track->GetITSsignal();
206     Bool_t isSA=kTRUE;
207     Double_t momITS=mom;
208     ULong_t trStatus=track->GetStatus();
209     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
210     UChar_t clumap=track->GetITSClusterMap();
211     Int_t nPointsForPid=0;
212     for(Int_t i=2; i<6; i++){
213       if(clumap&(1<<i)) ++nPointsForPid;
214     }
215
216     if(nPointsForPid<3) { // track not to be used for combined PID purposes
217       track->ResetStatus(AliESDtrack::kITSpid);
218       return;
219     }
220
221     Double_t p[10];
222
223     Bool_t mismatch=kTRUE, heavy=kTRUE;
224     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
225       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
226       Double_t bethe=fITSResponse.Bethe(momITS,mass);
227       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
228       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
229         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
230       } else {
231         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
232         mismatch=kFALSE;
233       }
234
235       // Check for particles heavier than (AliPID::kSPECIES - 1)
236       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
237
238     }
239
240     if (mismatch)
241        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
242
243     track->SetITSpid(p);
244
245     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
246   }
247   else {  // Likelihood method
248     Double_t condprobfun[AliPID::kSPECIES];
249     Double_t qclu[4];
250     track->GetITSdEdxSamples(qclu);
251     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
252     track->SetITSpid(condprobfun);
253   }
254
255 }
256 //_________________________________________________________________________
257 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
258 {
259   //
260   //   TOF PID using gaussian response
261   //
262
263   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
264   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
265   if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
266
267   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
268   Float_t timezero = fTOFResponse.GetT0bin(ibin);
269
270   Double_t time[AliPID::kSPECIES];
271   track->GetIntegratedTimes(time);
272
273   Double_t sigma[AliPID::kSPECIES];
274   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
275     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
276   }
277
278   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
279            Form("Expected TOF signals [ps]: %f %f %f %f %f",
280                   time[AliPID::kElectron],
281                   time[AliPID::kMuon],
282                   time[AliPID::kPion],
283                   time[AliPID::kKaon],
284                   time[AliPID::kProton]));
285
286   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
287            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
288                   sigma[AliPID::kElectron],
289                   sigma[AliPID::kMuon],
290                   sigma[AliPID::kPion],
291                   sigma[AliPID::kKaon],
292                   sigma[AliPID::kProton]
293                   ));
294
295   Double_t tof = track->GetTOFsignal() - timezero;
296
297   Double_t p[AliPID::kSPECIES];
298 //   Bool_t mismatch = kTRUE;
299   Bool_t heavy = kTRUE;
300   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
301     Double_t sig = sigma[j];
302     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
303         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
304     } else
305       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
306
307     // Check the mismatching
308     
309 //     Double_t mass = AliPID::ParticleMass(j);
310 //     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
311 //     if (p[j]>pm) mismatch = kFALSE;
312
313     // Check for particles heavier than (AliPID::kSPECIES - 1)
314     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
315
316   }
317
318   /*
319     if (mismatch)
320     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
321   */
322
323   track->SetTOFpid(p);
324
325   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
326   
327   // kTOFmismatch flas is not set because deprecated from 18/02/2013
328   //  if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
329   //  else track->ResetStatus(AliESDtrack::kTOFmismatch);
330 }
331 //_________________________________________________________________________
332 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
333 {
334   //
335   // Method to recalculate the TRD PID probabilities
336   //
337   Double_t prob[AliPID::kSPECIES];
338   GetComputeTRDProbability(track, AliPID::kSPECIES, prob);
339   track->SetTRDpid(prob);
340 }
341 //_________________________________________________________________________
342 void AliESDpid::CombinePID(AliESDtrack *track) const
343 {
344   //
345   // Combine the information of various detectors
346   // to determine the Particle Identification
347   //
348   Int_t ns=AliPID::kSPECIES;
349   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
350
351   if (track->IsOn(AliESDtrack::kITSpid)) {
352     Double_t d[10];
353     track->GetITSpid(d);
354     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
355   }
356
357   if (track->IsOn(AliESDtrack::kTPCpid)) {
358     Double_t d[10];
359     track->GetTPCpid(d);
360     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
361   }
362
363   if (track->IsOn(AliESDtrack::kTRDpid)) {
364     Double_t d[10];
365     track->GetTRDpid(d);
366     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
367   }
368
369   if (track->IsOn(AliESDtrack::kTOFpid)) {
370     Double_t d[10];
371     track->GetTOFpid(d);
372     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
373   }
374
375   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
376     Double_t d[10];
377     track->GetHMPIDpid(d);
378     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
379   }
380
381   track->SetESDpid(p);
382 }
383 //_________________________________________________________________________
384 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
385   //
386   // Check pid matching of TOF with TPC as reference
387   //
388     Bool_t status = kFALSE;
389     
390     Double_t exptimes[5];
391     track->GetIntegratedTimes(exptimes);
392     
393     Float_t p = track->P();
394     
395     Float_t dedx = track->GetTPCsignal();
396     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
397     
398     Double_t ptpc[3];
399     track->GetInnerPxPyPz(ptpc);
400     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
401     
402     for(Int_t i=0;i < 5;i++){
403         AliPID::EParticleType type=AliPID::EParticleType(i);
404         
405         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
406         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
407             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
408             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
409             
410             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
411                 status = kTRUE;
412             }
413         }
414     }
415     
416     // for nuclei
417     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
418     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
419     
420     
421     return status;
422 }
423
424 //_________________________________________________________________________
425 Float_t AliESDpid::GetSignalDeltaTOFold(const AliVParticle *track, AliPID::EParticleType type, Bool_t ratio/*=kFALSE*/) const
426 {
427   //
428   // TOF signal - expected
429   //
430   AliVTrack *vtrack=(AliVTrack*)track;
431   
432   const Double_t expTime = fTOFResponse.GetExpectedSignal(vtrack,type);
433   Double_t tofTime = 99999;
434   if (fTuneMConData && ((fTuneMConDataMask & kDetTOF) == kDetTOF) ) tofTime = (Double_t)this->GetTOFsignalTunedOnData((AliVTrack*)vtrack);
435   else tofTime=vtrack->GetTOFsignal();
436   tofTime = tofTime  - fTOFResponse.GetStartTime(vtrack->P());
437   Double_t delta=-9999.;
438
439   if (!ratio) delta=tofTime-expTime;
440   else if (expTime>1.e-20) delta=tofTime/expTime;
441   
442   return delta;
443 }
444
445 //_________________________________________________________________________
446 Float_t AliESDpid::GetNumberOfSigmasTOFold(const AliVParticle *track, AliPID::EParticleType type) const
447 {
448   //
449   // Number of sigma implementation for the TOF
450   //
451
452   AliVTrack *vtrack=(AliVTrack*)track;
453   Double_t tofTime = 99999;
454   if (fTuneMConData && ((fTuneMConDataMask & kDetTOF) == kDetTOF) ) tofTime = (Double_t)this->GetTOFsignalTunedOnData((AliVTrack*)vtrack);
455   else tofTime=vtrack->GetTOFsignal();
456   
457   Double_t expTime = fTOFResponse.GetExpectedSignal(vtrack,type);
458   return (tofTime - fTOFResponse.GetStartTime(vtrack->P()) - expTime)/fTOFResponse.GetExpectedSigma(vtrack->P(),expTime,AliPID::ParticleMassZ(type));
459 }