]>
Commit | Line | Data |
---|---|---|
1 | /************************************************************************** | |
2 | * Copyright(c) 2005-2007, 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 | // ITS PID method # 1 | |
20 | // Implementation of the ITS PID class | |
21 | // Very naive one... Should be made better by the detector experts... | |
22 | // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch | |
23 | //----------------------------------------------------------------- | |
24 | #include "AliITSpidESD.h" | |
25 | #include "AliITSpidESD1.h" | |
26 | #include "AliESDEvent.h" | |
27 | #include "AliESDtrack.h" | |
28 | ||
29 | ClassImp(AliITSpidESD1) | |
30 | ||
31 | AliITSpidESD1::AliITSpidESD1(): AliITSpidESD(), | |
32 | fMIP(0), | |
33 | fRes(0), | |
34 | fRange(0) | |
35 | { | |
36 | //Default constructor | |
37 | } | |
38 | //_________________________________________________________________________ | |
39 | AliITSpidESD1::AliITSpidESD1(Double_t *param): AliITSpidESD(), | |
40 | fMIP(param[0]), | |
41 | fRes(param[1]), | |
42 | fRange(param[2]) | |
43 | { | |
44 | // | |
45 | // The main constructor | |
46 | // | |
47 | } | |
48 | ||
49 | ||
50 | //_________________________________________________________________________ | |
51 | Int_t AliITSpidESD1::MakePID(AliESDEvent *event) | |
52 | { | |
53 | // | |
54 | // This function calculates the "detector response" PID probabilities | |
55 | // | |
56 | Int_t ntrk=event->GetNumberOfTracks(); | |
57 | for (Int_t i=0; i<ntrk; i++) { | |
58 | AliESDtrack *t=event->GetTrack(i); | |
59 | if ((t->GetStatus()&AliESDtrack::kITSin )==0) | |
60 | if ((t->GetStatus()&AliESDtrack::kITSout)==0) continue; | |
61 | Double_t mom=t->GetP(); | |
62 | Double_t dedx=t->GetITSsignal()/fMIP; | |
63 | Double_t p[10]; | |
64 | Bool_t mismatch=kTRUE, heavy=kTRUE; | |
65 | for (Int_t j=0; j<AliPID::kSPECIES; j++) { | |
66 | Double_t mass=AliPID::ParticleMass(j);//GeV/c^2 | |
67 | Double_t bethe=AliITSpidESD::Bethe(mom,mass); | |
68 | Double_t sigma=fRes*bethe; | |
69 | if (TMath::Abs(dedx-bethe) > fRange*sigma) { | |
70 | p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma; | |
71 | } else { | |
72 | p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma; | |
73 | mismatch=kFALSE; | |
74 | } | |
75 | ||
76 | // Check for particles heavier than (AliPID::kSPECIES - 1) | |
77 | if (dedx < (bethe + fRange*sigma)) heavy=kFALSE; | |
78 | ||
79 | } | |
80 | ||
81 | if (mismatch) | |
82 | for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1/AliPID::kSPECIES; | |
83 | ||
84 | t->SetITSpid(p); | |
85 | ||
86 | if (heavy) t->ResetStatus(AliESDtrack::kITSpid); | |
87 | ||
88 | } | |
89 | return 0; | |
90 | } |