]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaCorrectionSelector.cxx
o) adding log tags to all files
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaCorrectionSelector.cxx
1 /* $Id$ */
2
3 #include "AlidNdEtaCorrectionSelector.h"
4
5 #include <TStyle.h>
6 #include <TSystem.h>
7 #include <TCanvas.h>
8 #include <TParticle.h>
9 #include <TParticlePDG.h>
10
11 #include <AliLog.h>
12 #include <AliGenEventHeader.h>
13 #include <AliTracker.h>
14
15 #include "esdTrackCuts/AliESDtrackCuts.h"
16 #include "dNdEtaCorrection.h"
17
18 ClassImp(AlidNdEtaCorrectionSelector)
19
20 AlidNdEtaCorrectionSelector::AlidNdEtaCorrectionSelector(TTree *) :
21   AliSelector(),
22   fEsdTrackCuts(0),
23   fdNdEtaCorrection(0),
24   fdNdEtaCorrectionFinal(0)
25 {
26   //
27   // Constructor. Initialization of pointers
28   //
29 }
30
31 AlidNdEtaCorrectionSelector::~AlidNdEtaCorrectionSelector()
32 {
33   //
34   // Destructor
35   //
36
37   // histograms are in the output list and deleted when the output
38   // list is deleted by the TSelector dtor
39 }
40
41 void AlidNdEtaCorrectionSelector::Begin(TTree * tree)
42 {
43   // The Begin() function is called at the start of the query.
44   // When running with PROOF Begin() is only called on the client.
45   // The tree argument is deprecated (on PROOF 0 is passed).
46
47   AliSelector::Begin(tree);
48 }
49
50 void AlidNdEtaCorrectionSelector::SlaveBegin(TTree * tree)
51 {
52   // The SlaveBegin() function is called after the Begin() function.
53   // When running with PROOF SlaveBegin() is called on each slave server.
54   // The tree argument is deprecated (on PROOF 0 is passed).
55
56   AliSelector::SlaveBegin(tree);
57
58   fdNdEtaCorrection = new dNdEtaCorrection();
59
60   if (fChain)
61     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fChain->GetUserInfo()->FindObject("AliESDtrackCuts"));
62
63   if (!fEsdTrackCuts)
64     AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from user info");
65 }
66
67 Bool_t AlidNdEtaCorrectionSelector::Notify()
68 {
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
74   // user if needed.
75
76   if (AliSelector::Notify() == kFALSE)
77     return kFALSE;
78
79   return kTRUE;
80 }
81
82 Bool_t AlidNdEtaCorrectionSelector::Process(Long64_t entry)
83 {
84   // The Process() function is called for each entry in the tree (or possibly
85   // keyed object in the case of PROOF) to be processed. The entry argument
86   // specifies which entry in the currently loaded tree is to be processed.
87   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
88   // to read either all or the required parts of the data. When processing
89   // keyed objects with PROOF, the object is already loaded and is available
90   // via the fObject pointer.
91   //
92   // This function should contain the "body" of the analysis. It can contain
93   // simple or elaborate selection criteria, run algorithms on the data
94   // of the event and typically fill histograms.
95
96   // WARNING when a selector is used with a TChain, you must use
97   //  the pointer to the current TTree to call GetEntry(entry).
98   //  The entry is always the local entry number in the current tree.
99   //  Assuming that fChain is the pointer to the TChain being processed,
100   //  use fChain->GetTree()->GetEntry(entry).
101
102   if (AliSelector::Process(entry) == kFALSE)
103     return kFALSE;
104
105   // check prerequesites
106   if (!fESD)
107   {
108     AliDebug(AliLog::kError, "ESD branch not available");
109     return kFALSE;
110   }
111
112   AliHeader* header = GetHeader();
113   if (!header)
114   {
115     AliDebug(AliLog::kError, "Header not available");
116     return kFALSE;
117   }
118
119   if (!fEsdTrackCuts)
120   {
121     AliDebug(AliLog::kError, "fESDTrackCuts not available");
122     return kFALSE;
123   }
124
125   // ########################################################
126   // get the EDS vertex
127   const AliESDVertex* vtxESD = fESD->GetVertex();
128
129   // the vertex should be reconstructed
130   if (strcmp(vtxESD->GetName(),"default")==0)
131     return kTRUE;
132
133   Double_t vtx_res[3];
134   vtx_res[0] = vtxESD->GetXRes();
135   vtx_res[1] = vtxESD->GetYRes();
136   vtx_res[2] = vtxESD->GetZRes();
137
138   // the resolution should be reasonable???
139   if (vtx_res[2]==0 || vtx_res[2]>0.1)
140     return kTRUE;
141
142   // ########################################################
143   // get the MC vertex
144   AliGenEventHeader* genHeader = header->GenEventHeader();
145
146   TArrayF vtxMC(3);
147   genHeader->PrimaryVertex(vtxMC);
148
149   // ########################################################
150   // loop over mc particles
151   TTree* particleTree = GetKinematics();
152   TParticle* particle = 0;
153   particleTree->SetBranchAddress("Particles", &particle);
154
155   Int_t nPrim  = header->GetNprimary();
156   Int_t nTotal = header->GetNtrack();
157
158   for (Int_t i_mc = nTotal - nPrim; i_mc < nTotal; ++i_mc)
159   {
160     particleTree->GetEntry(i_mc);
161
162     if (!particle)
163       continue;
164
165     if (IsPrimaryCharged(particle, nPrim) == kFALSE)
166       continue;
167
168     fdNdEtaCorrection->FillGene(vtxMC[2], particle->Eta());
169   }// end of mc particle
170
171   // ########################################################
172   // loop over esd tracks
173   Int_t nTracks = fESD->GetNumberOfTracks();
174   for (Int_t t=0; t<nTracks; t++)
175   {
176     AliESDtrack* esdTrack = fESD->GetTrack(t);
177
178     // cut the esd track?
179     if (!fEsdTrackCuts->AcceptTrack(esdTrack))
180       continue;
181
182     // using the eta of the mc particle
183     Int_t label = TMath::Abs(esdTrack->GetLabel());
184     if (label == 0)
185     {
186       AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", t));
187       continue;
188     }
189     particleTree->GetEntry(nTotal - nPrim + label);
190
191     fdNdEtaCorrection->FillMeas(vtxMC[2], particle->Eta());
192
193   } // end of track loop
194
195   return kTRUE;
196 }
197
198 void AlidNdEtaCorrectionSelector::SlaveTerminate()
199 {
200   // The SlaveTerminate() function is called after all entries or objects
201   // have been processed. When running with PROOF SlaveTerminate() is called
202   // on each slave server.
203
204   AliSelector::SlaveTerminate();
205
206   // Add the histograms to the output on each slave server
207   if (!fOutput)
208   {
209     AliDebug(AliLog::kError, "ERROR: Output list not initialized");
210     return;
211   }
212
213   fOutput->Add(fdNdEtaCorrection->GetGeneratedHistogram());
214   fOutput->Add(fdNdEtaCorrection->GetMeasuredHistogram());
215 }
216
217 void AlidNdEtaCorrectionSelector::Terminate()
218 {
219   // The Terminate() function is the last function to be called during
220   // a query. It always runs on the client, it can be used to present
221   // the results graphically or save the results to file.
222
223   AliSelector::Terminate();
224
225   fdNdEtaCorrectionFinal = new dNdEtaCorrection();
226   TH2F* measuredHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("etaVsVtx_meas"));
227   TH2F* generatedHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("etaVsVtx_gene"));
228   if (!measuredHistogram || !generatedHistogram)
229   {
230      AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p %p", (void*) generatedHistogram, (void*) measuredHistogram));
231     return;
232   }
233   fdNdEtaCorrectionFinal->SetGeneratedHistogram(generatedHistogram);
234   fdNdEtaCorrectionFinal->SetMeasuredHistogram(measuredHistogram);
235
236   fdNdEtaCorrectionFinal->Finish();
237
238   TFile* fout = new TFile("correction_map.root","RECREATE");
239
240   fEsdTrackCuts->SaveHistograms("esd_track_cuts");
241   fdNdEtaCorrectionFinal->SaveHistograms();
242
243   fout->Write();
244   fout->Close();
245
246   fdNdEtaCorrectionFinal->DrawHistograms();
247 }