merging post-CVS changes into SVN:
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaTask.cxx
1 /* $Id$ */
2
3 #include "AlidNdEtaTask.h"
4
5 #include <TStyle.h>
6 #include <TSystem.h>
7 #include <TCanvas.h>
8 #include <TVector3.h>
9 #include <TChain.h>
10 #include <TFile.h>
11 #include <TH1F.h>
12 #include <TH2F.h>
13 #include <TH3F.h>
14 #include <TParticle.h>
15 #include <TRandom.h>
16 #include <TNtuple.h>
17 #include <TObjString.h>
18 #include <TF1.h>
19
20 #include <AliLog.h>
21 #include <AliESDVertex.h>
22 #include <AliESDEvent.h>
23 #include <AliStack.h>
24 #include <AliHeader.h>
25 #include <AliGenEventHeader.h>
26 #include <AliMultiplicity.h>
27 #include <AliAnalysisManager.h>
28 #include <AliMCEventHandler.h>
29 #include <AliMCEvent.h>
30 #include <AliESDInputHandler.h>
31
32 #include "esdTrackCuts/AliESDtrackCuts.h"
33 #include "AliPWG0Helper.h"
34 #include "AliCorrection.h"
35 #include "AliCorrectionMatrix3D.h"
36 #include "dNdEta/dNdEtaAnalysis.h"
37
38 ClassImp(AlidNdEtaTask)
39
40 AlidNdEtaTask::AlidNdEtaTask(const char* opt) :
41   AliAnalysisTask("AlidNdEtaTask", ""),
42   fESD(0),
43   fOutput(0),
44   fOption(opt),
45   fAnalysisMode(kTPC),
46   fReadMC(kFALSE),
47   fEsdTrackCuts(0),
48   fdNdEtaAnalysisESD(0),
49   fMult(0),
50   fMultVtx(0),
51   fEvents(0),
52   fdNdEtaAnalysis(0),
53   fdNdEtaAnalysisTr(0),
54   fdNdEtaAnalysisTrVtx(0),
55   fdNdEtaAnalysisTracks(0),
56   fVertex(0),
57   fPartPt(0)
58 {
59   //
60   // Constructor. Initialization of pointers
61   //
62
63   // Define input and output slots here
64   DefineInput(0, TChain::Class());
65   DefineOutput(0, TList::Class());
66 }
67
68 AlidNdEtaTask::~AlidNdEtaTask()
69 {
70   //
71   // Destructor
72   //
73
74   // histograms are in the output list and deleted when the output
75   // list is deleted by the TSelector dtor
76
77   if (fOutput) {
78     delete fOutput;
79     fOutput = 0;
80   }
81 }
82
83 //________________________________________________________________________
84 void AlidNdEtaTask::ConnectInputData(Option_t *)
85 {
86   // Connect ESD
87   // Called once
88
89   Printf("AlidNdEtaTask::ConnectInputData called");
90
91   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
92   if (!tree) {
93     Printf("ERROR: Could not read tree from input slot 0");
94   } else {
95     // Disable all branches and enable only the needed ones
96     //tree->SetBranchStatus("*", 0);
97
98     tree->SetBranchStatus("fTriggerMask", 1);
99     tree->SetBranchStatus("fSPDVertex*", 1);
100
101     if (fAnalysisMode == kSPD)
102       tree->SetBranchStatus("fSPDMult*", 1);
103
104     if (fAnalysisMode == kTPC) {
105       AliESDtrackCuts::EnableNeededBranches(tree);
106       tree->SetBranchStatus("fTracks.fLabel", 1);
107     }
108
109     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
110
111     if (!esdH) {
112       Printf("ERROR: Could not get ESDInputHandler");
113     } else
114       fESD = esdH->GetEvent();
115   }
116 }
117
118 void AlidNdEtaTask::CreateOutputObjects()
119 {
120   // create result objects and add to output list
121
122   fOutput = new TList;
123   fOutput->SetOwner();
124
125   TString detector;
126   if (fAnalysisMode == kTPC)
127   {
128     detector = "TPC";
129   }
130   else if (fAnalysisMode == kSPD)
131     detector = "SPD";
132
133   fdNdEtaAnalysisESD = new dNdEtaAnalysis("fdNdEtaAnalysisESD", "fdNdEtaAnalysisESD", detector);
134   fOutput->Add(fdNdEtaAnalysisESD);
135
136   fMult = new TH1F("fMult", "fMult;Ntracks;Count", 201, -0.5, 200.5);
137   fOutput->Add(fMult);
138
139   fMultVtx = new TH1F("fMultVtx", "fMultVtx;Ntracks;Count", 201, -0.5, 200.5);
140   fOutput->Add(fMultVtx);
141
142   for (Int_t i=0; i<3; ++i)
143   {
144     fPartEta[i] = new TH1F(Form("dndeta_check_%d", i), Form("dndeta_check_%d", i), 60, -6, 6);
145     fPartEta[i]->Sumw2();
146     fOutput->Add(fPartEta[i]);
147   }
148
149   fEvents = new TH1F("dndeta_check_vertex", "dndeta_check_vertex", 40, -20, 20);
150   fOutput->Add(fEvents);
151
152   if (fReadMC)
153   {
154     fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta", detector);
155     fOutput->Add(fdNdEtaAnalysis);
156
157     fdNdEtaAnalysisTr = new dNdEtaAnalysis("dndetaTr", "dndetaTr", detector);
158     fOutput->Add(fdNdEtaAnalysisTr);
159
160     fdNdEtaAnalysisTrVtx = new dNdEtaAnalysis("dndetaTrVtx", "dndetaTrVtx", detector);
161     fOutput->Add(fdNdEtaAnalysisTrVtx);
162
163     fdNdEtaAnalysisTracks = new dNdEtaAnalysis("dndetaTracks", "dndetaTracks", detector);
164     fOutput->Add(fdNdEtaAnalysisTracks);
165
166     fVertex = new TH3F("vertex_check", "vertex_check", 50, -50, 50, 50, -50, 50, 50, -50, 50);
167     fOutput->Add(fVertex);
168
169     fPartPt =  new TH1F("dndeta_check_pt", "dndeta_check_pt", 1000, 0, 10);
170     fPartPt->Sumw2();
171     fOutput->Add(fPartPt);
172   }
173 }
174
175 void AlidNdEtaTask::Exec(Option_t*)
176 {
177   // process the event
178
179   // Check prerequisites
180   if (!fESD)
181   {
182     AliDebug(AliLog::kError, "ESD branch not available");
183     return;
184   }
185
186   // trigger definition
187   Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD->GetTriggerMask(), AliPWG0Helper::kMB2);
188
189   Bool_t eventVertex = AliPWG0Helper::IsVertexReconstructed(fESD->GetVertex());
190
191   // get the ESD vertex
192   const AliESDVertex* vtxESD = fESD->GetVertex();
193   Double_t vtx[3];
194   vtxESD->GetXYZ(vtx);
195
196   // post the data already here
197   PostData(0, fOutput);
198
199   // create list of (label, eta, pt) tuples
200   Int_t inputCount = 0;
201   Int_t* labelArr = 0;
202   Float_t* etaArr = 0;
203   Float_t* ptArr = 0;
204   if (fAnalysisMode == kSPD)
205   {
206     // get tracklets
207     const AliMultiplicity* mult = fESD->GetMultiplicity();
208     if (!mult)
209     {
210       AliDebug(AliLog::kError, "AliMultiplicity not available");
211       return;
212     }
213
214     labelArr = new Int_t[mult->GetNumberOfTracklets()];
215     etaArr = new Float_t[mult->GetNumberOfTracklets()];
216     ptArr = new Float_t[mult->GetNumberOfTracklets()];
217
218     // get multiplicity from ITS tracklets
219     for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
220     {
221       //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
222
223       // This removes non-tracklets in PDC06 data. Very bad solution. New solution is implemented for newer data. Keeping this for compatibility.
224       if (mult->GetDeltaPhi(i) < -1000)
225         continue;
226
227       etaArr[inputCount] = mult->GetEta(i);
228       labelArr[inputCount] = mult->GetLabel(i);
229       ptArr[inputCount] = 0; // no pt for tracklets
230       ++inputCount;
231     }
232   }
233   else if (fAnalysisMode == kTPC)
234   {
235     if (!fEsdTrackCuts)
236     {
237       AliDebug(AliLog::kError, "fESDTrackCuts not available");
238       return;
239     }
240
241     // get multiplicity from ESD tracks
242     TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
243     Int_t nGoodTracks = list->GetEntries();
244
245     labelArr = new Int_t[nGoodTracks];
246     etaArr = new Float_t[nGoodTracks];
247     ptArr = new Float_t[nGoodTracks];
248
249     // loop over esd tracks
250     for (Int_t i=0; i<nGoodTracks; i++)
251     {
252       AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(i));
253       if (!esdTrack)
254       {
255         AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", i));
256         continue;
257       }
258
259       etaArr[inputCount] = esdTrack->Eta();
260       labelArr[inputCount] = TMath::Abs(esdTrack->GetLabel());
261       ptArr[inputCount] = esdTrack->Pt();
262       ++inputCount;
263
264       Printf("%f", esdTrack->Pt());
265     }
266   }
267   else
268     return;
269
270   // Processing of ESD information (always)
271   if (eventTriggered)
272   {
273     // control hist
274     fMult->Fill(inputCount);
275     if (eventVertex)
276     {
277       // control hist
278       fMultVtx->Fill(inputCount);
279
280       for (Int_t i=0; i<inputCount; ++i)
281       {
282         Float_t eta = etaArr[i];
283         Float_t pt  = ptArr[i];
284
285         fdNdEtaAnalysisESD->FillTrack(vtx[2], eta, pt);
286
287         if (TMath::Abs(vtx[2]) < 20)
288         {
289           fPartEta[0]->Fill(eta);
290
291           if (vtx[2] < 0)
292             fPartEta[1]->Fill(eta);
293           else
294             fPartEta[2]->Fill(eta);
295         }
296       }
297
298       // for event count per vertex
299       fdNdEtaAnalysisESD->FillEvent(vtx[2], inputCount);
300
301       // control hists
302       fEvents->Fill(vtx[2]);
303     }
304   }
305
306   if (fReadMC)   // Processing of MC information (optional)
307   {
308     AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
309     if (!eventHandler) {
310       Printf("ERROR: Could not retrieve MC event handler");
311       return;
312     }
313
314     AliMCEvent* mcEvent = eventHandler->MCEvent();
315     if (!mcEvent) {
316       Printf("ERROR: Could not retrieve MC event");
317       return;
318     }
319
320     AliStack* stack = mcEvent->Stack();
321     if (!stack)
322     {
323       AliDebug(AliLog::kError, "Stack not available");
324       return;
325     }
326
327     AliHeader* header = mcEvent->Header();
328     if (!header)
329     {
330       AliDebug(AliLog::kError, "Header not available");
331       return;
332     }
333
334     // get the MC vertex
335     AliGenEventHeader* genHeader = header->GenEventHeader();
336     TArrayF vtxMC(3);
337     genHeader->PrimaryVertex(vtxMC);
338
339     // loop over mc particles
340     Int_t nPrim  = stack->GetNprimary();
341
342     Int_t nAcceptedParticles = 0;
343
344     for (Int_t iMc = 0; iMc < nPrim; ++iMc)
345     {
346       //Printf("Getting particle %d", iMc);
347       TParticle* particle = stack->Particle(iMc);
348
349       if (!particle)
350         continue;
351
352       if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
353         continue;
354
355       AliDebug(AliLog::kDebug+1, Form("Accepted primary %d, unique ID: %d", iMc, particle->GetUniqueID()));
356       ++nAcceptedParticles;
357
358       Float_t eta = particle->Eta();
359       Float_t pt = particle->Pt();
360
361       fdNdEtaAnalysis->FillTrack(vtxMC[2], eta, pt);
362       fVertex->Fill(particle->Vx(), particle->Vy(), particle->Vz());
363
364       if (eventTriggered)
365       {
366         fdNdEtaAnalysisTr->FillTrack(vtxMC[2], eta, pt);
367         if (eventVertex)
368           fdNdEtaAnalysisTrVtx->FillTrack(vtxMC[2], eta, pt);
369       }
370
371       if (TMath::Abs(eta) < 0.8)
372         fPartPt->Fill(pt);
373     }
374
375     fdNdEtaAnalysis->FillEvent(vtxMC[2], nAcceptedParticles);
376     if (eventTriggered)
377     {
378       fdNdEtaAnalysisTr->FillEvent(vtxMC[2], nAcceptedParticles);
379       if (eventVertex)
380         fdNdEtaAnalysisTrVtx->FillEvent(vtxMC[2], nAcceptedParticles);
381     }
382
383     if (eventTriggered && eventVertex)
384     {
385       // from tracks is only done for triggered and vertex reconstructed events
386
387       for (Int_t i=0; i<inputCount; ++i)
388       {
389         Int_t label = labelArr[i];
390
391         if (label < 0)
392           continue;
393
394         //Printf("Getting particle of track %d", label);
395         TParticle* particle = stack->Particle(label);
396
397         if (!particle)
398         {
399           AliDebug(AliLog::kError, Form("ERROR: Could not retrieve particle %d.", label));
400           continue;
401         }
402
403         fdNdEtaAnalysisTracks->FillTrack(vtxMC[2], particle->Eta(), particle->Pt());
404       } // end of track loop
405
406       // for event count per vertex
407       fdNdEtaAnalysisTracks->FillEvent(vtxMC[2], inputCount);
408     }
409   }
410
411   delete[] etaArr;
412   delete[] labelArr;
413   delete[] ptArr;
414 }
415
416 void AlidNdEtaTask::Terminate(Option_t *)
417 {
418   // The Terminate() function is the last function to be called during
419   // a query. It always runs on the client, it can be used to present
420   // the results graphically or save the results to file.
421
422   fOutput = dynamic_cast<TList*> (GetOutputData(0));
423   if (!fOutput) {
424     Printf("ERROR: fOutput not available");
425     return;
426   }
427
428   fdNdEtaAnalysisESD = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("fdNdEtaAnalysisESD"));
429   fMult = dynamic_cast<TH1F*> (fOutput->FindObject("fMult"));
430   fMultVtx = dynamic_cast<TH1F*> (fOutput->FindObject("fMultVtx"));
431   fEvents = dynamic_cast<TH1F*> (fOutput->FindObject("dndeta_check_vertex"));
432   for (Int_t i=0; i<3; ++i)
433     fPartEta[i] = dynamic_cast<TH1F*> (fOutput->FindObject(Form("dndeta_check_%d", i)));
434
435   if (!fdNdEtaAnalysisESD)
436   {
437     AliDebug(AliLog::kError, "ERROR: fdNdEtaAnalysisESD not available");
438     return;
439   }
440
441   if (fMult && fMultVtx)
442   {
443     new TCanvas;
444     fMult->Draw();
445     fMultVtx->SetLineColor(2);
446     fMultVtx->DrawCopy("SAME");
447   }
448
449   if (fPartEta[0])
450   {
451     Int_t events1 = (Int_t) fEvents->Integral(fEvents->GetXaxis()->FindBin(-19.9), fEvents->GetXaxis()->FindBin(-0.001));
452     Int_t events2 = (Int_t) fEvents->Integral(fEvents->GetXaxis()->FindBin(0.001), fEvents->GetXaxis()->FindBin(19.9));
453
454     fPartEta[0]->Scale(1.0 / (events1 + events2));
455     fPartEta[1]->Scale(1.0 / events1);
456     fPartEta[2]->Scale(1.0 / events2);
457
458     for (Int_t i=0; i<3; ++i)
459       fPartEta[i]->Scale(1.0 / fPartEta[i]->GetBinWidth(1));
460
461     new TCanvas("control", "control", 500, 500);
462     for (Int_t i=0; i<3; ++i)
463     {
464       fPartEta[i]->SetLineColor(i+1);
465       fPartEta[i]->Draw((i==0) ? "" : "SAME");
466     }
467   }
468
469   if (fEvents)
470   {
471     new TCanvas("control3", "control3", 500, 500);
472     fEvents->DrawCopy();
473   }
474
475   TFile* fout = new TFile("analysis_esd_raw.root", "RECREATE");
476
477   if (fdNdEtaAnalysisESD)
478     fdNdEtaAnalysisESD->SaveHistograms();
479
480   if (fEsdTrackCuts)
481     fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");
482
483   if (fMult)
484     fMult->Write();
485
486   if (fMultVtx)
487     fMultVtx->Write();
488
489   for (Int_t i=0; i<3; ++i)
490     if (fPartEta[i])
491       fPartEta[i]->Write();
492
493   if (fEvents)
494     fEvents->Write();
495
496   fout->Write();
497   fout->Close();
498
499   Printf("Writting result to analysis_esd_raw.root");
500
501   if (fReadMC)
502   {
503     fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
504     fdNdEtaAnalysisTr = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndetaTr"));
505     fdNdEtaAnalysisTrVtx = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndetaTrVtx"));
506     fdNdEtaAnalysisTracks = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndetaTracks"));
507     fPartPt = dynamic_cast<TH1F*> (fOutput->FindObject("dndeta_check_pt"));
508
509     if (!fdNdEtaAnalysis || !fdNdEtaAnalysisTr || !fdNdEtaAnalysisTrVtx || !fPartPt)
510     {
511       AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p %p", (void*) fdNdEtaAnalysis, (void*) fPartPt));
512       return;
513     }
514
515     fdNdEtaAnalysis->Finish(0, -1, AlidNdEtaCorrection::kNone);
516     fdNdEtaAnalysisTr->Finish(0, -1, AlidNdEtaCorrection::kNone);
517     fdNdEtaAnalysisTrVtx->Finish(0, -1, AlidNdEtaCorrection::kNone);
518     fdNdEtaAnalysisTracks->Finish(0, -1, AlidNdEtaCorrection::kNone);
519
520     Int_t events = (Int_t) fdNdEtaAnalysis->GetData()->GetEventCorrection()->GetMeasuredHistogram()->Integral();
521     fPartPt->Scale(1.0/events);
522     fPartPt->Scale(1.0/fPartPt->GetBinWidth(1));
523
524     TFile* fout = new TFile("analysis_mc.root","RECREATE");
525
526     fdNdEtaAnalysis->SaveHistograms();
527     fdNdEtaAnalysisTr->SaveHistograms();
528     fdNdEtaAnalysisTrVtx->SaveHistograms();
529     fdNdEtaAnalysisTracks->SaveHistograms();
530     fPartPt->Write();
531
532     fout->Write();
533     fout->Close();
534
535     Printf("Writting result to analysis_mc.root");
536
537     if (fPartPt)
538     {
539       new TCanvas("control2", "control2", 500, 500);
540       fPartPt->DrawCopy();
541     }
542   }
543 }