adding nsigma study
[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 <TH3F.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 "AlidNdEtaCorrection.h"
24
25 ClassImp(AlidNdEtaSystematicsSelector)
26
27 AlidNdEtaSystematicsSelector::AlidNdEtaSystematicsSelector() :
28   AliSelectorRL(),
29   fSecondaries(0),
30   fSigmaVertex(0),
31   fEsdTrackCuts(0)
32 {
33   //
34   // Constructor. Initialization of pointers
35   //
36
37   for (Int_t i=0; i<4; ++i)
38     fdNdEtaCorrection[i] = 0;
39 }
40
41 AlidNdEtaSystematicsSelector::~AlidNdEtaSystematicsSelector()
42 {
43   //
44   // Destructor
45   //
46 }
47
48 void AlidNdEtaSystematicsSelector::Begin(TTree* tree)
49 {
50   // Begin function
51
52   AliSelectorRL::Begin(tree);
53
54   ReadUserObjects(tree);
55 }
56
57 void AlidNdEtaSystematicsSelector::ReadUserObjects(TTree* tree)
58 {
59   // read the user objects, called from slavebegin and begin
60
61   if (!fEsdTrackCuts && fInput)
62     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
63
64   if (!fEsdTrackCuts && tree)
65     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
66
67   if (!fEsdTrackCuts)
68      AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
69 }
70
71 void AlidNdEtaSystematicsSelector::SlaveBegin(TTree* tree)
72 {
73   // The SlaveBegin() function is called after the Begin() function.
74   // When running with PROOF SlaveBegin() is called on each slave server.
75   // The tree argument is deprecated (on PROOF 0 is passed).
76
77   AliSelector::SlaveBegin(tree);
78
79   ReadUserObjects(tree);
80
81   TString option(GetOption());
82
83   printf("Running AlidNdEtaSystematicsSelector with options %s\n", option.Data());
84
85   if (option.Contains("secondaries"))
86   {
87     fSecondaries = new TH3F("fSecondaries", "fSecondaries;NacceptedTracks;CratioSecondaries;p_{T}", 2000, -0.5, 205.5, 16, 0.45, 2.05, 10, 0, 10);
88   }
89
90   if (option.Contains("particle-composition"))
91   {
92     for (Int_t i=0; i<4; ++i)
93     {
94       TString name;
95       name.Form("correction_%d", i);
96       fdNdEtaCorrection[i] = new AlidNdEtaCorrection(name, name);
97     }
98   }
99
100   if (option.Contains("sigma-vertex"))
101   {
102     fSigmaVertex = new TH1F("fSigmaVertex", "fSigmaVertex;Nsigma2vertex;NacceptedTracks", 10, 0.25, 5.25);
103     printf("WARNING: sigma-vertex analysis enabled. This will produce weird results in the AliESDtrackCuts histograms\n");
104   }
105 }
106
107 Bool_t AlidNdEtaSystematicsSelector::Process(Long64_t entry)
108 {
109   // The Process() function is called for each entry in the tree (or possibly
110   // keyed object in the case of PROOF) to be processed. The entry argument
111   // specifies which entry in the currently loaded tree is to be processed.
112   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
113   // to read either all or the required parts of the data. When processing
114   // keyed objects with PROOF, the object is already loaded and is available
115   // via the fObject pointer.
116   //
117   // This function should contain the "body" of the analysis. It can contain
118   // simple or elaborate selection criteria, run algorithms on the data
119   // of the event and typically fill histograms.
120
121   // WARNING when a selector is used with a TChain, you must use
122   //  the pointer to the current TTree to call GetEntry(entry).
123   //  The entry is always the local entry number in the current tree.
124   //  Assuming that fTree is the pointer to the TChain being processed,
125   //  use fTree->GetTree()->GetEntry(entry).
126
127   if (AliSelectorRL::Process(entry) == kFALSE)
128     return kFALSE;
129
130   // Check prerequisites
131   if (!fESD)
132   {
133     AliDebug(AliLog::kError, "ESD branch not available");
134     return kFALSE;
135   }
136
137   if (!fEsdTrackCuts)
138   {
139     AliDebug(AliLog::kError, "fESDTrackCuts not available");
140     return kFALSE;
141   }
142
143   AliStack* stack = GetStack();
144   if (!stack)
145   {
146     AliDebug(AliLog::kError, "Stack not available");
147     return kFALSE;
148   }
149
150   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
151
152   if (fdNdEtaCorrection[0])
153     FillCorrectionMaps(list);
154
155   if (fSecondaries)
156     FillSecondaries(list);
157
158   if (fSigmaVertex)
159     FillSigmaVertex();
160
161   delete list;
162   list = 0;
163
164   return kTRUE;
165 }
166
167 void AlidNdEtaSystematicsSelector::FillCorrectionMaps(TObjArray* listOfTracks)
168 {
169   // fills the correction maps for different particle species
170
171   AliStack* stack = GetStack();
172   AliHeader* header = GetHeader();
173
174   Bool_t vertexReconstructed = AliPWG0Helper::IsVertexReconstructed(fESD);
175   Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD);
176
177   for (Int_t i=0; i<4; ++i)
178   {
179     fdNdEtaCorrection[i]->IncreaseEventCount();
180     if (eventTriggered)
181       fdNdEtaCorrection[i]->IncreaseTriggeredEventCount();
182   }
183
184   // get the MC vertex
185   AliGenEventHeader* genHeader = header->GenEventHeader();
186
187   TArrayF vtxMC(3);
188   genHeader->PrimaryVertex(vtxMC);
189
190   // loop over mc particles
191   Int_t nPrim  = stack->GetNprimary();
192
193   for (Int_t iMc = 0; iMc < nPrim; ++iMc)
194   {
195     TParticle* particle = stack->Particle(iMc);
196
197     if (!particle)
198     {
199       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
200       continue;
201     }
202
203     if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
204       continue;
205
206     Float_t eta = particle->Eta();
207     Float_t pt = particle->Pt();
208
209     Int_t id = -1;
210     switch (TMath::Abs(particle->GetPdgCode()))
211     {
212       case 211: id = 0; break;
213       case 321: id = 1; break;
214       case 2212: id = 2; break;
215       default: id = 3; break;
216     }
217
218     if (vertexReconstructed)
219       fdNdEtaCorrection[id]->FillParticle(vtxMC[2], eta, pt);
220
221     fdNdEtaCorrection[id]->FillParticleAllEvents(eta, pt);
222     if (eventTriggered)
223       fdNdEtaCorrection[id]->FillParticleWhenEventTriggered(eta, pt);
224   }// end of mc particle
225
226   // loop over esd tracks
227   TIterator* iter = listOfTracks->MakeIterator();
228   TObject* obj = 0;
229   while ((obj = iter->Next()) != 0)
230   {
231     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
232     if (!esdTrack)
233       continue;
234
235     // using the properties of the mc particle
236     Int_t label = TMath::Abs(esdTrack->GetLabel());
237     if (label == 0)
238     {
239       AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", label));
240       continue;
241     }
242
243     TParticle* particle = stack->Particle(label);
244     if (!particle)
245     {
246       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
247       continue;
248     }
249
250     Int_t id = -1;
251     switch (TMath::Abs(particle->GetPdgCode()))
252     {
253       case 211:  id = 0; break;
254       case 321:  id = 1; break;
255       case 2212: id = 2; break;
256       default:   id = 3; break;
257     }
258
259     if (vertexReconstructed)
260       fdNdEtaCorrection[id]->FillParticleWhenMeasuredTrack(vtxMC[2], particle->Eta(), particle->Pt());
261   } // end of track loop
262
263   delete iter;
264   iter = 0;
265 }
266
267 void AlidNdEtaSystematicsSelector::FillSecondaries(TObjArray* listOfTracks)
268 {
269   // fills the secondary histograms
270
271   AliStack* stack = GetStack();
272
273   TH1* nPrimaries = new TH1F("secondaries_primaries", "secondaries_primaries", fSecondaries->GetZaxis()->GetNbins(), fSecondaries->GetZaxis()->GetXmin(), fSecondaries->GetZaxis()->GetXmax());
274   TH1* nSecondaries = new TH1F("secondaries_secondaries", "secondaries_secondaries", fSecondaries->GetZaxis()->GetNbins(), fSecondaries->GetZaxis()->GetXmin(), fSecondaries->GetZaxis()->GetXmax());
275
276   TIterator* iter = listOfTracks->MakeIterator();
277   TObject* obj = 0;
278   while ((obj = iter->Next()) != 0)
279   {
280     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
281     if (!esdTrack)
282       continue;
283
284     Int_t label = TMath::Abs(esdTrack->GetLabel());
285     if (label == 0)
286     {
287       AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", label));
288       continue;
289     }
290
291     TParticle* particle = stack->Particle(label);
292     if (!particle)
293     {
294       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
295       continue;
296     }
297
298     Int_t nPrim  = stack->GetNprimary();
299     if (label < nPrim)
300       nPrimaries->Fill(particle->Pt());
301     else
302       nSecondaries->Fill(particle->Pt());
303   }
304
305   for (Int_t i=1; i<=nPrimaries->GetNbinsX(); ++i)
306   {
307     Int_t primaries = (Int_t) nPrimaries->GetBinContent(i);
308     Int_t secondaries = (Int_t) nSecondaries->GetBinContent(i);
309
310     if (primaries + secondaries > 0)
311     {
312       AliDebug(AliLog::kDebug, Form("The ratio between primaries and secondaries is %d:%d = %f", primaries, secondaries, ((secondaries > 0) ? (Double_t) primaries / secondaries : -1)));
313
314       for (Double_t factor = 0.5; factor < 2.01; factor += 0.1)
315       {
316         Double_t nTracks = (Double_t) primaries + (Double_t) secondaries * factor;
317         fSecondaries->Fill(nTracks, factor, nPrimaries->GetBinCenter(i));
318         //if (secondaries > 0) printf("We fill: %f %f %f\n", nTracks, factor, nPrimaries->GetBinCenter(i));
319       }
320     }
321   }
322
323   delete nPrimaries;
324   nPrimaries = 0;
325
326   delete nSecondaries;
327   nSecondaries = 0;
328
329   delete iter;
330   iter = 0;
331 }
332
333 void AlidNdEtaSystematicsSelector::FillSigmaVertex()
334 {
335   // fills the fSigmaVertex histogram
336
337   // save the old value
338   Float_t oldSigmaVertex = fEsdTrackCuts->GetMinNsigmaToVertex();
339
340   // set to maximum
341   fEsdTrackCuts->SetMinNsigmaToVertex(5);
342
343   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
344
345   TIterator* iter = list->MakeIterator();
346   TObject* obj = 0;
347   while ((obj = iter->Next()) != 0)
348   {
349     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
350     if (!esdTrack)
351       continue;
352
353     Float_t sigma2Vertex = fEsdTrackCuts->GetSigmaToVertex(esdTrack);
354
355     for (Double_t nSigma = 0.5; nSigma < 5.1; nSigma += 0.5)
356     {
357       if (sigma2Vertex < nSigma)
358         fSigmaVertex->Fill(nSigma);
359     }
360   }
361
362   delete iter;
363   iter = 0;
364
365   delete list;
366   list = 0;
367
368   // set back the old value
369   fEsdTrackCuts->SetMinNsigmaToVertex(oldSigmaVertex);
370 }
371
372 void AlidNdEtaSystematicsSelector::SlaveTerminate()
373 {
374   // The SlaveTerminate() function is called after all entries or objects
375   // have been processed. When running with PROOF SlaveTerminate() is called
376   // on each slave server.
377
378   AliSelector::SlaveTerminate();
379
380   // Add the histograms to the output on each slave server
381   if (!fOutput)
382   {
383     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
384     return;
385   }
386
387   if (fSecondaries)
388     fOutput->Add(fSecondaries);
389
390   for (Int_t i=0; i<4; ++i)
391     if (fdNdEtaCorrection[i])
392       fOutput->Add(fdNdEtaCorrection[i]);
393
394   if (fSigmaVertex)
395     fOutput->Add(fSigmaVertex);
396 }
397
398 void AlidNdEtaSystematicsSelector::Terminate()
399 {
400   // The Terminate() function is the last function to be called during
401   // a query. It always runs on the client, it can be used to present
402   // the results graphically or save the results to file.
403
404   AliSelector::Terminate();
405
406   fSecondaries = dynamic_cast<TH3F*> (fOutput->FindObject("fSecondaries"));
407   for (Int_t i=0; i<4; ++i)
408     fdNdEtaCorrection[i] = dynamic_cast<AlidNdEtaCorrection*> (fOutput->FindObject(Form("correction_%d", i)));
409   fSigmaVertex = dynamic_cast<TH1F*> (fOutput->FindObject("fSigmaVertex"));
410
411   TFile* fout = TFile::Open("systematics.root", "RECREATE");
412
413   if (fEsdTrackCuts)
414     fEsdTrackCuts->SaveHistograms("esd_track_cuts");
415
416   if (fSecondaries)
417     fSecondaries->Write();
418
419   if (fSigmaVertex)
420     fSigmaVertex->Write();
421
422   for (Int_t i=0; i<4; ++i)
423     if (fdNdEtaCorrection[i])
424       fdNdEtaCorrection[i]->SaveHistograms();
425
426   fout->Write();
427   fout->Close();
428
429   if (fSecondaries)
430   {
431     new TCanvas;
432     fSecondaries->Draw();
433   }
434
435   if (fSigmaVertex)
436   {
437     new TCanvas;
438     fSigmaVertex->Draw();
439   }
440 }