]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaAnalysisESDSelector.cxx
implemented analysis using 3d information (vtx_z, eta, pt)
[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::IsEventTriggered(fESD) == kFALSE)
127     return kTRUE;
128
129   if (AliPWG0Helper::IsVertexReconstructed(fESD) == kFALSE)
130     return kTRUE;
131
132   // ########################################################
133   // get the EDS vertex
134   const AliESDVertex* vtxESD = fESD->GetVertex();
135   Double_t vtx[3];
136   vtxESD->GetXYZ(vtx);
137
138   // get number of "good" tracks
139   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
140   Int_t nGoodTracks = list->GetEntries();
141
142   Float_t vertexRecoCorr = fdNdEtaCorrection->GetVertexRecoCorrection(vtx[2], nGoodTracks);
143   if (vertexRecoCorr <= 0)
144   {
145     AliDebug(AliLog::kError, Form("INFO: Skipping event because vertexRecoCorr is <= 0 (%f)", vertexRecoCorr));
146     delete list;
147     return kTRUE;
148   }
149
150   // loop over esd tracks
151   for (Int_t t=0; t<nGoodTracks; t++)
152   {
153     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(t));
154     if (!esdTrack)
155     {
156       AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
157       continue;
158     }
159
160     Double_t p[3];
161     esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
162     TVector3 vector(p);
163
164     Float_t theta = vector.Theta();
165     Float_t eta   = -TMath::Log(TMath::Tan(theta/2.));
166     Float_t pt = vector.Pt();
167
168     Float_t track2particleCorr = fdNdEtaCorrection->GetTrack2ParticleCorrection(vtx[2], eta, pt);
169
170     Float_t weight = vertexRecoCorr * track2particleCorr;
171     if (weight <= 0)
172     {
173       AliDebug(AliLog::kError, Form("INFO: Skipping track because weight is <= 0 (track %d, weight %f)", t, weight));
174       continue;
175     }
176
177     fdNdEtaAnalysis->FillTrack(vtx[2], eta, pt, weight);
178   } // end of track loop
179
180   delete list;
181   list = 0;
182
183   // for event count per vertex
184   fdNdEtaAnalysis->FillEvent(vtx[2], vertexRecoCorr);
185
186   return kTRUE;
187 }
188
189 void AlidNdEtaAnalysisESDSelector::SlaveTerminate()
190 {
191   // The SlaveTerminate() function is called after all entries or objects
192   // have been processed. When running with PROOF SlaveTerminate() is called
193   // on each slave server.
194
195   AliSelector::SlaveTerminate();
196
197   // Add the histograms to the output on each slave server
198   if (!fOutput)
199   {
200     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
201     return;
202   }
203
204   fOutput->Add(fdNdEtaAnalysis);
205 }
206
207 void AlidNdEtaAnalysisESDSelector::Terminate()
208 {
209   // The Terminate() function is the last function to be called during
210   // a query. It always runs on the client, it can be used to present
211   // the results graphically or save the results to file.
212
213   AliSelector::Terminate();
214
215   fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
216
217   if (!fdNdEtaAnalysis)
218   {
219     AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fdNdEtaAnalysis));
220     return;
221   }
222
223   TFile* fout = new TFile("analysis_esd.root","RECREATE");
224
225   if (fdNdEtaAnalysis)
226     fdNdEtaAnalysis->SaveHistograms();
227
228   if (fEsdTrackCuts)
229     fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");
230
231   fout->Write();
232   fout->Close();
233 }