3 #include "AlidNdEtaCorrectionSelector.h"
9 #include <TParticlePDG.h>
13 #include <TSelector.h>
17 #include <AliTracker.h>
18 #include <AliESDVertex.h>
20 #include <AliESDtrack.h>
21 #include <AliRunLoader.h>
24 #include <AliHeader.h>
25 #include <AliGenEventHeader.h>
26 #include <../PYTHIA6/AliGenPythiaEventHeader.h>
27 #include <../EVGEN/AliGenCocktailEventHeader.h>
29 #include "esdTrackCuts/AliESDtrackCuts.h"
30 #include "dNdEta/AlidNdEtaCorrection.h"
31 #include "AliPWG0Helper.h"
32 #include "AliPWG0depHelper.h"
34 ClassImp(AlidNdEtaCorrectionSelector)
36 AlidNdEtaCorrectionSelector::AlidNdEtaCorrectionSelector() :
49 // Constructor. Initialization of pointers
53 AlidNdEtaCorrectionSelector::~AlidNdEtaCorrectionSelector()
59 // histograms are in the output list and deleted when the output
60 // list is deleted by the TSelector dtor
63 Bool_t AlidNdEtaCorrectionSelector::SignOK(TParticlePDG* particle)
65 // returns if a particle with this sign should be counted
66 // this is determined by the value of fSignMode, which should have the same sign
68 // if fSignMode is 0 all particles are counted
75 printf("WARNING: not counting a particle that does not have a pdg particle\n");
79 Double_t charge = particle->Charge();
92 void AlidNdEtaCorrectionSelector::ReadUserObjects(TTree* tree)
94 // read the user objects, called from slavebegin and begin
96 if (!fEsdTrackCuts && fInput)
97 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
99 if (!fEsdTrackCuts && tree)
100 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
103 AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
106 void AlidNdEtaCorrectionSelector::Begin(TTree * tree)
108 // The Begin() function is called at the start of the query.
109 // When running with PROOF Begin() is only called on the client.
110 // The tree argument is deprecated (on PROOF 0 is passed).
112 AliSelectorRL::Begin(tree);
114 ReadUserObjects(tree);
116 TString option = GetOption();
117 AliInfo(Form("Called with option %s.", option.Data()));
119 if (option.Contains("only-positive"))
121 AliInfo("Processing only positive particles.");
124 else if (option.Contains("only-negative"))
126 AliInfo("Processing only negative particles.");
131 void AlidNdEtaCorrectionSelector::SlaveBegin(TTree * tree)
133 // The SlaveBegin() function is called after the Begin() function.
134 // When running with PROOF SlaveBegin() is called on each slave server.
135 // The tree argument is deprecated (on PROOF 0 is passed).
137 AliSelectorRL::SlaveBegin(tree);
139 ReadUserObjects(tree);
141 fdNdEtaCorrection = new AlidNdEtaCorrection("dndeta_correction", "dndeta_correction");
143 fPIDParticles = new TH1F("pid_particles", "PID of generated primary particles", 10001, -5000.5, 5000.5);
144 fPIDTracks = new TH1F("pid_tracks", "MC PID of reconstructed tracks", 10001, -5000.5, 5000.5);
146 fClustersITSPos = new TH1F("clusters_its_pos", "clusters_its_pos", 7, -0.5, 6.5);
147 fClustersTPCPos = new TH1F("clusters_tpc_pos", "clusters_tpc_pos", 160, -0.5, 159.5);
149 fClustersITSNeg = new TH1F("clusters_its_neg", "clusters_its_neg", 7, -0.5, 6.5);
150 fClustersTPCNeg = new TH1F("clusters_tpc_neg", "clusters_tpc_neg", 160, -0.5, 159.5);
153 void AlidNdEtaCorrectionSelector::Init(TTree* tree)
155 // read the user objects
157 AliSelectorRL::Init(tree);
159 // Enable only the needed branches
162 tree->SetBranchStatus("*", 0);
163 tree->SetBranchStatus("fTriggerMask", 1);
164 tree->SetBranchStatus("fSPDVertex*", 1);
165 tree->SetBranchStatus("fTracks.fLabel", 1);
166 tree->SetBranchStatus("fTracks.fITSncls", 1);
167 tree->SetBranchStatus("fTracks.fTPCncls", 1);
169 AliESDtrackCuts::EnableNeededBranches(tree);
173 Bool_t AlidNdEtaCorrectionSelector::Process(Long64_t entry)
175 // The Process() function is called for each entry in the tree (or possibly
176 // keyed object in the case of PROOF) to be processed. The entry argument
177 // specifies which entry in the currently loaded tree is to be processed.
178 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
179 // to read either all or the required parts of the data. When processing
180 // keyed objects with PROOF, the object is already loaded and is available
181 // via the fObject pointer.
183 // This function should contain the "body" of the analysis. It can contain
184 // simple or elaborate selection criteria, run algorithms on the data
185 // of the event and typically fill histograms.
187 // WARNING when a selector is used with a TChain, you must use
188 // the pointer to the current TTree to call GetEntry(entry).
189 // The entry is always the local entry number in the current tree.
190 // Assuming that fTree is the pointer to the TChain being processed,
191 // use fTree->GetTree()->GetEntry(entry).
193 AliDebug(AliLog::kDebug+1,"Processing event ...\n");
196 if (AliSelectorRL::Process(entry) == kFALSE)
199 // check prerequesites
202 AliDebug(AliLog::kError, "ESD branch not available");
206 AliHeader* header = GetHeader();
209 AliDebug(AliLog::kError, "Header not available");
213 AliStack* stack = GetStack();
216 AliDebug(AliLog::kError, "Stack not available");
222 AliDebug(AliLog::kError, "fESDTrackCuts not available");
226 Bool_t vertexReconstructed = AliPWG0Helper::IsVertexReconstructed(fESD);
228 Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD);
231 AliGenEventHeader* genHeader = header->GenEventHeader();
234 genHeader->PrimaryVertex(vtxMC);
236 // loop over mc particles
237 Int_t nPrim = stack->GetNprimary();
239 for (Int_t iMc = 0; iMc < nPrim; ++iMc)
241 AliDebug(AliLog::kDebug+1, Form("MC Loop: Processing particle %d.", iMc));
243 TParticle* particle = stack->Particle(iMc);
247 AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
251 if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
254 if (SignOK(particle->GetPDG()) == kFALSE)
257 Float_t eta = particle->Eta();
258 Float_t pt = particle->Pt();
260 if (vertexReconstructed) {
261 fdNdEtaCorrection->FillParticle(vtxMC[2], eta, pt);
263 if (pt > 0.1 && pt < 0.2)
264 fPIDParticles->Fill(particle->GetPdgCode());
266 }// end of mc particle
268 // ########################################################
269 // loop over esd tracks
270 Int_t nTracks = fESD->GetNumberOfTracks();
272 // count the number of "good" tracks as parameter for vertex reconstruction efficiency
273 Int_t nGoodTracks = 0;
274 for (Int_t t=0; t<nTracks; t++)
276 AliDebug(AliLog::kDebug+1, Form("ESD Loop: Processing track %d.", t));
278 AliESDtrack* esdTrack = fESD->GetTrack(t);
280 // cut the esd track?
281 if (!fEsdTrackCuts->AcceptTrack(esdTrack))
286 // using the properties of the mc particle
287 Int_t label = TMath::Abs(esdTrack->GetLabel());
290 AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", t));
294 TParticle* particle = stack->Particle(label);
297 AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
301 if (SignOK(particle->GetPDG()) == kFALSE)
304 if (vertexReconstructed)
306 fdNdEtaCorrection->FillParticleWhenMeasuredTrack(vtxMC[2], particle->Eta(), particle->Pt());
307 if (particle->Pt() > 0.1 && particle->Pt() < 0.2)
309 fPIDTracks->Fill(particle->GetPdgCode());
310 if (particle->GetPDG()->Charge() > 0)
312 fClustersITSPos->Fill(esdTrack->GetITSclusters(0));
313 fClustersTPCPos->Fill(esdTrack->GetTPCclusters(0));
317 fClustersITSNeg->Fill(esdTrack->GetITSclusters(0));
318 fClustersTPCNeg->Fill(esdTrack->GetTPCclusters(0));
322 } // end of track loop
324 // stuff regarding the vertex reco correction and trigger bias correction
325 if (eventTriggered) {
326 fdNdEtaCorrection->FillEventWithTrigger(vtxMC[2], nGoodTracks);
327 if (vertexReconstructed)
328 fdNdEtaCorrection->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
331 // getting process information
332 Int_t processtype = AliPWG0depHelper::GetPythiaEventProcessType(header);
333 AliDebug(AliLog::kDebug+1,Form(" Found pythia procces type %d", processtype));
336 AliDebug(AliLog::kError, Form("Unkown Pythia process type %d.", processtype));
338 fdNdEtaCorrection->FillEventAll(vtxMC[2], nGoodTracks, "INEL");
340 if (processtype!=92 && processtype!=93)
341 fdNdEtaCorrection->FillEventAll(vtxMC[2], nGoodTracks, "NSD");
343 if (processtype!=92 && processtype!=93 && processtype!=94)
344 fdNdEtaCorrection->FillEventAll(vtxMC[2], nGoodTracks, "ND");
349 void AlidNdEtaCorrectionSelector::SlaveTerminate()
351 // The SlaveTerminate() function is called after all entries or objects
352 // have been processed. When running with PROOF SlaveTerminate() is called
353 // on each slave server.
355 AliSelectorRL::SlaveTerminate();
357 // Add the histograms to the output on each slave server
360 AliDebug(AliLog::kError, "ERROR: Output list not initialized");
364 fOutput->Add(fdNdEtaCorrection);
367 void AlidNdEtaCorrectionSelector::Terminate()
369 // The Terminate() function is the last function to be called during
370 // a query. It always runs on the client, it can be used to present
371 // the results graphically or save the results to file.
373 AliSelectorRL::Terminate();
375 fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fOutput->FindObject("dndeta_correction"));
376 if (!fdNdEtaCorrection)
378 AliDebug(AliLog::kError, "Could not read object from output list");
382 fdNdEtaCorrection->Finish();
384 TFile* fout = new TFile(Form("correction_map%s.root", GetOption()), "RECREATE");
387 fEsdTrackCuts->SaveHistograms("esd_track_cuts");
388 fdNdEtaCorrection->SaveHistograms();
393 fdNdEtaCorrection->DrawHistograms();
395 if (fPIDParticles && fPIDTracks)
397 new TCanvas("pidcanvas", "pidcanvas", 500, 500);
399 fPIDParticles->Draw();
400 fPIDTracks->SetLineColor(2);
401 fPIDTracks->Draw("SAME");
403 TDatabasePDG* pdgDB = new TDatabasePDG;
405 for (Int_t i=0; i <= fPIDParticles->GetNbinsX()+1; ++i)
406 if (fPIDParticles->GetBinContent(i) > 0)
407 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));
413 if (fClustersITSPos && fClustersITSNeg && fClustersTPCPos && fClustersTPCNeg)
415 TCanvas* canvas = new TCanvas("clusters", "clusters", 1000, 500);
416 canvas->Divide(2, 1);
419 fClustersITSPos->Draw();
420 fClustersITSNeg->SetLineColor(kRed);
421 fClustersITSNeg->Draw("SAME");
424 fClustersTPCPos->Draw();
425 fClustersTPCNeg->SetLineColor(kRed);
426 fClustersTPCNeg->Draw("SAME");