]>
Commit | Line | Data |
---|---|---|
00514d01 | 1 | // $Id$ |
91f4b7c5 | 2 | // |
980821ba | 3 | // General QA task (S.Aiola). |
91f4b7c5 | 4 | // |
5 | // | |
00514d01 | 6 | |
7 | #include <TChain.h> | |
8 | #include <TClonesArray.h> | |
9 | #include <TH1F.h> | |
10 | #include <TH2F.h> | |
11 | #include <TList.h> | |
12 | #include <TLorentzVector.h> | |
13 | #include <TParticle.h> | |
14 | ||
15 | #include "AliAnalysisManager.h" | |
16 | #include "AliCentrality.h" | |
17 | #include "AliVCluster.h" | |
18 | #include "AliESDtrack.h" | |
19 | #include "AliEmcalJet.h" | |
20 | #include "AliAODTrack.h" | |
21 | #include "AliESDtrack.h" | |
00514d01 | 22 | #include "AliVEventHandler.h" |
23 | #include "AliPicoTrack.h" | |
24 | ||
25 | #include "AliAnalysisTaskSAQA.h" | |
26 | ||
27 | ClassImp(AliAnalysisTaskSAQA) | |
28 | ||
29 | //________________________________________________________________________ | |
30 | AliAnalysisTaskSAQA::AliAnalysisTaskSAQA() : | |
31 | AliAnalysisTaskSE("AliAnalysisTaskSAQA"), | |
32 | fOutput(0), | |
33 | fTracksName("Tracks"), | |
34 | fCaloName("CaloClusters"), | |
00514d01 | 35 | fTrgClusName("ClustersL1GAMMAFEE"), |
36 | fTracks(0), | |
37 | fCaloClusters(0), | |
38 | fJets(0), | |
39 | fTrgClusters(0), | |
40 | fCent(0), | |
41 | fHistCentrality(0), | |
42 | fHistTracksCent(0), | |
43 | fHistClusCent(0), | |
44 | fHistTracksPt(0), | |
45 | fHistClustersEnergy(0), | |
46 | fHistEPcorrelation(0), | |
00514d01 | 47 | fHistTrPhiEta(0), |
48 | fHistClusPhiEta(0), | |
00514d01 | 49 | fHistMaxTrgCluster(0), |
c554a987 | 50 | fHistClusPhiCorr(0), |
51 | fHistTracksPhiCorr(0), | |
00514d01 | 52 | Ptbins(100), |
53 | Ptlow(0), | |
54 | Ptup(50), | |
55 | Ebins(100), | |
56 | Elow(0), | |
57 | Eup(50) | |
58 | { | |
59 | // Default constructor. | |
60 | ||
61 | for (Int_t i = 0; i < 5; i++) { | |
62 | fHistTrackPhi[i] = 0; | |
63 | fHistTrackEta[i] = 0; | |
64 | } | |
65 | ||
91f4b7c5 | 66 | // Output slot #1 writes into a TH1 container |
67 | DefineOutput(1, TList::Class()); | |
00514d01 | 68 | } |
69 | ||
70 | //________________________________________________________________________ | |
71 | AliAnalysisTaskSAQA::AliAnalysisTaskSAQA(const char *name) : | |
72 | AliAnalysisTaskSE(name), | |
73 | fOutput(0), | |
74 | fTracksName("Tracks"), | |
75 | fCaloName("CaloClusters"), | |
00514d01 | 76 | fTrgClusName("ClustersL1GAMMAFEE"), |
77 | fTracks(0), | |
78 | fCaloClusters(0), | |
79 | fJets(0), | |
80 | fTrgClusters(0), | |
81 | fCent(0), | |
82 | fHistCentrality(0), | |
83 | fHistTracksCent(0), | |
84 | fHistClusCent(0), | |
85 | fHistTracksPt(0), | |
86 | fHistClustersEnergy(0), | |
87 | fHistEPcorrelation(0), | |
00514d01 | 88 | fHistTrPhiEta(0), |
89 | fHistClusPhiEta(0), | |
00514d01 | 90 | fHistMaxTrgCluster(0), |
c554a987 | 91 | fHistClusPhiCorr(0), |
92 | fHistTracksPhiCorr(0), | |
00514d01 | 93 | Ptbins(100), |
94 | Ptlow(0), | |
95 | Ptup(50), | |
96 | Ebins(100), | |
97 | Elow(0), | |
98 | Eup(50) | |
99 | { | |
100 | // Standard constructor. | |
101 | ||
102 | for (Int_t i = 0; i < 5; i++) { | |
103 | fHistTrackPhi[i] = 0; | |
104 | fHistTrackEta[i] = 0; | |
105 | } | |
106 | ||
91f4b7c5 | 107 | // Output slot #1 writes into a TH1 container |
108 | DefineOutput(1, TList::Class()); | |
00514d01 | 109 | } |
110 | ||
111 | //________________________________________________________________________ | |
112 | AliAnalysisTaskSAQA::~AliAnalysisTaskSAQA() | |
113 | { | |
114 | // Destructor | |
115 | } | |
116 | ||
117 | //________________________________________________________________________ | |
118 | void AliAnalysisTaskSAQA::UserCreateOutputObjects() | |
119 | { | |
120 | // Create histograms | |
121 | ||
122 | fOutput = new TList(); | |
123 | fOutput->SetOwner(); // IMPORTANT! | |
124 | ||
125 | fHistCentrality = new TH1F("fHistCentrality","Event centrality distribution", Ebins, 0, 100); | |
126 | fHistCentrality->GetXaxis()->SetTitle("Centrality (%)"); | |
127 | fHistCentrality->GetYaxis()->SetTitle("counts"); | |
128 | fOutput->Add(fHistCentrality); | |
129 | ||
130 | fHistTracksCent = new TH2F("fHistTracksCent","Tracks vs. centrality", Ebins, 0, 100, Ebins, 0, 4000); | |
131 | fHistTracksCent->GetXaxis()->SetTitle("Centrality (%)"); | |
132 | fHistTracksCent->GetYaxis()->SetTitle("No. of tracks"); | |
133 | fOutput->Add(fHistTracksCent); | |
134 | ||
135 | fHistClusCent = new TH2F("fHistClusCent","Clusters vs. centrality", Ebins, 0, 100, Ebins, 0, 2000); | |
136 | fHistClusCent->GetXaxis()->SetTitle("Centrality (%)"); | |
137 | fHistClusCent->GetYaxis()->SetTitle("No. of clusters"); | |
138 | fOutput->Add(fHistClusCent); | |
139 | ||
140 | fHistTracksPt = new TH1F("fHistTracksPt","P_{T} spectrum of reconstructed tracks", Ptbins, Ptlow, Ptup); | |
141 | fHistTracksPt->GetXaxis()->SetTitle("P_{T} [GeV/c]"); | |
142 | fHistTracksPt->GetYaxis()->SetTitle("counts"); | |
143 | fOutput->Add(fHistTracksPt); | |
144 | ||
145 | fHistClustersEnergy = new TH1F("fHistClustersEnergy","Energy spectrum of clusters", Ebins, Elow, Eup); | |
146 | fHistClustersEnergy->GetXaxis()->SetTitle("E [GeV]"); | |
147 | fHistClustersEnergy->GetYaxis()->SetTitle("counts"); | |
148 | fOutput->Add(fHistClustersEnergy); | |
149 | ||
150 | fHistEPcorrelation = new TH2F("fHistEPcorrelation","Energy-momentum correlation", Ptbins, Ptlow, Ptup, Ebins, Elow, Eup); | |
151 | fHistEPcorrelation->GetXaxis()->SetTitle("P [GeV/c]"); | |
152 | fHistEPcorrelation->GetYaxis()->SetTitle("E [GeV]"); | |
153 | fOutput->Add(fHistEPcorrelation); | |
00514d01 | 154 | |
155 | fHistTrPhiEta = new TH2F("fHistTrPhiEta","Phi-Eta distribution of tracks", 20, -2, 2, 32, 0, 6.4); | |
156 | fHistTrPhiEta->GetXaxis()->SetTitle("Eta"); | |
157 | fHistTrPhiEta->GetYaxis()->SetTitle("Phi"); | |
158 | fOutput->Add(fHistTrPhiEta); | |
159 | ||
160 | fHistClusPhiEta = new TH2F("fHistClusPhiEta","Phi-Eta distribution of clusters", 20, -2, 2, 32, 0, 6.4); | |
161 | fHistClusPhiEta->GetXaxis()->SetTitle("Eta"); | |
162 | fHistClusPhiEta->GetYaxis()->SetTitle("Phi"); | |
163 | fOutput->Add(fHistClusPhiEta); | |
164 | ||
00514d01 | 165 | fHistMaxTrgCluster = new TH1F("fHistMaxTrgCluster","Energy distribution of max trigger clusters", Ebins, Elow, Eup); |
166 | fHistMaxTrgCluster->GetXaxis()->SetTitle("E [GeV]"); | |
167 | fHistMaxTrgCluster->GetYaxis()->SetTitle("counts"); | |
168 | fOutput->Add(fHistMaxTrgCluster); | |
c554a987 | 169 | |
170 | fHistTracksPhiCorr = new TH1F("fHistTracksPhiCorr", "fHistTracksPhiCorr", 128, -1.6, 4.8); | |
171 | fHistTracksPhiCorr->GetXaxis()->SetTitle("#Delta#phi"); | |
172 | fHistTracksPhiCorr->GetYaxis()->SetTitle("counts"); | |
173 | fOutput->Add(fHistTracksPhiCorr); | |
174 | ||
175 | fHistClusPhiCorr = new TH1F("fHistClusPhiCorr", "fHistClusPhiCorr", 128, -1.6, 4.8); | |
176 | fHistClusPhiCorr->GetXaxis()->SetTitle("#Delta#phi"); | |
177 | fHistClusPhiCorr->GetYaxis()->SetTitle("counts"); | |
178 | fOutput->Add(fHistClusPhiCorr); | |
00514d01 | 179 | |
180 | for (Int_t i = 0; i < 5; i++) { | |
181 | TString histnamephi("fHistTrackPhi_"); | |
182 | histnamephi += i; | |
183 | fHistTrackPhi[i] = new TH1F(histnamephi.Data(),histnamephi.Data(), 128, 0, 6.4); | |
184 | fHistTrackPhi[i]->GetXaxis()->SetTitle("Phi"); | |
185 | fOutput->Add(fHistTrackPhi[i]); | |
186 | ||
187 | TString histnameeta("fHistTrackEta_"); | |
188 | histnameeta += i; | |
189 | fHistTrackEta[i] = new TH1F(histnameeta.Data(),histnameeta.Data(), 100, -2, 2); | |
190 | fHistTrackEta[i]->GetXaxis()->SetTitle("Eta"); | |
191 | fOutput->Add(fHistTrackEta[i]); | |
192 | } | |
193 | ||
194 | fHistTrackPhi[0]->SetLineColor(kRed); | |
195 | fHistTrackEta[0]->SetLineColor(kRed); | |
196 | fHistTrackPhi[1]->SetLineColor(kBlue); | |
197 | fHistTrackEta[1]->SetLineColor(kBlue); | |
198 | fHistTrackPhi[2]->SetLineColor(kGreen); | |
199 | fHistTrackEta[2]->SetLineColor(kGreen); | |
200 | fHistTrackPhi[3]->SetLineColor(kOrange); | |
201 | fHistTrackEta[3]->SetLineColor(kOrange); | |
202 | fHistTrackPhi[4]->SetLineColor(kBlack); | |
203 | fHistTrackEta[4]->SetLineColor(kBlack); | |
204 | ||
205 | PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram | |
206 | } | |
207 | ||
208 | void AliAnalysisTaskSAQA::RetrieveEventObjects() | |
209 | { | |
210 | fCaloClusters = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fCaloName)); | |
211 | if (!fCaloClusters) { | |
212 | AliWarning(Form("Could not retrieve clusters!")); | |
213 | } | |
214 | ||
215 | fTracks = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTracksName)); | |
216 | if (!fTracks) { | |
217 | AliWarning(Form("Could not retrieve tracks!")); | |
218 | } | |
219 | ||
00514d01 | 220 | if (strcmp(fTrgClusName,"")) { |
221 | fTrgClusters = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTrgClusName)); | |
222 | if (!fTrgClusters) { | |
223 | AliWarning(Form("Could not retrieve trigger clusters!")); | |
224 | } | |
225 | } | |
226 | ||
227 | fCent = InputEvent()->GetCentrality(); | |
228 | } | |
229 | ||
230 | AliVTrack* AliAnalysisTaskSAQA::GetTrack(const Int_t i) const | |
231 | { | |
232 | if (fTracks) | |
233 | return dynamic_cast<AliVTrack*>(fTracks->At(i)); | |
234 | else | |
235 | return 0; | |
236 | } | |
237 | ||
238 | Int_t AliAnalysisTaskSAQA::GetNumberOfTracks() const | |
239 | { | |
240 | if (fTracks) | |
241 | return fTracks->GetEntriesFast(); | |
242 | else | |
243 | return 0; | |
244 | } | |
245 | ||
246 | AliVCluster* AliAnalysisTaskSAQA::GetCaloCluster(const Int_t i) const | |
247 | { | |
248 | if (fCaloClusters) | |
249 | return dynamic_cast<AliVCluster*>(fCaloClusters->At(i)); | |
250 | else | |
251 | return 0; | |
252 | } | |
253 | ||
254 | Int_t AliAnalysisTaskSAQA::GetNumberOfCaloClusters() const | |
255 | { | |
256 | if (fCaloClusters) | |
257 | return fCaloClusters->GetEntriesFast(); | |
258 | else | |
259 | return 0; | |
260 | } | |
261 | ||
00514d01 | 262 | AliVCluster* AliAnalysisTaskSAQA::GetTrgCluster(const Int_t i) const |
263 | { | |
264 | if (fTrgClusters) | |
265 | return dynamic_cast<AliVCluster*>(fTrgClusters->At(i)); | |
266 | else | |
267 | return 0; | |
268 | } | |
269 | ||
270 | Int_t AliAnalysisTaskSAQA::GetNumberOfTrgClusters() const | |
271 | { | |
272 | if (fTrgClusters) | |
273 | return fTrgClusters->GetEntriesFast(); | |
274 | else | |
275 | return 0; | |
276 | } | |
277 | ||
278 | void AliAnalysisTaskSAQA::FillHistograms() | |
279 | { | |
91f4b7c5 | 280 | Float_t cent = 100; |
00514d01 | 281 | |
91f4b7c5 | 282 | if (fCent) |
283 | cent = fCent->GetCentralityPercentile("V0M"); | |
284 | else | |
285 | AliWarning("Centrality not available!"); | |
286 | ||
287 | fHistCentrality->Fill(cent); | |
288 | fHistTracksCent->Fill(cent, GetNumberOfTracks()); | |
289 | fHistClusCent->Fill(cent, GetNumberOfCaloClusters()); | |
00514d01 | 290 | |
291 | // Cluster loop | |
292 | Int_t nclusters = GetNumberOfCaloClusters(); | |
293 | //cout << nclusters << " clusters" << endl; | |
294 | for (Int_t iClusters = 0; iClusters < nclusters; iClusters++) { | |
295 | AliVCluster* cluster = GetCaloCluster(iClusters); | |
296 | if (!cluster) { | |
297 | printf("ERROR: Could not receive cluster %d\n", iClusters); | |
298 | continue; | |
299 | } | |
300 | ||
301 | if (!(cluster->IsEMCAL())) continue; | |
302 | ||
303 | fHistClustersEnergy->Fill(cluster->E()); | |
304 | ||
305 | Float_t pos[3]; | |
306 | cluster->GetPosition(pos); | |
307 | TVector3 clusVec(pos); | |
308 | fHistClusPhiEta->Fill(clusVec.Eta(), clusVec.Phi()); | |
c554a987 | 309 | |
310 | for(Int_t ic2 = iClusters+1; ic2 < nclusters; ic2++) { | |
311 | AliVCluster* cluster2 = GetCaloCluster(ic2); | |
312 | if (!cluster2) { | |
313 | printf("ERROR: Could not receive cluster %d\n", ic2); | |
314 | continue; | |
315 | } | |
316 | ||
317 | if (!(cluster2->IsEMCAL())) continue; | |
318 | ||
319 | Float_t pos2[3]; | |
320 | cluster2->GetPosition(pos2); | |
321 | TVector3 clusVec2(pos2); | |
322 | ||
323 | Float_t dphi = clusVec.Phi() - clusVec2.Phi(); | |
324 | if (dphi < -1.6) dphi += TMath::Pi() * 2; | |
325 | if (dphi > 4.8) dphi -= TMath::Pi() * 2; | |
326 | fHistClusPhiCorr->Fill(dphi); | |
327 | } | |
00514d01 | 328 | } //cluster loop |
329 | ||
00514d01 | 330 | // Track loop |
331 | Int_t ntracks = GetNumberOfTracks(); | |
332 | //cout << ntracks << " tracks" << endl; | |
333 | for(Int_t i = 0; i < ntracks; i++) { | |
334 | //cout << "track n. " << i << endl; | |
335 | AliVTrack* track = GetTrack(i); // pointer to reconstructed to track | |
336 | if(!track) { | |
337 | AliError(Form("ERROR: Could not retrieve esdtrack %d",i)); | |
338 | continue; | |
339 | } | |
340 | ||
341 | if (!AcceptTrack(track)) continue; | |
342 | ||
c554a987 | 343 | for(Int_t it2 = i+1; it2 < ntracks; it2++) { |
344 | AliVTrack* track2 = GetTrack(it2); | |
345 | if(!track2) { | |
346 | AliError(Form("ERROR: Could not retrieve track %d", it2)); | |
347 | continue; | |
348 | } | |
349 | ||
350 | if (!AcceptTrack(track2)) continue; | |
351 | ||
352 | Float_t dphi = track->Phi() - track2->Phi(); | |
353 | if (dphi < -1.6) dphi += TMath::Pi() * 2; | |
354 | if (dphi > 4.8) dphi -= TMath::Pi() * 2; | |
355 | fHistTracksPhiCorr->Fill(dphi); | |
356 | } // second track loop | |
357 | ||
00514d01 | 358 | fHistTracksPt->Fill(track->Pt()); |
359 | ||
360 | Int_t clId = track->GetEMCALcluster(); | |
c554a987 | 361 | if (clId > -1 && clId < nclusters) { |
00514d01 | 362 | AliVCluster* cluster = GetCaloCluster(clId); |
363 | if (cluster) | |
364 | fHistEPcorrelation->Fill(track->P(),cluster->E()); | |
365 | } | |
366 | ||
c554a987 | 367 | Int_t label = track->GetLabel(); |
00514d01 | 368 | |
c554a987 | 369 | fHistTrPhiEta->Fill(track->Eta(), track->Phi()); |
00514d01 | 370 | |
c554a987 | 371 | fHistTrackEta[4]->Fill(track->Eta()); |
372 | fHistTrackPhi[4]->Fill(track->Phi()); | |
00514d01 | 373 | |
374 | if (label >= 0 && label < 4) { | |
c554a987 | 375 | fHistTrackEta[label]->Fill(track->Eta()); |
376 | fHistTrackPhi[label]->Fill(track->Phi()); | |
00514d01 | 377 | } |
378 | else { | |
379 | AliWarning(Form("Track label %d not recognized!",label)); | |
380 | } | |
00514d01 | 381 | } |
382 | ||
00514d01 | 383 | Int_t ntrgclusters = GetNumberOfTrgClusters(); |
384 | Float_t maxe = 0; | |
385 | //cout << ntrgclusters << " clusters" << endl; | |
386 | for (Int_t iClusters = 0; iClusters < ntrgclusters; iClusters++) { | |
387 | AliVCluster* cluster = GetTrgCluster(iClusters); | |
388 | if (!cluster) { | |
389 | printf("ERROR: Could not receive cluster %d\n", iClusters); | |
390 | continue; | |
391 | } | |
392 | ||
393 | if (!(cluster->IsEMCAL())) continue; | |
394 | ||
395 | if (cluster->E() > maxe) | |
396 | maxe = cluster->E(); | |
397 | ||
398 | } //cluster loop | |
399 | fHistMaxTrgCluster->Fill(maxe); | |
400 | ||
401 | } | |
402 | ||
403 | //________________________________________________________________________ | |
404 | Bool_t AliAnalysisTaskSAQA::AcceptTrack(AliVTrack* /*track*/) | |
405 | { | |
406 | return kTRUE; | |
407 | } | |
408 | ||
409 | //________________________________________________________________________ | |
410 | void AliAnalysisTaskSAQA::UserExec(Option_t *) | |
411 | { | |
412 | // Main loop, called for each event. | |
413 | // Add jets to event if not yet there | |
414 | ||
415 | RetrieveEventObjects(); | |
416 | ||
417 | FillHistograms(); | |
418 | ||
419 | // information for this iteration of the UserExec in the container | |
420 | PostData(1, fOutput); | |
421 | } | |
422 | ||
423 | //________________________________________________________________________ | |
424 | void AliAnalysisTaskSAQA::Terminate(Option_t *) | |
425 | { | |
426 | // Called once at the end of the analysis. | |
427 | } |