]>
Commit | Line | Data |
---|---|---|
8c6a71ab | 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 | //----------------------------------------------------------------- | |
17 | // Implementation of the TPC PID class | |
18 | // Very naive one... Should be made better by the detector experts... | |
19 | // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch | |
20 | //----------------------------------------------------------------- | |
21 | ||
22 | #include "AliTPCpidESD.h" | |
af885e0f | 23 | #include "AliESDEvent.h" |
8c6a71ab | 24 | #include "AliESDtrack.h" |
25 | ||
26 | ClassImp(AliTPCpidESD) | |
27 | ||
28 | //_________________________________________________________________________ | |
179c6296 | 29 | AliTPCpidESD::AliTPCpidESD(Double_t *param): |
30 | fMIP(0.), | |
31 | fRes(0.), | |
32 | fRange(0.) | |
8c6a71ab | 33 | { |
34 | // | |
35 | // The main constructor | |
36 | // | |
37 | fMIP=param[0]; | |
38 | fRes=param[1]; | |
39 | fRange=param[2]; | |
40 | } | |
41 | ||
42 | Double_t AliTPCpidESD::Bethe(Double_t bg) { | |
43 | // | |
44 | // This is the Bethe-Bloch function normalised to 1 at the minimum | |
45 | // | |
46 | Double_t bg2=bg*bg; | |
cf6a95b7 | 47 | Double_t bethe; |
48 | if (bg<3.5e1) | |
49 | bethe=(1.+ bg2)/bg2*(log(5940*bg2) - bg2/(1.+ bg2)); | |
50 | else // Density effect ( approximately :) | |
51 | bethe=1.15*(1.+ bg2)/bg2*(log(3.5*5940*bg) - bg2/(1.+ bg2)); | |
8c6a71ab | 52 | return bethe/11.091; |
53 | } | |
54 | ||
55 | //_________________________________________________________________________ | |
af885e0f | 56 | Int_t AliTPCpidESD::MakePID(AliESDEvent *event) |
8c6a71ab | 57 | { |
58 | // | |
59 | // This function calculates the "detector response" PID probabilities | |
60 | // | |
8c6a71ab | 61 | Int_t ntrk=event->GetNumberOfTracks(); |
62 | for (Int_t i=0; i<ntrk; i++) { | |
63 | AliESDtrack *t=event->GetTrack(i); | |
64 | if ((t->GetStatus()&AliESDtrack::kTPCin )==0) | |
65 | if ((t->GetStatus()&AliESDtrack::kTPCout)==0) continue; | |
8c6a71ab | 66 | Double_t p[10]; |
7a8614f3 | 67 | Double_t mom=t->GetP(); |
68 | const AliExternalTrackParam *in=t->GetInnerParam(); | |
69 | if (in) mom=in->GetP(); | |
70 | Double_t dedx=t->GetTPCsignal()/fMIP; | |
71 | Bool_t mismatch=kTRUE, heavy=kTRUE; | |
72 | for (Int_t j=0; j<AliPID::kSPECIES; j++) { | |
304864ab | 73 | Double_t mass=AliPID::ParticleMass(j); |
8c6a71ab | 74 | Double_t bethe=Bethe(mom/mass); |
75 | Double_t sigma=fRes*bethe; | |
76 | if (TMath::Abs(dedx-bethe) > fRange*sigma) { | |
431be10d | 77 | p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma; |
7a8614f3 | 78 | } else { |
79 | p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma; | |
80 | mismatch=kFALSE; | |
8c6a71ab | 81 | } |
7a8614f3 | 82 | |
83 | // Check for particles heavier than (AliPID::kSPECIES - 1) | |
84 | if (dedx < (bethe + fRange*sigma)) heavy=kFALSE; | |
85 | ||
8c6a71ab | 86 | } |
7a8614f3 | 87 | |
88 | if (mismatch) | |
89 | for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1/AliPID::kSPECIES; | |
90 | ||
8c6a71ab | 91 | t->SetTPCpid(p); |
7a8614f3 | 92 | |
93 | if (heavy) t->ResetStatus(AliESDtrack::kTPCpid); | |
94 | ||
8c6a71ab | 95 | } |
96 | return 0; | |
97 | } |