]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliESDpid.cxx
better just warning than error in case of mirrors requesting when putting to local...
[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     Double_t mom = t->GetTPCmomentum();
79     dedx = t->GetTPCsignal();
80     track->SetTPCsignalTunedOnData(dedx);
81     if(dedx < 20) return dedx;
82
83     AliPID::EParticleType type = AliPID::kPion;
84
85     AliMCEventHandler* eventHandler=dynamic_cast<AliMCEventHandler*>(fEventHandler);
86     if (eventHandler) {
87         AliMCEvent* mcEvent = eventHandler->MCEvent();
88         if(mcEvent){
89             Bool_t kGood = kTRUE;
90             AliMCParticle *MCpart = (AliMCParticle *) mcEvent->GetTrack(TMath::Abs(t->GetLabel()));
91             TParticle *part = MCpart->Particle();
92             
93             Int_t iS = TMath::Abs(part->GetPdgCode());
94
95             if(iS==AliPID::ParticleCode(AliPID::kElectron)){
96                 type = AliPID::kElectron;
97             }
98             else if(iS==AliPID::ParticleCode(AliPID::kMuon)){
99                 type = AliPID::kMuon;
100             }
101             else if(iS==AliPID::ParticleCode(AliPID::kPion)){
102                 type = AliPID::kPion;
103             }
104             else if(iS==AliPID::ParticleCode(AliPID::kKaon)){
105                 type = AliPID::kKaon;
106             }
107             else if(iS==AliPID::ParticleCode(AliPID::kProton)){
108                 type = AliPID::kProton;
109             }
110             else if(iS==AliPID::ParticleCode(AliPID::kDeuteron)){ // d
111                 type = AliPID::kDeuteron;
112             }
113             else if(iS==AliPID::ParticleCode(AliPID::kTriton)){ // t
114                 type = AliPID::kTriton;
115             }
116             else if(iS==AliPID::ParticleCode(AliPID::kHe3)){ // 3He
117                 type = AliPID::kHe3;
118             }
119             else if(iS==AliPID::ParticleCode(AliPID::kAlpha)){ // 4He
120                 type = AliPID::kAlpha;
121             }
122             else
123                 kGood = kFALSE;
124
125             if(kGood){
126                 Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
127                 Double_t sigma=fTPCResponse.GetExpectedSigma(mom,t->GetTPCsignalN(),type);
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 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
139 {
140   //
141   //  TPC pid using bethe-bloch and gaussian response
142   //
143   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
144     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
145
146     Double_t mom = track->GetP();
147     const AliExternalTrackParam *in=track->GetInnerParam();
148     if (in) mom = in->GetP();
149
150     Double_t p[AliPID::kSPECIES];
151     Double_t dedx=track->GetTPCsignal(); 
152     Bool_t mismatch=kTRUE, heavy=kTRUE;
153
154     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
155       AliPID::EParticleType type=AliPID::EParticleType(j);
156       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
157       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
158       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
159         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
160       } else {
161         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
162         mismatch=kFALSE;
163       }
164
165       // Check for particles heavier than (AliPID::kSPECIES - 1)
166       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
167
168     }
169
170     if (mismatch)
171        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
172
173     track->SetTPCpid(p);
174
175     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
176
177 }
178 //_________________________________________________________________________
179 void AliESDpid::MakeITSPID(AliESDtrack *track) const
180 {
181   //
182   // ITS PID
183   // Two options, depending on fITSPIDmethod:
184   //  1) Truncated mean method
185   //  2) Likelihood, using charges measured in all 4 layers and 
186   //     Landau+gaus response functions
187   //
188
189   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
190       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
191
192   Double_t mom=track->GetP();  
193   if (fITSPIDmethod == kITSTruncMean) {
194     Double_t dedx=track->GetITSsignal();
195     Bool_t isSA=kTRUE;
196     Double_t momITS=mom;
197     ULong_t trStatus=track->GetStatus();
198     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
199     UChar_t clumap=track->GetITSClusterMap();
200     Int_t nPointsForPid=0;
201     for(Int_t i=2; i<6; i++){
202       if(clumap&(1<<i)) ++nPointsForPid;
203     }
204
205     if(nPointsForPid<3) { // track not to be used for combined PID purposes
206       track->ResetStatus(AliESDtrack::kITSpid);
207       return;
208     }
209
210     Double_t p[10];
211
212     Bool_t mismatch=kTRUE, heavy=kTRUE;
213     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
214       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
215       Double_t bethe=fITSResponse.Bethe(momITS,mass);
216       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
217       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
218         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
219       } else {
220         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
221         mismatch=kFALSE;
222       }
223
224       // Check for particles heavier than (AliPID::kSPECIES - 1)
225       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
226
227     }
228
229     if (mismatch)
230        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
231
232     track->SetITSpid(p);
233
234     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
235   }
236   else {  // Likelihood method
237     Double_t condprobfun[AliPID::kSPECIES];
238     Double_t qclu[4];
239     track->GetITSdEdxSamples(qclu);
240     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
241     track->SetITSpid(condprobfun);
242   }
243
244 }
245 //_________________________________________________________________________
246 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
247 {
248   //
249   //   TOF PID using gaussian response
250   //
251
252   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
253   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
254   if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
255
256   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
257   Float_t timezero = fTOFResponse.GetT0bin(ibin);
258
259   Double_t time[AliPID::kSPECIES];
260   track->GetIntegratedTimes(time);
261
262   Double_t sigma[AliPID::kSPECIES];
263   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
264     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
265   }
266
267   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
268            Form("Expected TOF signals [ps]: %f %f %f %f %f",
269                   time[AliPID::kElectron],
270                   time[AliPID::kMuon],
271                   time[AliPID::kPion],
272                   time[AliPID::kKaon],
273                   time[AliPID::kProton]));
274
275   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
276            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
277                   sigma[AliPID::kElectron],
278                   sigma[AliPID::kMuon],
279                   sigma[AliPID::kPion],
280                   sigma[AliPID::kKaon],
281                   sigma[AliPID::kProton]
282                   ));
283
284   Double_t tof = track->GetTOFsignal() - timezero;
285
286   Double_t p[AliPID::kSPECIES];
287   Bool_t mismatch = kTRUE, heavy = kTRUE;
288   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
289     Double_t sig = sigma[j];
290     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
291         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
292     } else
293       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
294
295     // Check the mismatching
296     
297     Double_t mass = AliPID::ParticleMass(j);
298     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
299     if (p[j]>pm) mismatch = kFALSE;
300
301     // Check for particles heavier than (AliPID::kSPECIES - 1)
302     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
303
304   }
305
306   /*
307     if (mismatch)
308     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
309   */
310
311   track->SetTOFpid(p);
312
313   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
314   
315   if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
316   else track->ResetStatus(AliESDtrack::kTOFmismatch);
317 }
318 //_________________________________________________________________________
319 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
320 {
321   //
322   // Method to recalculate the TRD PID probabilities
323   //
324   Double_t prob[AliPID::kSPECIES];
325   ComputeTRDProbability(track, AliPID::kSPECIES, prob);
326   track->SetTRDpid(prob);
327 }
328 //_________________________________________________________________________
329 void AliESDpid::CombinePID(AliESDtrack *track) const
330 {
331   //
332   // Combine the information of various detectors
333   // to determine the Particle Identification
334   //
335   Int_t ns=AliPID::kSPECIES;
336   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
337
338   if (track->IsOn(AliESDtrack::kITSpid)) {
339     Double_t d[10];
340     track->GetITSpid(d);
341     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
342   }
343
344   if (track->IsOn(AliESDtrack::kTPCpid)) {
345     Double_t d[10];
346     track->GetTPCpid(d);
347     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
348   }
349
350   if (track->IsOn(AliESDtrack::kTRDpid)) {
351     Double_t d[10];
352     track->GetTRDpid(d);
353     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
354   }
355
356   if (track->IsOn(AliESDtrack::kTOFpid)) {
357     Double_t d[10];
358     track->GetTOFpid(d);
359     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
360   }
361
362   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
363     Double_t d[10];
364     track->GetHMPIDpid(d);
365     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
366   }
367
368   track->SetESDpid(p);
369 }
370 //_________________________________________________________________________
371 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
372   //
373   // Check pid matching of TOF with TPC as reference
374   //
375     Bool_t status = kFALSE;
376     
377     Double_t exptimes[5];
378     track->GetIntegratedTimes(exptimes);
379     
380     Float_t p = track->P();
381     
382     Float_t dedx = track->GetTPCsignal();
383     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
384     
385     Double_t ptpc[3];
386     track->GetInnerPxPyPz(ptpc);
387     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
388     
389     for(Int_t i=0;i < 5;i++){
390         AliPID::EParticleType type=AliPID::EParticleType(i);
391         
392         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
393         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
394             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
395             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
396             
397             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
398                 status = kTRUE;
399             }
400         }
401     }
402     
403     // for nuclei
404     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
405     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
406     
407     
408     return status;
409 }
410
411 Float_t AliESDpid::NumberOfSigmasTOF(const AliVParticle *track, AliPID::EParticleType type, const Float_t /*timeZeroTOF*/) const
412 {
413   //
414   // Number of sigma implementation for the TOF
415   //
416   
417   AliVTrack *vtrack=(AliVTrack*)track;
418   // look for cached value first
419   if (vtrack->GetDetectorPID()){
420     return vtrack->GetDetectorPID()->GetNumberOfSigmas(kTOF, type);
421   }
422   
423   Double_t expTime = fTOFResponse.GetExpectedSignal(vtrack,type);
424   return (vtrack->GetTOFsignal() - fTOFResponse.GetStartTime(vtrack->P()) - expTime)/fTOFResponse.GetExpectedSigma(vtrack->P(),expTime,AliPID::ParticleMassZ(type));
425 }