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