]>
Commit | Line | Data |
---|---|---|
0f67a57c | 1 | /* $Id$ */ |
2 | ||
3 | #include "AlidNdEtaCorrectionTask.h" | |
4 | ||
5 | #include <TCanvas.h> | |
6 | #include <TChain.h> | |
7 | #include <TFile.h> | |
8 | #include <TH1F.h> | |
9 | #include <TParticle.h> | |
10 | ||
11 | #include <AliLog.h> | |
12 | #include <AliESDVertex.h> | |
13 | #include <AliESDEvent.h> | |
14 | #include <AliStack.h> | |
15 | #include <AliHeader.h> | |
16 | #include <AliGenEventHeader.h> | |
17 | #include <AliMultiplicity.h> | |
18 | #include <AliAnalysisManager.h> | |
19 | #include <AliMCEventHandler.h> | |
20 | #include <AliMCEvent.h> | |
21 | #include <AliESDInputHandler.h> | |
22 | ||
23 | #include "esdTrackCuts/AliESDtrackCuts.h" | |
24 | #include "AliPWG0Helper.h" | |
25 | //#include "AliCorrection.h" | |
26 | //#include "AliCorrectionMatrix3D.h" | |
27 | #include "dNdEta/dNdEtaAnalysis.h" | |
28 | #include "dNdEta/AlidNdEtaCorrection.h" | |
29 | ||
30 | ClassImp(AlidNdEtaCorrectionTask) | |
31 | ||
32 | AlidNdEtaCorrectionTask::AlidNdEtaCorrectionTask(const char* opt) : | |
33 | AliAnalysisTask("AlidNdEtaCorrectionTask", ""), | |
34 | fESD(0), | |
35 | fOutput(0), | |
36 | fOption(opt), | |
37 | fAnalysisMode(kTPC), | |
38 | fSignMode(0), | |
39 | fEsdTrackCuts(0), | |
40 | fdNdEtaCorrection(0), | |
41 | fdNdEtaAnalysisMC(0), | |
42 | fdNdEtaAnalysisESD(0), | |
43 | fPIDParticles(0), | |
44 | fPIDTracks(0) | |
45 | { | |
46 | // | |
47 | // Constructor. Initialization of pointers | |
48 | // | |
49 | ||
50 | // Define input and output slots here | |
51 | DefineInput(0, TChain::Class()); | |
52 | DefineOutput(0, TList::Class()); | |
53 | } | |
54 | ||
55 | AlidNdEtaCorrectionTask::~AlidNdEtaCorrectionTask() | |
56 | { | |
57 | // | |
58 | // Destructor | |
59 | // | |
60 | ||
61 | // histograms are in the output list and deleted when the output | |
62 | // list is deleted by the TSelector dtor | |
63 | ||
64 | if (fOutput) { | |
65 | delete fOutput; | |
66 | fOutput = 0; | |
67 | } | |
68 | } | |
69 | ||
70 | //________________________________________________________________________ | |
71 | void AlidNdEtaCorrectionTask::ConnectInputData(Option_t *) | |
72 | { | |
73 | // Connect ESD | |
74 | // Called once | |
75 | ||
76 | Printf("AlidNdEtaCorrectionTask::ConnectInputData called"); | |
77 | ||
78 | TTree* tree = dynamic_cast<TTree*> (GetInputData(0)); | |
79 | if (!tree) { | |
80 | Printf("ERROR: Could not read tree from input slot 0"); | |
81 | } else { | |
82 | // Disable all branches and enable only the needed ones | |
83 | tree->SetBranchStatus("*", 0); | |
84 | ||
85 | tree->SetBranchStatus("fTriggerMask", 1); | |
86 | tree->SetBranchStatus("fSPDVertex*", 1); | |
87 | ||
88 | if (fAnalysisMode == kSPD) | |
89 | tree->SetBranchStatus("fSPDMult*", 1); | |
90 | ||
91 | if (fAnalysisMode == kTPC) { | |
92 | AliESDtrackCuts::EnableNeededBranches(tree); | |
93 | tree->SetBranchStatus("fTracks.fLabel", 1); | |
94 | } | |
95 | ||
96 | tree->SetCacheSize(0); | |
97 | ||
98 | AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
99 | ||
100 | if (!esdH) { | |
101 | Printf("ERROR: Could not get ESDInputHandler"); | |
102 | } else | |
103 | fESD = esdH->GetEvent(); | |
104 | } | |
105 | } | |
106 | ||
107 | void AlidNdEtaCorrectionTask::CreateOutputObjects() | |
108 | { | |
109 | // create result objects and add to output list | |
110 | ||
111 | if (fOption.Contains("only-positive")) | |
112 | { | |
113 | Printf("INFO: Processing only positive particles."); | |
114 | fSignMode = 1; | |
115 | } | |
116 | else if (fOption.Contains("only-negative")) | |
117 | { | |
118 | Printf("INFO: Processing only negative particles."); | |
119 | fSignMode = -1; | |
120 | } | |
121 | ||
122 | TString detector; | |
123 | if (fAnalysisMode == kTPC) | |
124 | { | |
125 | detector = "TPC"; | |
126 | } | |
127 | else if (fAnalysisMode == kSPD) | |
128 | detector = "SPD"; | |
129 | ||
130 | fOutput = new TList; | |
131 | fOutput->SetOwner(); | |
132 | ||
133 | fdNdEtaCorrection = new AlidNdEtaCorrection("dndeta_correction", "dndeta_correction", detector); | |
134 | fOutput->Add(fdNdEtaCorrection); | |
135 | ||
136 | fPIDParticles = new TH1F("pid_particles", "PID of generated primary particles", 10001, -5000.5, 5000.5); | |
137 | fOutput->Add(fPIDParticles); | |
138 | ||
139 | fPIDTracks = new TH1F("pid_tracks", "MC PID of reconstructed tracks", 10001, -5000.5, 5000.5); | |
140 | fOutput->Add(fPIDTracks); | |
141 | ||
142 | fdNdEtaAnalysisMC = new dNdEtaAnalysis("dndetaMC", "dndetaMC", detector); | |
143 | fOutput->Add(fdNdEtaAnalysisMC); | |
144 | ||
145 | fdNdEtaAnalysisESD = new dNdEtaAnalysis("dndetaESD", "dndetaESD", detector); | |
146 | fOutput->Add(fdNdEtaAnalysisESD); | |
147 | } | |
148 | ||
149 | void AlidNdEtaCorrectionTask::Exec(Option_t*) | |
150 | { | |
151 | // process the event | |
152 | ||
153 | // Check prerequisites | |
154 | if (!fESD) | |
155 | { | |
156 | AliDebug(AliLog::kError, "ESD branch not available"); | |
157 | return; | |
158 | } | |
159 | ||
160 | // trigger definition | |
161 | Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD->GetTriggerMask(), AliPWG0Helper::kMB2); | |
162 | ||
163 | Bool_t eventVertex = AliPWG0Helper::IsVertexReconstructed(fESD->GetVertex()); | |
164 | ||
165 | // post the data already here | |
166 | PostData(0, fOutput); | |
167 | ||
168 | // create list of (label, eta, pt) tuples | |
169 | Int_t inputCount = 0; | |
170 | Int_t* labelArr = 0; | |
171 | Float_t* etaArr = 0; | |
172 | Float_t* ptArr = 0; | |
173 | if (fAnalysisMode == kSPD) | |
174 | { | |
175 | // get tracklets | |
176 | const AliMultiplicity* mult = fESD->GetMultiplicity(); | |
177 | if (!mult) | |
178 | { | |
179 | AliDebug(AliLog::kError, "AliMultiplicity not available"); | |
180 | return; | |
181 | } | |
182 | ||
183 | labelArr = new Int_t[mult->GetNumberOfTracklets()]; | |
184 | etaArr = new Float_t[mult->GetNumberOfTracklets()]; | |
185 | ptArr = new Float_t[mult->GetNumberOfTracklets()]; | |
186 | ||
187 | // get multiplicity from ITS tracklets | |
188 | for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i) | |
189 | { | |
190 | //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i)); | |
191 | ||
192 | // This removes non-tracklets in PDC06 data. Very bad solution. New solution is implemented for newer data. Keeping this for compatibility. | |
193 | if (mult->GetDeltaPhi(i) < -1000) | |
194 | continue; | |
195 | ||
196 | etaArr[inputCount] = mult->GetEta(i); | |
197 | labelArr[inputCount] = mult->GetLabel(i); | |
198 | ptArr[inputCount] = 0; // no pt for tracklets | |
199 | ++inputCount; | |
200 | } | |
201 | } | |
202 | else if (fAnalysisMode == kTPC) | |
203 | { | |
204 | if (!fEsdTrackCuts) | |
205 | { | |
206 | AliDebug(AliLog::kError, "fESDTrackCuts not available"); | |
207 | return; | |
208 | } | |
209 | ||
210 | // get multiplicity from ESD tracks | |
211 | TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD); | |
212 | Int_t nGoodTracks = list->GetEntries(); | |
213 | ||
214 | labelArr = new Int_t[nGoodTracks]; | |
215 | etaArr = new Float_t[nGoodTracks]; | |
216 | ptArr = new Float_t[nGoodTracks]; | |
217 | ||
218 | // loop over esd tracks | |
219 | for (Int_t i=0; i<nGoodTracks; i++) | |
220 | { | |
221 | AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(i)); | |
222 | if (!esdTrack) | |
223 | { | |
224 | AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", i)); | |
225 | continue; | |
226 | } | |
227 | ||
228 | etaArr[inputCount] = esdTrack->Eta(); | |
229 | labelArr[inputCount] = TMath::Abs(esdTrack->GetLabel()); | |
230 | ptArr[inputCount] = esdTrack->Pt(); | |
231 | ++inputCount; | |
232 | } | |
233 | } | |
234 | else | |
235 | return; | |
236 | ||
237 | AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()); | |
238 | if (!eventHandler) { | |
239 | Printf("ERROR: Could not retrieve MC event handler"); | |
240 | return; | |
241 | } | |
242 | ||
243 | AliMCEvent* mcEvent = eventHandler->MCEvent(); | |
244 | if (!mcEvent) { | |
245 | Printf("ERROR: Could not retrieve MC event"); | |
246 | return; | |
247 | } | |
248 | ||
249 | AliStack* stack = mcEvent->Stack(); | |
250 | if (!stack) | |
251 | { | |
252 | AliDebug(AliLog::kError, "Stack not available"); | |
253 | return; | |
254 | } | |
255 | ||
256 | AliHeader* header = mcEvent->Header(); | |
257 | if (!header) | |
258 | { | |
259 | AliDebug(AliLog::kError, "Header not available"); | |
260 | return; | |
261 | } | |
262 | ||
263 | // get process type | |
264 | Int_t processType = AliPWG0Helper::GetPythiaEventProcessType(header); | |
265 | AliDebug(AliLog::kDebug+1, Form("Found pythia procces type %d", processType)); | |
266 | ||
267 | if (processType<0) | |
268 | AliDebug(AliLog::kError, Form("Unknown Pythia process type %d.", processType)); | |
269 | ||
270 | // get the MC vertex | |
271 | AliGenEventHeader* genHeader = header->GenEventHeader(); | |
272 | TArrayF vtxMC(3); | |
273 | genHeader->PrimaryVertex(vtxMC); | |
274 | ||
275 | // loop over mc particles | |
276 | Int_t nPrim = stack->GetNprimary(); | |
277 | ||
278 | for (Int_t iMc = 0; iMc < nPrim; ++iMc) | |
279 | { | |
280 | //Printf("Getting particle %d", iMc); | |
281 | TParticle* particle = stack->Particle(iMc); | |
282 | ||
283 | if (!particle) | |
284 | continue; | |
285 | ||
286 | if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE) | |
287 | continue; | |
288 | ||
289 | if (SignOK(particle->GetPDG()) == kFALSE) | |
290 | continue; | |
291 | ||
292 | Float_t eta = particle->Eta(); | |
293 | Float_t pt = particle->Pt(); | |
294 | ||
295 | fdNdEtaCorrection->FillMCParticle(vtxMC[2], eta, pt, eventTriggered, eventVertex, processType); | |
296 | ||
297 | if (eventTriggered) | |
298 | if (eventVertex) | |
299 | fdNdEtaAnalysisMC->FillTrack(vtxMC[2], eta, pt); | |
300 | } | |
301 | ||
302 | for (Int_t i=0; i<inputCount; ++i) | |
303 | { | |
304 | Int_t label = labelArr[i]; | |
305 | ||
306 | if (label < 0) | |
307 | { | |
308 | AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", label)); | |
309 | continue; | |
310 | } | |
311 | ||
312 | TParticle* particle = stack->Particle(label); | |
313 | if (!particle) | |
314 | { | |
315 | AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label)); | |
316 | continue; | |
317 | } | |
318 | ||
319 | if (SignOK(particle->GetPDG()) == kFALSE) | |
320 | continue; | |
321 | ||
322 | if (eventTriggered && eventVertex) | |
323 | { | |
324 | fdNdEtaCorrection->FillTrackedParticle(vtxMC[2], particle->Eta(), particle->Pt()); | |
325 | fdNdEtaAnalysisESD->FillTrack(vtxMC[2], particle->Eta(), particle->Pt()); | |
326 | if (particle->Pt() > 0.1 && particle->Pt() < 0.2) | |
327 | { | |
328 | fPIDTracks->Fill(particle->GetPdgCode()); | |
329 | } | |
330 | } | |
331 | } | |
332 | ||
333 | if (eventTriggered && eventVertex) | |
334 | fdNdEtaAnalysisMC->FillEvent(vtxMC[2], inputCount); | |
335 | ||
336 | // stuff regarding the vertex reco correction and trigger bias correction | |
337 | fdNdEtaCorrection->FillEvent(vtxMC[2], inputCount, eventTriggered, eventVertex, processType); | |
338 | if (eventTriggered) { | |
339 | if (eventVertex) | |
340 | { | |
341 | fdNdEtaAnalysisESD->FillEvent(vtxMC[2], inputCount); | |
342 | } | |
343 | } | |
344 | ||
345 | delete[] etaArr; | |
346 | delete[] labelArr; | |
347 | delete[] ptArr; | |
348 | } | |
349 | ||
350 | void AlidNdEtaCorrectionTask::Terminate(Option_t *) | |
351 | { | |
352 | // The Terminate() function is the last function to be called during | |
353 | // a query. It always runs on the client, it can be used to present | |
354 | // the results graphically or save the results to file. | |
355 | ||
356 | fOutput = dynamic_cast<TList*> (GetOutputData(0)); | |
357 | if (!fOutput) { | |
358 | Printf("ERROR: fOutput not available"); | |
359 | return; | |
360 | } | |
361 | ||
362 | fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fOutput->FindObject("dndeta_correction")); | |
363 | fdNdEtaAnalysisMC = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndetaMC")); | |
364 | fdNdEtaAnalysisESD = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndetaESD")); | |
365 | if (!fdNdEtaCorrection || !fdNdEtaAnalysisMC || !fdNdEtaAnalysisESD) | |
366 | { | |
367 | AliDebug(AliLog::kError, "Could not read object from output list"); | |
368 | return; | |
369 | } | |
370 | ||
371 | fdNdEtaCorrection->Finish(); | |
372 | ||
373 | TString fileName; | |
374 | fileName.Form("correction_map%s.root", fOption.Data()); | |
375 | TFile* fout = new TFile(fileName, "RECREATE"); | |
376 | ||
377 | if (fEsdTrackCuts) | |
378 | fEsdTrackCuts->SaveHistograms("esd_track_cuts"); | |
379 | fdNdEtaCorrection->SaveHistograms(); | |
380 | fdNdEtaAnalysisMC->SaveHistograms(); | |
381 | fdNdEtaAnalysisESD->SaveHistograms(); | |
382 | ||
383 | fout->Write(); | |
384 | fout->Close(); | |
385 | ||
386 | fdNdEtaCorrection->DrawHistograms(); | |
387 | ||
388 | if (fPIDParticles && fPIDTracks) | |
389 | { | |
390 | new TCanvas("pidcanvas", "pidcanvas", 500, 500); | |
391 | ||
392 | fPIDParticles->Draw(); | |
393 | fPIDTracks->SetLineColor(2); | |
394 | fPIDTracks->Draw("SAME"); | |
395 | ||
396 | TDatabasePDG* pdgDB = new TDatabasePDG; | |
397 | ||
398 | for (Int_t i=0; i <= fPIDParticles->GetNbinsX()+1; ++i) | |
399 | if (fPIDParticles->GetBinContent(i) > 0) | |
400 | printf("PDG = %d (%s): generated: %d, reconstructed: %d, ratio: %f\n", (Int_t) fPIDParticles->GetBinCenter(i), pdgDB->GetParticle((Int_t) fPIDParticles->GetBinCenter(i))->GetName(), (Int_t) fPIDParticles->GetBinContent(i), (Int_t) fPIDTracks->GetBinContent(i), ((fPIDTracks->GetBinContent(i) > 0) ? fPIDParticles->GetBinContent(i) / fPIDTracks->GetBinContent(i) : -1)); | |
401 | ||
402 | delete pdgDB; | |
403 | pdgDB = 0; | |
404 | } | |
405 | ||
406 | Printf("Writting result to %s", fileName.Data()); | |
407 | } | |
408 | ||
409 | Bool_t AlidNdEtaCorrectionTask::SignOK(TParticlePDG* particle) | |
410 | { | |
411 | // returns if a particle with this sign should be counted | |
412 | // this is determined by the value of fSignMode, which should have the same sign | |
413 | // as the charge | |
414 | // if fSignMode is 0 all particles are counted | |
415 | ||
416 | if (fSignMode == 0) | |
417 | return kTRUE; | |
418 | ||
419 | if (!particle) | |
420 | { | |
421 | printf("WARNING: not counting a particle that does not have a pdg particle\n"); | |
422 | return kFALSE; | |
423 | } | |
424 | ||
425 | Double_t charge = particle->Charge(); | |
426 | ||
427 | if (fSignMode > 0) | |
428 | if (charge < 0) | |
429 | return kFALSE; | |
430 | ||
431 | if (fSignMode < 0) | |
432 | if (charge > 0) | |
433 | return kFALSE; | |
434 | ||
435 | return kTRUE; | |
436 | } |