]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDpid.cxx
AliESDHeader: AliTriggerConfiguration and more trigger scalers added
[u/mrichter/AliRoot.git] / STEER / 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 "AliLog.h"
30 #include "AliPID.h"
31 #include "AliTOFHeader.h"
32 #include "AliESDpid.h"
33 #include "AliESDEvent.h"
34 #include "AliESDtrack.h"
35
36 ClassImp(AliESDpid)
37
38 Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF) const {
39   //
40   //  Calculate probabilities for all detectors, except if TPConly==kTRUE
41   //  and combine PID
42   //  
43   //   Option TPConly==kTRUE is used during reconstruction, 
44   //  because ITS tracking uses TPC pid
45   //  HMPID and TRD pid are done in detector reconstructors
46   //
47
48   /*
49   Float_t timeZeroTOF = 0;
50   if (subtractT0) 
51     timeZeroTOF = event->GetT0();
52   */
53   Int_t nTrk=event->GetNumberOfTracks();
54   for (Int_t iTrk=0; iTrk<nTrk; iTrk++) {  
55     AliESDtrack *track=event->GetTrack(iTrk);
56     MakeTPCPID(track);
57     if (!TPConly) {
58       MakeITSPID(track);
59       MakeTOFPID(track, timeZeroTOF);
60       //MakeHMPIDPID(track);
61       //MakeTRDPID(track);
62     }
63     CombinePID(track);
64   }
65   return 0;
66 }
67 //_________________________________________________________________________
68 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
69 {
70   //
71   //  TPC pid using bethe-bloch and gaussian response
72   //
73   if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
74     if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
75
76     Double_t mom = track->GetP();
77     const AliExternalTrackParam *in=track->GetInnerParam();
78     if (in) mom = in->GetP();
79
80     Double_t p[AliPID::kSPECIES];
81     Double_t dedx=track->GetTPCsignal(); 
82     Bool_t mismatch=kTRUE, heavy=kTRUE;
83
84     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
85       AliPID::EParticleType type=AliPID::EParticleType(j);
86       Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type); 
87       Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
88       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
89         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
90       } else {
91         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
92         mismatch=kFALSE;
93       }
94
95       // Check for particles heavier than (AliPID::kSPECIES - 1)
96       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
97
98     }
99
100     if (mismatch)
101        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
102
103     track->SetTPCpid(p);
104
105     if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
106
107 }
108 //_________________________________________________________________________
109 void AliESDpid::MakeITSPID(AliESDtrack *track) const
110 {
111   //
112   // ITS PID
113   // Two options, depending on fITSPIDmethod:
114   //  1) Truncated mean method
115   //  2) Likelihood, using charges measured in all 4 layers and 
116   //     Landau+gaus response functions
117   //
118
119   if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
120       (track->GetStatus()&AliESDtrack::kITSout)==0) return;
121
122   Double_t mom=track->GetP();  
123   if (fITSPIDmethod == kITSTruncMean) {
124     Double_t dedx=track->GetITSsignal();
125     Bool_t isSA=kTRUE;
126     Double_t momITS=mom;
127     ULong_t trStatus=track->GetStatus();
128     if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
129     UChar_t clumap=track->GetITSClusterMap();
130     Int_t nPointsForPid=0;
131     for(Int_t i=2; i<6; i++){
132       if(clumap&(1<<i)) ++nPointsForPid;
133     }
134
135     if(nPointsForPid<3) { // track not to be used for combined PID purposes
136       track->ResetStatus(AliESDtrack::kITSpid);
137       return;
138     }
139
140     Double_t p[10];
141
142     Bool_t mismatch=kTRUE, heavy=kTRUE;
143     for (Int_t j=0; j<AliPID::kSPECIES; j++) {
144       Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
145       Double_t bethe=fITSResponse.Bethe(momITS,mass);
146       Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
147       if (TMath::Abs(dedx-bethe) > fRange*sigma) {
148         p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
149       } else {
150         p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
151         mismatch=kFALSE;
152       }
153
154       // Check for particles heavier than (AliPID::kSPECIES - 1)
155       if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
156
157     }
158
159     if (mismatch)
160        for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
161
162     track->SetITSpid(p);
163
164     if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
165   }
166   else {  // Likelihood method
167     Double_t condprobfun[AliPID::kSPECIES];
168     Double_t qclu[4];
169     track->GetITSdEdxSamples(qclu);
170     fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
171     track->SetITSpid(condprobfun);
172   }
173
174 }
175 //_________________________________________________________________________
176 void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
177 {
178   //
179   //   TOF PID using gaussian response
180   //
181
182   if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
183   if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
184
185   Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
186   Float_t timezero = fTOFResponse.GetT0bin(ibin);
187
188   Double_t time[AliPID::kSPECIESN];
189   track->GetIntegratedTimes(time);
190
191   Double_t sigma[AliPID::kSPECIES];
192   for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
193     sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
194   }
195
196   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
197            Form("Expected TOF signals [ps]: %f %f %f %f %f",
198                   time[AliPID::kElectron],
199                   time[AliPID::kMuon],
200                   time[AliPID::kPion],
201                   time[AliPID::kKaon],
202                   time[AliPID::kProton]));
203
204   AliDebugGeneral("AliESDpid::MakeTOFPID",2,
205            Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
206                   sigma[AliPID::kElectron],
207                   sigma[AliPID::kMuon],
208                   sigma[AliPID::kPion],
209                   sigma[AliPID::kKaon],
210                   sigma[AliPID::kProton]
211                   ));
212
213   Double_t tof = track->GetTOFsignal() - timezero;
214
215   Double_t p[AliPID::kSPECIES];
216   Bool_t mismatch = kTRUE, heavy = kTRUE;
217   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
218     Double_t sig = sigma[j];
219     if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
220         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
221     } else
222       p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
223
224     // Check the mismatching
225     Double_t mass = AliPID::ParticleMass(j);
226     Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
227     if (p[j]>pm) mismatch = kFALSE;
228
229     // Check for particles heavier than (AliPID::kSPECIES - 1)
230     if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
231
232   }
233
234   if (mismatch)
235     for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
236
237   track->SetTOFpid(p);
238
239   if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);    
240   if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
241   
242 }
243 //_________________________________________________________________________
244 void AliESDpid::MakeTRDPID(AliESDtrack *track) const
245 {
246   //
247   // Method to recalculate the TRD PID probabilities
248   //
249   if((track->GetStatus()&AliESDtrack::kTRDout)==0) return;
250   Double_t prob[AliPID::kSPECIES]; Float_t mom[6];
251   Double_t dedx[48];  // Allocate space for the maximum number of TRD slices
252   for(Int_t ilayer = 0; ilayer < 6; ilayer++){
253     mom[ilayer] = track->GetTRDmomentum(ilayer);
254     for(Int_t islice = 0; islice < track->GetNumberOfTRDslices(); islice++){
255       dedx[ilayer*track->GetNumberOfTRDslices()+islice] = track->GetTRDslice(ilayer, islice);
256     }
257   }
258   fTRDResponse.GetResponse(track->GetNumberOfTRDslices(), dedx, mom, prob);
259   track->SetTRDpid(prob);
260 }
261 //_________________________________________________________________________
262 void AliESDpid::CombinePID(AliESDtrack *track) const
263 {
264   //
265   // Combine the information of various detectors
266   // to determine the Particle Identification
267   //
268   Int_t ns=AliPID::kSPECIES;
269   Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
270
271   if (track->IsOn(AliESDtrack::kITSpid)) {
272     Double_t d[10];
273     track->GetITSpid(d);
274     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
275   }
276
277   if (track->IsOn(AliESDtrack::kTPCpid)) {
278     Double_t d[10];
279     track->GetTPCpid(d);
280     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
281   }
282
283   if (track->IsOn(AliESDtrack::kTRDpid)) {
284     Double_t d[10];
285     track->GetTRDpid(d);
286     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
287   }
288
289   if (track->IsOn(AliESDtrack::kTOFpid)) {
290     Double_t d[10];
291     track->GetTOFpid(d);
292     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
293   }
294
295   if (track->IsOn(AliESDtrack::kHMPIDpid)) {
296     Double_t d[10];
297     track->GetHMPIDpid(d);
298     for (Int_t j=0; j<ns; j++) p[j]*=d[j];
299   }
300
301   track->SetESDpid(p);
302 }
303 //_________________________________________________________________________
304 Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
305     Bool_t status = kFALSE;
306     
307     Double_t exptimes[5];
308     track->GetIntegratedTimes(exptimes);
309     
310     Float_t p = track->P();
311     
312     Float_t dedx = track->GetTPCsignal();
313     Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
314     
315     Double_t ptpc[3];
316     track->GetInnerPxPyPz(ptpc);
317     Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
318     
319     for(Int_t i=0;i < 5;i++){
320         AliPID::EParticleType type=AliPID::EParticleType(i);
321         
322         Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
323         if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
324             Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
325             Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
326             
327             if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
328                 status = kTRUE;
329             }
330         }
331     }
332     
333     // for nuclei
334     Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
335     if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
336     
337     
338     return status;
339 }
340 //_________________________________________________________________________
341 void AliESDpid::SetTOFResponse(AliVEvent *vevent,EStartTimeType_t option){
342   //
343   // Set TOF response function
344   // Input option for event_time used
345   //
346
347     AliESDEvent *event=(AliESDEvent*)vevent;
348   
349     Float_t t0spread = 0.; //event->GetEventTimeSpread();
350     if(t0spread < 10) t0spread = 80;
351
352     // T0 from TOF algorithm
353
354     Bool_t flagT0TOF=kFALSE;
355     Bool_t flagT0T0=kFALSE;
356     Float_t *startTime = new Float_t[fTOFResponse.GetNmomBins()];
357     Float_t *startTimeRes = new Float_t[fTOFResponse.GetNmomBins()];
358     Int_t *startTimeMask = new Int_t[fTOFResponse.GetNmomBins()];
359
360     // T0-TOF arrays
361     Float_t *estimatedT0event = new Float_t[fTOFResponse.GetNmomBins()];
362     Float_t *estimatedT0resolution = new Float_t[fTOFResponse.GetNmomBins()];
363     for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
364       estimatedT0event[i]=0.0;
365       estimatedT0resolution[i]=0.0;
366       startTimeMask[i] = 0;
367     }
368
369     Float_t resT0A=75,resT0C=65,resT0AC=55;
370     if(event->GetT0TOF()){ // check if T0 detector information is available
371         flagT0T0=kTRUE;
372     }
373
374
375     AliTOFHeader *tofHeader =(AliTOFHeader*)event->GetTOFHeader();
376
377     if(tofHeader){ // read global info and T0-TOF info from ESD
378       fTOFResponse.SetTimeResolution(tofHeader->GetTOFResolution());
379       t0spread = tofHeader->GetT0spread(); // read t0 sprad
380       if(t0spread < 10) t0spread = 80;
381
382       flagT0TOF=kTRUE;
383       for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){ // read T0-TOF default value
384         startTime[i]=tofHeader->GetDefaultEventTimeVal();
385         startTimeRes[i]=tofHeader->GetDefaultEventTimeRes();
386         if(startTimeRes[i] < 1.e-5) startTimeRes[i] = t0spread;
387       }
388
389       TArrayI *ibin=tofHeader->GetNvalues();
390       TArrayF *t0Bin=tofHeader->GetEventTimeValues();
391       TArrayF *t0ResBin=tofHeader->GetEventTimeRes();
392       for(Int_t j=0;j < tofHeader->GetNbins();j++){ // fill T0-TOF in p-bins
393         Int_t icurrent = (Int_t)ibin->GetAt(j);
394         startTime[icurrent]=t0Bin->GetAt(j);
395         startTimeRes[icurrent]=t0ResBin->GetAt(j);
396         if(startTimeRes[icurrent] < 1.e-5) startTimeRes[icurrent] = t0spread;
397      }
398     }
399
400     // for cut of 3 sigma on t0 spread
401     Float_t t0cut = 3 * t0spread;
402     if(t0cut < 500) t0cut = 500;
403
404     if(option == kFILL_T0){ // T0-FILL is used
405         for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
406           estimatedT0event[i]=0.0;
407           estimatedT0resolution[i]=t0spread;
408         }
409         fTOFResponse.SetT0event(estimatedT0event);
410         fTOFResponse.SetT0resolution(estimatedT0resolution);
411     }
412
413     if(option == kTOF_T0){ // T0-TOF is used when available (T0-FILL otherwise) from ESD
414         if(flagT0TOF){
415             fTOFResponse.SetT0event(startTime);
416             fTOFResponse.SetT0resolution(startTimeRes);
417             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
418               if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
419               fTOFResponse.SetT0binMask(i,startTimeMask[i]);
420             }
421         }
422         else{
423             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
424               estimatedT0event[i]=0.0;
425               estimatedT0resolution[i]=t0spread;
426               fTOFResponse.SetT0binMask(i,startTimeMask[i]);
427             }
428             fTOFResponse.SetT0event(estimatedT0event);
429             fTOFResponse.SetT0resolution(estimatedT0resolution);
430         }
431     }
432     else if(option == kBest_T0){ // T0-T0 or T0-TOF are used when available (T0-FILL otherwise) from ESD
433         Float_t t0AC=-10000;
434         Float_t t0A=-10000;
435         Float_t t0C=-10000;
436         if(flagT0T0){
437             t0AC= event->GetT0TOF()[0];
438             t0A= event->GetT0TOF()[1];
439             t0C= event->GetT0TOF()[2];
440         }
441
442         Float_t t0t0Best = 0;
443         Float_t t0t0BestRes = 9999;
444         Int_t t0used=0;
445         if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
446             t0t0Best = t0AC;
447             t0t0BestRes = resT0AC;
448             t0used=6;
449         }
450         else if(TMath::Abs(t0C) < t0cut){
451             t0t0Best = t0C;
452             t0t0BestRes = resT0C;
453             t0used=4;
454         }
455         else if(TMath::Abs(t0A) < t0cut){
456             t0t0Best = t0A;
457             t0t0BestRes = resT0A;
458             t0used=2;
459         }
460
461         if(flagT0TOF){ // if T0-TOF info is available
462             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
463                 if(t0t0BestRes < 999){
464                   if(startTimeRes[i] < t0spread){
465                     Double_t wtot = 1./startTimeRes[i]/startTimeRes[i] + 1./t0t0BestRes/t0t0BestRes;
466                     Double_t t0best = startTime[i]/startTimeRes[i]/startTimeRes[i] + t0t0Best/t0t0BestRes/t0t0BestRes;
467                     estimatedT0event[i]=t0best / wtot;
468                     estimatedT0resolution[i]=1./TMath::Sqrt(wtot);
469                     startTimeMask[i] = t0used+1;
470                   }
471                   else {
472                     estimatedT0event[i]=t0t0Best;
473                     estimatedT0resolution[i]=t0t0BestRes;
474                     startTimeMask[i] = t0used;
475                   }
476                 }
477                 else{
478                   estimatedT0event[i]=startTime[i];
479                   estimatedT0resolution[i]=startTimeRes[i];
480                   if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
481                 }
482                 fTOFResponse.SetT0binMask(i,startTimeMask[i]);
483             }
484             fTOFResponse.SetT0event(estimatedT0event);
485             fTOFResponse.SetT0resolution(estimatedT0resolution);
486         }
487         else{ // if no T0-TOF info is available
488             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
489               fTOFResponse.SetT0binMask(i,t0used);
490               if(t0t0BestRes < 999){
491                 estimatedT0event[i]=t0t0Best;
492                 estimatedT0resolution[i]=t0t0BestRes;
493               }
494               else{
495                 estimatedT0event[i]=0.0;
496                 estimatedT0resolution[i]=t0spread;
497               }
498             }
499             fTOFResponse.SetT0event(estimatedT0event);
500             fTOFResponse.SetT0resolution(estimatedT0resolution);
501         }
502     }
503
504     else if(option == kT0_T0){ // T0-T0 is used when available (T0-FILL otherwise) from ESD
505         Float_t t0AC=-10000;
506         Float_t t0A=-10000;
507         Float_t t0C=-10000;
508         if(flagT0T0){
509             t0AC= event->GetT0TOF()[0];
510             t0A= event->GetT0TOF()[1];
511             t0C= event->GetT0TOF()[2];
512         }
513
514         if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
515             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
516               estimatedT0event[i]=t0AC;
517               estimatedT0resolution[i]=resT0AC;
518               fTOFResponse.SetT0binMask(i,6);
519             }
520         }
521         else if(TMath::Abs(t0C) < t0cut){
522             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
523               estimatedT0event[i]=t0C;
524               estimatedT0resolution[i]=resT0C;
525               fTOFResponse.SetT0binMask(i,4);
526             }
527         }
528         else if(TMath::Abs(t0A) < t0cut){
529             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
530               estimatedT0event[i]=t0A;
531               estimatedT0resolution[i]=resT0A;
532               fTOFResponse.SetT0binMask(i,2);
533             }
534         }
535         else{
536             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
537               estimatedT0event[i]=0.0;
538               estimatedT0resolution[i]=t0spread;
539               fTOFResponse.SetT0binMask(i,0);
540             }
541         }
542         fTOFResponse.SetT0event(estimatedT0event);
543         fTOFResponse.SetT0resolution(estimatedT0resolution);
544     }
545     delete [] startTime;
546     delete [] startTimeRes;
547     delete [] startTimeMask;
548     delete [] estimatedT0event;
549     delete [] estimatedT0resolution;
550 }