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