]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliESDpid.cxx
Added basic locking for the analysis manager session and data handlers to protect...
[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
39 ClassImp(AliESDpid)
40
41 Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF) const {
42   //
43   //  Calculate probabilities for all detectors, except if TPConly==kTRUE
44   //  and combine PID
45   //  
46   //   Option TPConly==kTRUE is used during reconstruction, 
47   //  because ITS tracking uses TPC pid
48   //  HMPID and TRD pid are done in detector reconstructors
49   //
50
51   /*
52   Float_t timeZeroTOF = 0;
53   if (subtractT0) 
54     timeZeroTOF = event->GetT0();
55   */
56   Int_t nTrk=event->GetNumberOfTracks();
57   for (Int_t iTrk=0; iTrk<nTrk; iTrk++) {  
58     AliESDtrack *track=event->GetTrack(iTrk);
59     MakeTPCPID(track);
60     if (!TPConly) {
61       MakeITSPID(track);
62       MakeTOFPID(track, timeZeroTOF);
63       //MakeHMPIDPID(track);
64       //MakeTRDPID(track);
65     }
66     CombinePID(track);
67   }
68   return 0;
69 }
70 //_________________________________________________________________________
71 Float_t AliESDpid::GetTPCsignalTunedOnData(const AliVTrack *t) const {
72     AliESDtrack *track = (AliESDtrack *) t;
73     Float_t dedx = track->GetTPCsignalTunedOnData();
74
75     if(dedx > 0) return dedx;
76
77     Double_t mom = t->GetTPCmomentum();
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                 Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
126                 Double_t sigma=fTPCResponse.GetExpectedSigma(mom,t->GetTPCsignalN(),type);
127                 dedx = gRandom->Gaus(bethe,sigma);
128                 if(iS == AliPID::ParticleCode(AliPID::kHe3) || iS == AliPID::ParticleCode(AliPID::kAlpha)) dedx *= 5;
129             }
130         }
131     }
132
133     track->SetTPCsignalTunedOnData(dedx);
134     return dedx;
135 }
136 //_________________________________________________________________________
137 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
138 {
139   //
140   //  TPC pid using bethe-bloch and gaussian response
141   //
142   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
143     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
144
145     Double_t mom = track->GetP();
146     const AliExternalTrackParam *in=track->GetInnerParam();
147     if (in) mom = in->GetP();
148
149     Double_t p[AliPID::kSPECIES];
150     Double_t dedx=track->GetTPCsignal(); 
151     Bool_t mismatch=kTRUE, heavy=kTRUE;
152
153     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
154       AliPID::EParticleType type=AliPID::EParticleType(j);
155       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
156       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
157       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
158         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
159       } else {
160         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
161         mismatch=kFALSE;
162       }
163
164       // Check for particles heavier than (AliPID::kSPECIES - 1)
165       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
166
167     }
168
169     if (mismatch)
170        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
171
172     track->SetTPCpid(p);
173
174     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
175
176 }
177 //_________________________________________________________________________
178 void AliESDpid::MakeITSPID(AliESDtrack *track) const
179 {
180   //
181   // ITS PID
182   // Two options, depending on fITSPIDmethod:
183   //  1) Truncated mean method
184   //  2) Likelihood, using charges measured in all 4 layers and 
185   //     Landau+gaus response functions
186   //
187
188   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
189       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
190
191   Double_t mom=track->GetP();  
192   if (fITSPIDmethod == kITSTruncMean) {
193     Double_t dedx=track->GetITSsignal();
194     Bool_t isSA=kTRUE;
195     Double_t momITS=mom;
196     ULong_t trStatus=track->GetStatus();
197     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
198     UChar_t clumap=track->GetITSClusterMap();
199     Int_t nPointsForPid=0;
200     for(Int_t i=2; i<6; i++){
201       if(clumap&(1<<i)) ++nPointsForPid;
202     }
203
204     if(nPointsForPid<3) { // track not to be used for combined PID purposes
205       track->ResetStatus(AliESDtrack::kITSpid);
206       return;
207     }
208
209     Double_t p[10];
210
211     Bool_t mismatch=kTRUE, heavy=kTRUE;
212     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
213       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
214       Double_t bethe=fITSResponse.Bethe(momITS,mass);
215       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
216       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
217         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
218       } else {
219         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
220         mismatch=kFALSE;
221       }
222
223       // Check for particles heavier than (AliPID::kSPECIES - 1)
224       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
225
226     }
227
228     if (mismatch)
229        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
230
231     track->SetITSpid(p);
232
233     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
234   }
235   else {  // Likelihood method
236     Double_t condprobfun[AliPID::kSPECIES];
237     Double_t qclu[4];
238     track->GetITSdEdxSamples(qclu);
239     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
240     track->SetITSpid(condprobfun);
241   }
242
243 }
244 //_________________________________________________________________________
245 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
246 {
247   //
248   //   TOF PID using gaussian response
249   //
250
251   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
252   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
253   if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
254
255   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
256   Float_t timezero = fTOFResponse.GetT0bin(ibin);
257
258   Double_t time[AliPID::kSPECIESN];
259   track->GetIntegratedTimes(time);
260
261   Double_t sigma[AliPID::kSPECIES];
262   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
263     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
264   }
265
266   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
267            Form("Expected TOF signals [ps]: %f %f %f %f %f",
268                   time[AliPID::kElectron],
269                   time[AliPID::kMuon],
270                   time[AliPID::kPion],
271                   time[AliPID::kKaon],
272                   time[AliPID::kProton]));
273
274   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
275            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
276                   sigma[AliPID::kElectron],
277                   sigma[AliPID::kMuon],
278                   sigma[AliPID::kPion],
279                   sigma[AliPID::kKaon],
280                   sigma[AliPID::kProton]
281                   ));
282
283   Double_t tof = track->GetTOFsignal() - timezero;
284
285   Double_t p[AliPID::kSPECIES];
286   Bool_t mismatch = kTRUE, heavy = kTRUE;
287   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
288     Double_t sig = sigma[j];
289     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
290         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
291     } else
292       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
293
294     // Check the mismatching
295     
296     Double_t mass = AliPID::ParticleMass(j);
297     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
298     if (p[j]>pm) mismatch = kFALSE;
299
300     // Check for particles heavier than (AliPID::kSPECIES - 1)
301     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
302
303   }
304
305   /*
306     if (mismatch)
307     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
308   */
309
310   track->SetTOFpid(p);
311
312   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
313   
314   if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
315   else track->ResetStatus(AliESDtrack::kTOFmismatch);
316 }
317 //_________________________________________________________________________
318 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
319 {
320   //
321   // Method to recalculate the TRD PID probabilities
322   //
323   Double_t prob[AliPID::kSPECIES];
324   ComputeTRDProbability(track, AliPID::kSPECIES, prob);
325   track->SetTRDpid(prob);
326 }
327 //_________________________________________________________________________
328 void AliESDpid::CombinePID(AliESDtrack *track) const
329 {
330   //
331   // Combine the information of various detectors
332   // to determine the Particle Identification
333   //
334   Int_t ns=AliPID::kSPECIES;
335   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
336
337   if (track->IsOn(AliESDtrack::kITSpid)) {
338     Double_t d[10];
339     track->GetITSpid(d);
340     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
341   }
342
343   if (track->IsOn(AliESDtrack::kTPCpid)) {
344     Double_t d[10];
345     track->GetTPCpid(d);
346     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
347   }
348
349   if (track->IsOn(AliESDtrack::kTRDpid)) {
350     Double_t d[10];
351     track->GetTRDpid(d);
352     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
353   }
354
355   if (track->IsOn(AliESDtrack::kTOFpid)) {
356     Double_t d[10];
357     track->GetTOFpid(d);
358     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
359   }
360
361   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
362     Double_t d[10];
363     track->GetHMPIDpid(d);
364     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
365   }
366
367   track->SetESDpid(p);
368 }
369 //_________________________________________________________________________
370 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
371   //
372   // Check pid matching of TOF with TPC as reference
373   //
374     Bool_t status = kFALSE;
375     
376     Double_t exptimes[5];
377     track->GetIntegratedTimes(exptimes);
378     
379     Float_t p = track->P();
380     
381     Float_t dedx = track->GetTPCsignal();
382     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
383     
384     Double_t ptpc[3];
385     track->GetInnerPxPyPz(ptpc);
386     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
387     
388     for(Int_t i=0;i < 5;i++){
389         AliPID::EParticleType type=AliPID::EParticleType(i);
390         
391         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
392         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
393             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
394             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
395             
396             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
397                 status = kTRUE;
398             }
399         }
400     }
401     
402     // for nuclei
403     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
404     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
405     
406     
407     return status;
408 }