]>
Commit | Line | Data |
---|---|---|
809a4336 | 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 | #include <TH2F.h> | |
17 | #include <TList.h> | |
18 | #include <TMath.h> | |
19 | ||
722347d8 | 20 | #include "AliAODTrack.h" |
21 | #include "AliAODMCParticle.h" | |
809a4336 | 22 | #include "AliESDtrack.h" |
722347d8 | 23 | #include "AliLog.h" |
24 | #include "AliMCParticle.h" | |
809a4336 | 25 | #include "AliPID.h" |
26 | ||
27 | #include "AliHFEpidTOF.h" | |
28 | #include "AliHFEpidBase.h" | |
29 | ||
30 | ||
31 | ClassImp(AliHFEpidTOF) | |
32 | ||
33 | //___________________________________________________________________ | |
34 | AliHFEpidTOF::AliHFEpidTOF(const Char_t *name): | |
35 | AliHFEpidBase(name) | |
36 | , fPID(0x0) | |
37 | , fQAList(0x0) | |
38 | { | |
39 | // | |
40 | // Constructor | |
41 | // | |
42 | } | |
43 | //___________________________________________________________________ | |
44 | AliHFEpidTOF::AliHFEpidTOF(const AliHFEpidTOF &c): | |
45 | AliHFEpidBase("") | |
46 | , fPID(0x0) | |
47 | , fQAList(0x0) | |
48 | { | |
49 | // | |
50 | // Copy operator | |
51 | // | |
52 | ||
53 | c.Copy(*this); | |
54 | } | |
55 | //___________________________________________________________________ | |
56 | AliHFEpidTOF &AliHFEpidTOF::operator=(const AliHFEpidTOF &ref){ | |
57 | // | |
58 | // Assignment operator | |
59 | // | |
60 | ||
61 | if(this != &ref){ | |
62 | ref.Copy(*this); | |
63 | } | |
64 | ||
65 | return *this; | |
66 | } | |
67 | //___________________________________________________________________ | |
68 | AliHFEpidTOF::~AliHFEpidTOF(){ | |
69 | // | |
70 | // Destructor | |
71 | // | |
72 | if(fPID) delete fPID; | |
73 | if(fQAList){ | |
74 | fQAList->Delete(); | |
75 | delete fQAList; | |
76 | } | |
77 | } | |
78 | //___________________________________________________________________ | |
79 | void AliHFEpidTOF::Copy(TObject &ref) const { | |
80 | // | |
81 | // Performs the copying of the object | |
82 | // | |
83 | AliHFEpidTOF &target = dynamic_cast<AliHFEpidTOF &>(ref); | |
84 | ||
85 | target.fPID = fPID; | |
86 | target.fQAList = fQAList; | |
87 | ||
88 | AliHFEpidBase::Copy(ref); | |
89 | } | |
90 | //___________________________________________________________________ | |
91 | Bool_t AliHFEpidTOF::InitializePID(){ | |
92 | // | |
93 | // InitializePID: TOF experts have to implement code here | |
94 | // | |
95 | return kTRUE; | |
96 | } | |
97 | ||
98 | //___________________________________________________________________ | |
722347d8 | 99 | Int_t AliHFEpidTOF::IsSelected(AliHFEpidObject *vtrack) |
809a4336 | 100 | { |
101 | ||
102 | // | |
103 | // as of 22/05/2006 : | |
104 | // returns AliPID based on the ESD TOF PID decision | |
105 | // the ESD PID will be checked and if necessary improved | |
106 | // in the mean time. Best of luck | |
107 | // | |
108 | // returns 10 (== kUnknown)if PID can not be assigned | |
109 | // | |
722347d8 | 110 | if(vtrack->fAnalysisType == AliHFEpidObject::kESDanalysis){ |
111 | AliESDtrack *esdTrack = dynamic_cast<AliESDtrack *>(vtrack->fRecTrack); | |
112 | if(!esdTrack) return 0; | |
113 | AliMCParticle *mcTrack = dynamic_cast<AliMCParticle *>(vtrack->fMCtrack); | |
114 | return MakePIDesd(esdTrack, mcTrack); | |
115 | } else { | |
116 | AliAODTrack *aodTrack = dynamic_cast<AliAODTrack *>(vtrack->fRecTrack); | |
117 | if(!aodTrack) return 0; | |
118 | AliAODMCParticle *aodmc = dynamic_cast<AliAODMCParticle *>(vtrack->fMCtrack); | |
119 | return MakePIDaod(aodTrack, aodmc); | |
120 | } | |
121 | } | |
809a4336 | 122 | |
722347d8 | 123 | //___________________________________________________________________ |
124 | Int_t AliHFEpidTOF::MakePIDesd(AliESDtrack *track, AliMCParticle * /*mcTrack*/){ | |
ad75027f | 125 | Long_t status = 0; |
126 | status = track->GetStatus(); | |
127 | ||
128 | if(!(status & AliESDtrack::kTOFout)) return AliPID::kUnknown; | |
809a4336 | 129 | |
130 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(0.); | |
131 | ||
ad75027f | 132 | Double_t tItrackL = track->GetIntegratedLength(); |
133 | Double_t tTOFsignal = track->GetTOFsignal(); | |
134 | ||
135 | if(tItrackL > 0) | |
809a4336 | 136 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(1.); |
137 | ||
ad75027f | 138 | if(tTOFsignal > 0) |
809a4336 | 139 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(2.); |
140 | ||
141 | ||
ad75027f | 142 | if(tItrackL <=0 || tTOFsignal <=0) return AliPID::kUnknown; |
809a4336 | 143 | |
144 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(3.); | |
ad75027f | 145 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFsignal)))->Fill(tTOFsignal/1000.); |
146 | (dynamic_cast<TH1F *>(fQAList->At(kHistTOFlength)))->Fill(tItrackL); | |
809a4336 | 147 | // get the TOF pid probabilities |
ad75027f | 148 | Double_t tESDpid[5] = {0., 0., 0., 0., 0.}; |
75d81601 | 149 | Float_t tTOFpidSum = 0.; |
809a4336 | 150 | // find the largest PID probability |
ad75027f | 151 | track->GetTOFpid(tESDpid); |
152 | Double_t tMAXpid = 0.; | |
153 | Int_t tMAXindex = -1; | |
809a4336 | 154 | for(Int_t i=0; i<5; ++i){ |
75d81601 | 155 | tTOFpidSum += tESDpid[i]; |
ad75027f | 156 | if(tESDpid[i] > tMAXpid){ |
157 | tMAXpid = tESDpid[i]; | |
158 | tMAXindex = i; | |
809a4336 | 159 | } |
160 | } | |
161 | ||
ad75027f | 162 | Double_t p = track->GetOuterParam()->P(); |
163 | Double_t beta = (tItrackL/100.)/(TMath::C()*(tTOFsignal/1e12)); | |
809a4336 | 164 | |
75d81601 | 165 | if(TMath::Abs(tTOFpidSum - 1) > 0.01) return AliPID::kUnknown; |
809a4336 | 166 | else{ |
167 | // should be the same as AliPID flags | |
168 | ||
75d81601 | 169 | (dynamic_cast<TH2F *>(fQAList->At(kHistTOFpid0+tMAXindex)))->Fill(beta, p); |
170 | (dynamic_cast<TH2F *>(fQAList->At(kHistTOFpidBetavP)))->Fill(beta, p); | |
ad75027f | 171 | return tMAXindex; |
809a4336 | 172 | } |
173 | } | |
174 | ||
722347d8 | 175 | //___________________________________________________________________ |
176 | Int_t AliHFEpidTOF::MakePIDaod(AliAODTrack * /*aodTrack*/, AliAODMCParticle * /*mcTrack*/){ | |
177 | AliError("AOD PID not yet implemented"); | |
178 | return 0; | |
179 | } | |
180 | ||
809a4336 | 181 | //___________________________________________________________________ |
182 | void AliHFEpidTOF::AddQAhistograms(TList *qaList){ | |
183 | // | |
184 | // Create QA histograms for TOF PID | |
185 | // | |
186 | ||
187 | fQAList = new TList; | |
188 | fQAList->SetName("fTOFqaHistos"); | |
189 | fQAList->AddAt(new TH1F("hTOF_flags", "TOF flags;flags (see code for info);counts", 10, -0.25, 4.75), kHistTOFpidFlags); | |
75d81601 | 190 | fQAList->AddAt(new TH2F("fTOFbeta_v_P_no","beta -v- P; beta;momentum [GeV/c]", 120, 0, 1.2, 200, 0, 20), kHistTOFpidBetavP); |
809a4336 | 191 | fQAList->AddAt(new TH1F("hTOF_signal", "TOF signal; TOF signal [ns];counts", 1000, 12, 50), kHistTOFsignal); |
192 | fQAList->AddAt(new TH1F("hTOF_length", "TOF track length; length [cm];counts", 400, 300, 700), kHistTOFlength); | |
75d81601 | 193 | fQAList->AddAt(new TH2F("hTOFpid_electron", "TOF reco electron; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid0); |
194 | fQAList->AddAt(new TH2F("hTOFpid_muon", "TOF reco muon; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid1); | |
195 | fQAList->AddAt(new TH2F("hTOFpid_pion", "TOF reco pion; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid2); | |
196 | fQAList->AddAt(new TH2F("hTOFpid_kaon", "TOF reco kaon; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid3); | |
197 | fQAList->AddAt(new TH2F("hTOFpid_proton", "TOF reco proton; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid4); | |
809a4336 | 198 | |
199 | qaList->AddLast(fQAList); | |
200 | } |