1 #include "AlidNdEtaCorrectionSelector.h"
7 #include <TParticlePDG.h>
10 #include <AliGenEventHeader.h>
11 #include <AliTracker.h>
13 #include "../esdTrackCuts/AliESDtrackCuts.h"
14 #include "dNdEtaCorrection.h"
16 ClassImp(AlidNdEtaCorrectionSelector)
18 AlidNdEtaCorrectionSelector::AlidNdEtaCorrectionSelector(TTree *) :
22 fdNdEtaCorrectionFinal(0)
25 // Constructor. Initialization of pointers
29 AlidNdEtaCorrectionSelector::~AlidNdEtaCorrectionSelector()
35 // histograms are in the output list and deleted when the output
36 // list is deleted by the TSelector dtor
39 void AlidNdEtaCorrectionSelector::Begin(TTree * tree)
41 // The Begin() function is called at the start of the query.
42 // When running with PROOF Begin() is only called on the client.
43 // The tree argument is deprecated (on PROOF 0 is passed).
45 AliSelector::Begin(tree);
48 void AlidNdEtaCorrectionSelector::SlaveBegin(TTree * tree)
50 // The SlaveBegin() function is called after the Begin() function.
51 // When running with PROOF SlaveBegin() is called on each slave server.
52 // The tree argument is deprecated (on PROOF 0 is passed).
54 AliSelector::SlaveBegin(tree);
56 fdNdEtaCorrection = new dNdEtaCorrection();
59 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fChain->GetUserInfo()->FindObject("AliESDtrackCuts"));
62 printf("ERROR: Could not read EsdTrackCuts from user info\n");
64 AliLog::SetClassDebugLevel("ESDtrackQualityCuts",1);
67 Bool_t AlidNdEtaCorrectionSelector::Notify()
69 // The Notify() function is called when a new file is opened. This
70 // can be either for a new TTree in a TChain or when when a new TTree
71 // is started when using PROOF. Typically here the branch pointers
72 // will be retrieved. It is normaly not necessary to make changes
73 // to the generated code, but the routine can be extended by the
76 if (AliSelector::Notify() == kFALSE)
82 Bool_t AlidNdEtaCorrectionSelector::IsPrimary(const TParticle* aParticle, Int_t aTotalPrimaries)
85 // Returns if the given particle is a primary particle
86 // This function or a equivalent should be available in some common place of AliRoot
89 // if the particle has a daughter primary, we do not want to count it
90 if (aParticle->GetFirstDaughter() != -1 && aParticle->GetFirstDaughter() < aTotalPrimaries)
93 Int_t pdgCode = TMath::Abs(aParticle->GetPdgCode());
95 // skip quarks and gluon
96 if (pdgCode > 10 && pdgCode != 21)
102 Bool_t AlidNdEtaCorrectionSelector::Process(Long64_t entry)
104 // The Process() function is called for each entry in the tree (or possibly
105 // keyed object in the case of PROOF) to be processed. The entry argument
106 // specifies which entry in the currently loaded tree is to be processed.
107 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
108 // to read either all or the required parts of the data. When processing
109 // keyed objects with PROOF, the object is already loaded and is available
110 // via the fObject pointer.
112 // This function should contain the "body" of the analysis. It can contain
113 // simple or elaborate selection criteria, run algorithms on the data
114 // of the event and typically fill histograms.
116 // WARNING when a selector is used with a TChain, you must use
117 // the pointer to the current TTree to call GetEntry(entry).
118 // The entry is always the local entry number in the current tree.
119 // Assuming that fChain is the pointer to the TChain being processed,
120 // use fChain->GetTree()->GetEntry(entry).
122 if (AliSelector::Process(entry) == kFALSE)
125 if (!fESD || !fHeader)
128 // ########################################################
129 // get the EDS vertex
130 const AliESDVertex* vtxESD = fESD->GetVertex();
132 // the vertex should be reconstructed
133 if (strcmp(vtxESD->GetName(),"default")==0)
137 vtx_res[0] = vtxESD->GetXRes();
138 vtx_res[1] = vtxESD->GetYRes();
139 vtx_res[2] = vtxESD->GetZRes();
141 // the resolution should be reasonable???
142 if (vtx_res[2]==0 || vtx_res[2]>0.1)
145 // ########################################################
147 AliGenEventHeader* genHeader = fHeader->GenEventHeader();
150 genHeader->PrimaryVertex(vtxMC);
152 // ########################################################
153 // loop over mc particles
154 TTree* particleTree = GetKinematics();
155 TParticle* particle = 0;
156 particleTree->SetBranchAddress("Particles", &particle);
158 Int_t nPrim = fHeader->GetNprimary();
159 Int_t nTotal = fHeader->GetNtrack();
161 for (Int_t i_mc = nTotal - nPrim; i_mc < nTotal; ++i_mc)
163 particleTree->GetEntry(i_mc);
168 if (strcmp(particle->GetName(),"XXX") == 0)
170 printf("WARNING: There is a particle named XXX (%d).\n", i_mc);
174 TParticlePDG* pdgPart = particle->GetPDG();
176 if (strcmp(pdgPart->ParticleClass(),"Unknown") == 0)
178 printf("WARNING: There is a particle with an unknown particle class (%d pdg code %d).\n", i_mc, particle->GetPdgCode());
182 if (IsPrimary(particle, nPrim) == kFALSE)
185 if (pdgPart->Charge() == 0)
188 fdNdEtaCorrection->FillGene(vtxMC[2], particle->Eta());
190 }// end of mc particle
192 // ########################################################
193 // loop over esd tracks
194 Int_t nTracks = fESD->GetNumberOfTracks();
195 for (Int_t t=0; t<nTracks; t++)
197 AliESDtrack* esdTrack = fESD->GetTrack(t);
199 // cut the esd track?
200 if (!fEsdTrackCuts->AcceptTrack(esdTrack))
203 // using the eta of the mc particle
204 Int_t label = TMath::Abs(esdTrack->GetLabel());
207 printf("WARNING: cannot find corresponding mc part for track %d.", t);
210 particleTree->GetEntry(nTotal - nPrim + label);
212 fdNdEtaCorrection->FillMeas(vtxMC[2], particle->Eta());
214 } // end of track loop
219 void AlidNdEtaCorrectionSelector::SlaveTerminate()
221 // The SlaveTerminate() function is called after all entries or objects
222 // have been processed. When running with PROOF SlaveTerminate() is called
223 // on each slave server.
225 AliSelector::SlaveTerminate();
227 // Add the histograms to the output on each slave server
230 printf("ERROR: Output list not initialized\n");
234 fOutput->Add(fdNdEtaCorrection->GetGeneratedHistogram());
235 fOutput->Add(fdNdEtaCorrection->GetMeasuredHistogram());
238 void AlidNdEtaCorrectionSelector::Terminate()
240 // The Terminate() function is the last function to be called during
241 // a query. It always runs on the client, it can be used to present
242 // the results graphically or save the results to file.
244 AliSelector::Terminate();
246 fdNdEtaCorrectionFinal = new dNdEtaCorrection();
247 TH2F* measuredHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("etaVsVtx_meas"));
248 TH2F* generatedHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("etaVsVtx_gene"));
249 if (!measuredHistogram || !generatedHistogram)
251 printf("ERROR: Histograms not available %p %p\n", (void*) generatedHistogram, (void*) measuredHistogram);
254 fdNdEtaCorrectionFinal->SetGeneratedHistogram(generatedHistogram);
255 fdNdEtaCorrectionFinal->SetMeasuredHistogram(measuredHistogram);
257 fdNdEtaCorrectionFinal->Finish();
259 TFile* fout = new TFile("correction_map.root","RECREATE");
261 fEsdTrackCuts->SaveHistograms("esd_track_cuts");
262 fdNdEtaCorrectionFinal->SaveHistograms();
267 fdNdEtaCorrectionFinal->DrawHistograms();