]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaSystematicsSelector.cxx
Changes:
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaSystematicsSelector.cxx
1 /* $Id$ */
2
3 #include "AlidNdEtaSystematicsSelector.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 <TH2F.h>
12 #include <TH1F.h>
13 #include <TParticle.h>
14
15 #include <AliLog.h>
16 #include <AliESD.h>
17 #include <AliGenEventHeader.h>
18 #include <AliStack.h>
19 #include <AliHeader.h>
20
21 #include "esdTrackCuts/AliESDtrackCuts.h"
22 #include "AliPWG0Helper.h"
23 #include "dNdEta/AlidNdEtaCorrection.h"
24
25 ClassImp(AlidNdEtaSystematicsSelector)
26
27 AlidNdEtaSystematicsSelector::AlidNdEtaSystematicsSelector() :
28   AliSelectorRL(),
29   fSecondaries(0),
30   fSigmaVertex(0),
31   fEsdTrackCuts(0),
32   fPIDParticles(0),
33   fPIDTracks(0)
34 {
35   //
36   // Constructor. Initialization of pointers
37   //
38
39   for (Int_t i=0; i<4; ++i)
40     fdNdEtaCorrection[i] = 0;
41 }
42
43 AlidNdEtaSystematicsSelector::~AlidNdEtaSystematicsSelector()
44 {
45   //
46   // Destructor
47   //
48 }
49
50 void AlidNdEtaSystematicsSelector::Begin(TTree* tree)
51 {
52   // Begin function
53
54   AliSelectorRL::Begin(tree);
55
56   ReadUserObjects(tree);
57 }
58
59 void AlidNdEtaSystematicsSelector::ReadUserObjects(TTree* tree)
60 {
61   // read the user objects, called from slavebegin and begin
62
63   if (!fEsdTrackCuts && fInput)
64     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
65
66   if (!fEsdTrackCuts && tree)
67     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
68
69   if (!fEsdTrackCuts)
70      AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
71 }
72
73 void AlidNdEtaSystematicsSelector::SlaveBegin(TTree* tree)
74 {
75   // The SlaveBegin() function is called after the Begin() function.
76   // When running with PROOF SlaveBegin() is called on each slave server.
77   // The tree argument is deprecated (on PROOF 0 is passed).
78
79   AliSelector::SlaveBegin(tree);
80
81   ReadUserObjects(tree);
82
83   TString option(GetOption());
84
85   printf("Running AlidNdEtaSystematicsSelector with options %s\n", option.Data());
86
87   if (option.Contains("secondaries"))
88   {
89     fSecondaries = new TH2F("fSecondaries", "fSecondaries;Case;N", 8, 0.5, 8.5, 2001, -0.5, 2000.5);
90     fSecondaries->GetXaxis()->SetBinLabel(1, "All Primaries");
91     fSecondaries->GetXaxis()->SetBinLabel(2, "All Secondaries");
92     fSecondaries->GetXaxis()->SetBinLabel(3, "Primaries pT > 0.3 GeV/c");
93     fSecondaries->GetXaxis()->SetBinLabel(4, "Secondaries pT > 0.3 GeV/c");
94     fSecondaries->GetXaxis()->SetBinLabel(5, "Tracks from Primaries");
95     fSecondaries->GetXaxis()->SetBinLabel(6, "Tracks from Secondaries");
96     fSecondaries->GetXaxis()->SetBinLabel(7, "Accepted Tracks from Primaries");
97     fSecondaries->GetXaxis()->SetBinLabel(8, "Accepted Tracks from Secondaries");
98   }
99
100   if (option.Contains("particle-composition"))
101   {
102     for (Int_t i=0; i<4; ++i)
103     {
104       TString name;
105       name.Form("correction_%d", i);
106       fdNdEtaCorrection[i] = new AlidNdEtaCorrection(name, name);
107     }
108   }
109
110   if (option.Contains("sigma-vertex"))
111   {
112     fSigmaVertex = new TH1F("fSigmaVertex", "fSigmaVertex;Nsigma2vertex;NacceptedTracks", 10, 0.25, 5.25);
113     printf("WARNING: sigma-vertex analysis enabled. This will produce weird results in the AliESDtrackCuts histograms\n");
114   }
115
116   fPIDParticles = new TH1F("pid_particles", "PID of generated primary particles", 10001, -5000.5, 5000.5);
117   fPIDTracks = new TH1F("pid_tracks", "MC PID of reconstructed tracks", 10001, -5000.5, 5000.5);
118 }
119
120 Bool_t AlidNdEtaSystematicsSelector::Process(Long64_t entry)
121 {
122   // The Process() function is called for each entry in the tree (or possibly
123   // keyed object in the case of PROOF) to be processed. The entry argument
124   // specifies which entry in the currently loaded tree is to be processed.
125   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
126   // to read either all or the required parts of the data. When processing
127   // keyed objects with PROOF, the object is already loaded and is available
128   // via the fObject pointer.
129   //
130   // This function should contain the "body" of the analysis. It can contain
131   // simple or elaborate selection criteria, run algorithms on the data
132   // of the event and typically fill histograms.
133
134   // WARNING when a selector is used with a TChain, you must use
135   //  the pointer to the current TTree to call GetEntry(entry).
136   //  The entry is always the local entry number in the current tree.
137   //  Assuming that fTree is the pointer to the TChain being processed,
138   //  use fTree->GetTree()->GetEntry(entry).
139
140   if (AliSelectorRL::Process(entry) == kFALSE)
141     return kFALSE;
142
143   // Check prerequisites
144   if (!fESD)
145   {
146     AliDebug(AliLog::kError, "ESD branch not available");
147     return kFALSE;
148   }
149
150   if (!fEsdTrackCuts)
151   {
152     AliDebug(AliLog::kError, "fESDTrackCuts not available");
153     return kFALSE;
154   }
155
156   AliStack* stack = GetStack();
157   if (!stack)
158   {
159     AliDebug(AliLog::kError, "Stack not available");
160     return kFALSE;
161   }
162
163   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
164
165   if (fdNdEtaCorrection[0])
166     FillCorrectionMaps(list);
167
168   if (fSecondaries)
169     FillSecondaries();
170
171   if (fSigmaVertex)
172     FillSigmaVertex();
173
174   delete list;
175   list = 0;
176
177   return kTRUE;
178 }
179
180 void AlidNdEtaSystematicsSelector::FillCorrectionMaps(TObjArray* listOfTracks)
181 {
182   // fills the correction maps for different particle species
183
184   AliStack* stack = GetStack();
185   AliHeader* header = GetHeader();
186
187   Bool_t vertexReconstructed = AliPWG0Helper::IsVertexReconstructed(fESD);
188   Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD);
189
190   // get the MC vertex
191   AliGenEventHeader* genHeader = header->GenEventHeader();
192
193   TArrayF vtxMC(3);
194   genHeader->PrimaryVertex(vtxMC);
195
196   // loop over mc particles
197   Int_t nPrim  = stack->GetNprimary();
198
199   for (Int_t iMc = 0; iMc < nPrim; ++iMc)
200   {
201     TParticle* particle = stack->Particle(iMc);
202
203     if (!particle)
204     {
205       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
206       continue;
207     }
208
209     if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
210       continue;
211
212     Float_t eta = particle->Eta();
213     Float_t pt = particle->Pt();
214
215     Int_t id = -1;
216     switch (TMath::Abs(particle->GetPdgCode()))
217     {
218       case 211: id = 0; break;
219       case 321: id = 1; break;
220       case 2212: id = 2; break;
221       case 11:
222         {
223           /*if (pt < 0.1 && particle->GetMother(0) > -1)
224           {
225             TParticle* mother = stack->Particle(particle->GetMother(0));
226             printf("Mother pdg code is %d\n", mother->GetPdgCode());
227           } */
228           //particle->Dump();
229           //if (particle->GetUniqueID() == 1)
230
231           //printf("Found illegal particle mc id = %d file = %s event = %d\n", iMc,  fTree->GetCurrentFile()->GetName(), fTree->GetTree()->GetReadEntry());
232         }
233       default: id = 3; break;
234     }
235
236     if (vertexReconstructed)
237     {
238       fdNdEtaCorrection[id]->FillParticle(vtxMC[2], eta, pt);
239       //if (pt < 0.1)
240         fPIDParticles->Fill(particle->GetPdgCode());
241     }
242
243     fdNdEtaCorrection[id]->FillParticleAllEvents(eta, pt);
244     if (eventTriggered)
245       fdNdEtaCorrection[id]->FillParticleWhenEventTriggered(eta, pt);
246   }// end of mc particle
247
248   // loop over esd tracks
249   TIterator* iter = listOfTracks->MakeIterator();
250   TObject* obj = 0;
251   while ((obj = iter->Next()) != 0)
252   {
253     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
254     if (!esdTrack)
255       continue;
256
257     // using the properties of the mc particle
258     Int_t label = TMath::Abs(esdTrack->GetLabel());
259     TParticle* particle = stack->Particle(label);
260     if (!particle)
261     {
262       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
263       continue;
264     }
265
266     TParticle* mother = particle;
267     // find primary particle that created this particle
268     while (label >= nPrim)
269     {
270       //printf("Particle %d (pdg %d) is not a primary. Let's check its mother %d\n", label, mother->GetPdgCode(), mother->GetMother(0));
271
272       if (mother->GetMother(0) == -1)
273       {
274         AliDebug(AliLog::kError, Form("UNEXPECTED: Could not find mother of secondary particle %d.", label));
275         mother = 0;
276         break;
277       }
278
279       label = mother->GetMother(0);
280
281       mother = stack->Particle(label);
282       if (!mother)
283       {
284         AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (find mother loop).", label));
285         break;
286       }
287     }
288
289     if (!mother)
290       continue;
291
292     //printf("We continue with particle %d (pdg %d)\n", label, mother->GetPdgCode());
293
294     Int_t id = -1;
295     switch (TMath::Abs(mother->GetPdgCode()))
296     {
297       case 211:  id = 0; break;
298       case 321:  id = 1; break;
299       case 2212: id = 2; break;
300       default:   id = 3; break;
301     }
302
303     if (vertexReconstructed)
304     {
305       fdNdEtaCorrection[id]->FillParticleWhenMeasuredTrack(vtxMC[2], particle->Eta(), particle->Pt());
306       //if (particle->Pt() < 0.1)
307         fPIDTracks->Fill(particle->GetPdgCode());
308     }
309   } // end of track loop
310
311   Int_t nGoodTracks = listOfTracks->GetEntries();
312
313   for (Int_t i=0; i<4; ++i)
314   {
315     fdNdEtaCorrection[i]->FillEvent(vtxMC[2], nGoodTracks);
316     if (eventTriggered)
317     {
318       fdNdEtaCorrection[i]->FillEventWithTrigger(vtxMC[2], nGoodTracks);
319       if (vertexReconstructed)
320         fdNdEtaCorrection[i]->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
321     }
322   }
323
324   delete iter;
325   iter = 0;
326 }
327
328 void AlidNdEtaSystematicsSelector::FillSecondaries()
329 {
330   // fills the secondary histograms
331
332   AliStack* stack = GetStack();
333
334   Int_t particlesPrimaries = 0;
335   Int_t particlesSecondaries = 0;
336   Int_t particlesPrimariesPtCut = 0;
337   Int_t particlesSecondariesPtCut = 0;
338
339   for (Int_t iParticle = 0; iParticle < stack->GetNtrack(); iParticle++)
340   {
341     TParticle* particle = stack->Particle(iParticle);
342     if (!particle)
343     {
344       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (particle loop).", iParticle));
345       continue;
346     }
347
348     if (TMath::Abs(particle->Eta()) > 0.9)
349       continue;
350
351     Bool_t isPrimary = kFALSE;
352     if (iParticle < stack->GetNprimary())
353     {
354       if (AliPWG0Helper::IsPrimaryCharged(particle, stack->GetNprimary()) == kFALSE)
355         continue;
356
357       isPrimary = kTRUE;
358     }
359
360     if (isPrimary)
361       particlesPrimaries++;
362     else
363       particlesSecondaries++;
364
365     if (particle->Pt() < 0.3)
366       continue;
367
368     if (isPrimary)
369       particlesPrimariesPtCut++;
370     else
371       particlesSecondariesPtCut++;
372   }
373
374   fSecondaries->Fill(1, particlesPrimaries);
375   fSecondaries->Fill(2, particlesSecondaries);
376   fSecondaries->Fill(3, particlesPrimariesPtCut);
377   fSecondaries->Fill(4, particlesSecondariesPtCut);
378
379   Int_t allTracksPrimaries = 0;
380   Int_t allTracksSecondaries = 0;
381   Int_t acceptedTracksPrimaries = 0;
382   Int_t acceptedTracksSecondaries = 0;
383
384   for (Int_t iTrack = 0; iTrack < fESD->GetNumberOfTracks(); iTrack++)
385   {
386     AliESDtrack* esdTrack = fESD->GetTrack(iTrack);
387     if (!esdTrack)
388     {
389       AliDebug(AliLog::kError, Form("UNEXPECTED: Could not get track %d.", iTrack));
390       continue;
391     }
392
393     Int_t label = TMath::Abs(esdTrack->GetLabel());
394     TParticle* particle = stack->Particle(label);
395     if (!particle)
396     {
397       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
398       continue;
399     }
400
401    if (label < stack->GetNprimary())
402       allTracksPrimaries++;
403     else
404       allTracksSecondaries++;
405
406     if (!fEsdTrackCuts->AcceptTrack(esdTrack))
407       continue;
408
409     if (label < stack->GetNprimary())
410       acceptedTracksPrimaries++;
411     else
412       acceptedTracksSecondaries++;
413   }
414
415   fSecondaries->Fill(5, allTracksPrimaries);
416   fSecondaries->Fill(6, allTracksSecondaries);
417   fSecondaries->Fill(7, acceptedTracksPrimaries);
418   fSecondaries->Fill(8, acceptedTracksSecondaries);
419
420   //printf("P = %d, S = %d, P_pt = %d, S_pt = %d, T_P = %d, T_S = %d, T_P_acc = %d, T_S_acc = %d\n", particlesPrimaries, particlesSecondaries, particlesPrimariesPtCut, particlesSecondariesPtCut, allTracksPrimaries, allTracksSecondaries, acceptedTracksPrimaries, acceptedTracksSecondaries);
421 }
422
423 void AlidNdEtaSystematicsSelector::FillSigmaVertex()
424 {
425   // fills the fSigmaVertex histogram
426
427   // save the old value
428   Float_t oldSigmaVertex = fEsdTrackCuts->GetMinNsigmaToVertex();
429
430   // set to maximum
431   fEsdTrackCuts->SetMinNsigmaToVertex(5);
432
433   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
434
435   TIterator* iter = list->MakeIterator();
436   TObject* obj = 0;
437   while ((obj = iter->Next()) != 0)
438   {
439     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
440     if (!esdTrack)
441       continue;
442
443     Float_t sigma2Vertex = fEsdTrackCuts->GetSigmaToVertex(esdTrack);
444
445     for (Double_t nSigma = 0.5; nSigma < 5.1; nSigma += 0.5)
446     {
447       if (sigma2Vertex < nSigma)
448         fSigmaVertex->Fill(nSigma);
449     }
450   }
451
452   delete iter;
453   iter = 0;
454
455   delete list;
456   list = 0;
457
458   // set back the old value
459   fEsdTrackCuts->SetMinNsigmaToVertex(oldSigmaVertex);
460 }
461
462 void AlidNdEtaSystematicsSelector::SlaveTerminate()
463 {
464   // The SlaveTerminate() function is called after all entries or objects
465   // have been processed. When running with PROOF SlaveTerminate() is called
466   // on each slave server.
467
468   AliSelector::SlaveTerminate();
469
470   // Add the histograms to the output on each slave server
471   if (!fOutput)
472   {
473     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
474     return;
475   }
476
477   if (fSecondaries)
478     fOutput->Add(fSecondaries);
479
480   for (Int_t i=0; i<4; ++i)
481     if (fdNdEtaCorrection[i])
482       fOutput->Add(fdNdEtaCorrection[i]);
483
484   if (fSigmaVertex)
485     fOutput->Add(fSigmaVertex);
486 }
487
488 void AlidNdEtaSystematicsSelector::Terminate()
489 {
490   // The Terminate() function is the last function to be called during
491   // a query. It always runs on the client, it can be used to present
492   // the results graphically or save the results to file.
493
494   AliSelector::Terminate();
495
496   fSecondaries = dynamic_cast<TH2F*> (fOutput->FindObject("fSecondaries"));
497   for (Int_t i=0; i<4; ++i)
498     fdNdEtaCorrection[i] = dynamic_cast<AlidNdEtaCorrection*> (fOutput->FindObject(Form("correction_%d", i)));
499   fSigmaVertex = dynamic_cast<TH1F*> (fOutput->FindObject("fSigmaVertex"));
500
501   if (fPIDParticles)
502   {
503     TDatabasePDG* pdgDB = new TDatabasePDG;
504
505     for (Int_t i=0; i <= fPIDParticles->GetNbinsX()+1; ++i)
506       if (fPIDParticles->GetBinContent(i) > 0)
507         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));
508
509     delete pdgDB;
510     pdgDB = 0;
511   }
512
513   TFile* fout = TFile::Open("systematics.root", "RECREATE");
514
515   if (fEsdTrackCuts)
516     fEsdTrackCuts->SaveHistograms("esd_track_cuts");
517
518   if (fSecondaries)
519     fSecondaries->Write();
520
521   if (fSigmaVertex)
522     fSigmaVertex->Write();
523
524   for (Int_t i=0; i<4; ++i)
525     if (fdNdEtaCorrection[i])
526       fdNdEtaCorrection[i]->SaveHistograms();
527
528   fout->Write();
529   fout->Close();
530
531   if (fSecondaries)
532   {
533     new TCanvas;
534     fSecondaries->Draw("COLZ");
535   }
536
537   if (fSigmaVertex)
538   {
539     new TCanvas;
540     fSigmaVertex->Draw();
541   }
542 }