]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaAnalysisESDSelector.cxx
minor changes to improve serializability
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaAnalysisESDSelector.cxx
1 /* $Id$ */
2
3 #include "AlidNdEtaAnalysisESDSelector.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
12 #include <AliLog.h>
13 #include <AliESDVertex.h>
14 #include <AliESD.h>
15
16 #include "esdTrackCuts/AliESDtrackCuts.h"
17 #include "dNdEta/dNdEtaAnalysis.h"
18 #include "dNdEta/AlidNdEtaCorrection.h"
19 #include "AliPWG0Helper.h"
20
21 ClassImp(AlidNdEtaAnalysisESDSelector)
22
23 AlidNdEtaAnalysisESDSelector::AlidNdEtaAnalysisESDSelector() :
24   AliSelector(),
25   fdNdEtaAnalysisMBVtx(0),
26   fdNdEtaAnalysisMB(0),
27   fdNdEtaAnalysis(0),
28   fEsdTrackCuts(0),
29   fdNdEtaCorrection(0)
30 {
31   //
32   // Constructor. Initialization of pointers
33   //
34
35   //AliLog::SetClassDebugLevel("AlidNdEtaAnalysisESDSelector", AliLog::kDebug);
36 }
37
38 AlidNdEtaAnalysisESDSelector::~AlidNdEtaAnalysisESDSelector()
39 {
40   //
41   // Destructor
42   //
43
44   // histograms are in the output list and deleted when the output
45   // list is deleted by the TSelector dtor
46 }
47
48 void AlidNdEtaAnalysisESDSelector::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   if (fInput)
57   {
58     printf("Printing input list:\n");
59     fInput->Print();
60   }
61
62   if (!fEsdTrackCuts && fInput)
63     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
64
65   if (!fEsdTrackCuts && tree)
66     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
67
68   if (!fEsdTrackCuts)
69      AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
70
71
72   if (!fdNdEtaCorrection && fInput)
73     fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fInput->FindObject("dndeta_correction"));
74
75   if (!fdNdEtaCorrection && fTree)
76     fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fTree->GetUserInfo()->FindObject("dndeta_correction"));
77
78   if (!fdNdEtaCorrection)
79      AliDebug(AliLog::kError, "ERROR: Could not read dndeta_correction from input list.");
80
81
82   fdNdEtaAnalysisMBVtx = new dNdEtaAnalysis("dndeta_mbvtx", "dndeta_mbvtx");
83   fdNdEtaAnalysisMB = new dNdEtaAnalysis("dndeta_mb", "dndeta_mb");
84   fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
85 }
86
87 void AlidNdEtaAnalysisESDSelector::Init(TTree* tree)
88 {
89   // read the user objects
90
91   AliSelector::Init(tree);
92 }
93
94 Bool_t AlidNdEtaAnalysisESDSelector::Process(Long64_t entry)
95 {
96   // The Process() function is called for each entry in the tree (or possibly
97   // keyed object in the case of PROOF) to be processed. The entry argument
98   // specifies which entry in the currently loaded tree is to be processed.
99   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
100   // to read either all or the required parts of the data. When processing
101   // keyed objects with PROOF, the object is already loaded and is available
102   // via the fObject pointer.
103   //
104   // This function should contain the "body" of the analysis. It can contain
105   // simple or elaborate selection criteria, run algorithms on the data
106   // of the event and typically fill histograms.
107
108   // WARNING when a selector is used with a TChain, you must use
109   //  the pointer to the current TTree to call GetEntry(entry).
110   //  The entry is always the local entry number in the current tree.
111   //  Assuming that fTree is the pointer to the TChain being processed,
112   //  use fTree->GetTree()->GetEntry(entry).
113
114   if (AliSelector::Process(entry) == kFALSE)
115     return kFALSE;
116
117   // Check prerequisites
118   if (!fESD)
119   {
120     AliDebug(AliLog::kError, "ESD branch not available");
121     return kFALSE;
122   }
123
124   if (!fEsdTrackCuts)
125   {
126     AliDebug(AliLog::kError, "fESDTrackCuts not available");
127     return kFALSE;
128   }
129
130   if (!fdNdEtaCorrection)
131   {
132     AliDebug(AliLog::kError, "fdNdEtaCorrection not available");
133     return kFALSE;
134   }
135
136   if (AliPWG0Helper::IsEventTriggered(fESD) == kFALSE)
137     return kTRUE;
138
139   if (AliPWG0Helper::IsVertexReconstructed(fESD) == kFALSE)
140     return kTRUE;
141
142   // ########################################################
143   // get the EDS vertex
144   const AliESDVertex* vtxESD = fESD->GetVertex();
145   Double_t vtx[3];
146   vtxESD->GetXYZ(vtx);
147
148   // get number of "good" tracks
149   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
150   Int_t nGoodTracks = list->GetEntries();
151
152   Float_t vertexRecoCorr = fdNdEtaCorrection->GetVertexRecoCorrection(vtx[2], nGoodTracks);
153   if (vertexRecoCorr <= 0)
154   {
155     AliDebug(AliLog::kError, Form("INFO: Skipping event because vertexRecoCorr is <= 0 (%f)", vertexRecoCorr));
156     delete list;
157     return kTRUE;
158   }
159
160   Float_t triggerCorr = fdNdEtaCorrection->GetTriggerCorrection(vtx[2], nGoodTracks);
161   if (triggerCorr <= 0)
162   {
163     AliDebug(AliLog::kError, Form("INFO: Skipping event because triggerCorr is <= 0 (%f)", triggerCorr));
164     delete list;
165     return kTRUE;
166   }
167
168   // loop over esd tracks
169   for (Int_t t=0; t<nGoodTracks; t++)
170   {
171     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(t));
172     if (!esdTrack)
173     {
174       AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
175       continue;
176     }
177
178     Double_t p[3];
179     esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
180     TVector3 vector(p);
181
182     Float_t theta = vector.Theta();
183     Float_t eta   = -TMath::Log(TMath::Tan(theta/2.));
184     Float_t pt = vector.Pt();
185
186     Float_t track2particleCorr = fdNdEtaCorrection->GetTrack2ParticleCorrection(vtx[2], eta, pt);
187
188     Float_t weight = track2particleCorr * vertexRecoCorr * triggerCorr;
189     if (weight <= 0)
190     {
191       AliDebug(AliLog::kError, Form("INFO: Skipping track because weight is <= 0 (track %d, weight %f) (vtx %f, eta %f, pt %f)", t, weight, vtx[2], eta, pt));
192       continue;
193     }
194
195     fdNdEtaAnalysisMBVtx->FillTrack(vtx[2], eta, pt, track2particleCorr);
196     fdNdEtaAnalysisMB->FillTrack(vtx[2], eta, pt, track2particleCorr * vertexRecoCorr);
197     fdNdEtaAnalysis->FillTrack(vtx[2], eta, pt, weight);
198   } // end of track loop
199
200   delete list;
201   list = 0;
202
203   // for event count per vertex
204   fdNdEtaAnalysisMBVtx->FillEvent(vtx[2], 1);
205   fdNdEtaAnalysisMB->FillEvent(vtx[2], vertexRecoCorr);
206   fdNdEtaAnalysis->FillEvent(vtx[2], vertexRecoCorr * triggerCorr);
207
208   return kTRUE;
209 }
210
211 void AlidNdEtaAnalysisESDSelector::SlaveTerminate()
212 {
213   // The SlaveTerminate() function is called after all entries or objects
214   // have been processed. When running with PROOF SlaveTerminate() is called
215   // on each slave server.
216
217   AliSelector::SlaveTerminate();
218
219   // Add the histograms to the output on each slave server
220   if (!fOutput)
221   {
222     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
223     return;
224   }
225
226   fOutput->Add(fdNdEtaAnalysis);
227 }
228
229 void AlidNdEtaAnalysisESDSelector::Terminate()
230 {
231   // The Terminate() function is the last function to be called during
232   // a query. It always runs on the client, it can be used to present
233   // the results graphically or save the results to file.
234
235   AliSelector::Terminate();
236
237   fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
238
239   if (!fdNdEtaAnalysis)
240   {
241     AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fdNdEtaAnalysis));
242     return;
243   }
244
245   fdNdEtaAnalysis->Finish(fdNdEtaCorrection, 0.3);
246   fdNdEtaAnalysisMB->Finish(fdNdEtaCorrection, 0.3);
247   fdNdEtaAnalysisMBVtx->Finish(fdNdEtaCorrection, 0.3);
248
249   TFile* fout = new TFile("analysis_esd.root","RECREATE");
250
251   if (fdNdEtaAnalysis)
252     fdNdEtaAnalysis->SaveHistograms();
253
254   if (fdNdEtaAnalysisMB)
255     fdNdEtaAnalysisMB->SaveHistograms();
256
257   if (fdNdEtaAnalysisMBVtx)
258     fdNdEtaAnalysisMBVtx->SaveHistograms();
259
260   if (fEsdTrackCuts)
261     fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");
262
263   fout->Write();
264   fout->Close();
265 }