]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaAnalysisESDSelector.cxx
Redundant lines removed.
[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 #include <TH1F.h>
12
13 #include <AliLog.h>
14 #include <AliESDVertex.h>
15 #include <AliESD.h>
16
17 #include "esdTrackCuts/AliESDtrackCuts.h"
18 #include "dNdEta/dNdEtaAnalysis.h"
19 #include "AliPWG0Helper.h"
20
21 #include "AliGenEventHeader.h"
22 #include "AliHeader.h"
23 #include "AliStack.h"
24 #include "TParticle.h"
25
26 ClassImp(AlidNdEtaAnalysisESDSelector)
27
28 AlidNdEtaAnalysisESDSelector::AlidNdEtaAnalysisESDSelector() :
29   AliSelector(),
30   fdNdEtaAnalysis(0),
31   fMult(0),
32   fEsdTrackCuts(0)
33 {
34   //
35   // Constructor. Initialization of pointers
36   //
37
38   AliLog::SetClassDebugLevel("AlidNdEtaAnalysisESDSelector", AliLog::kDebug);
39 }
40
41 AlidNdEtaAnalysisESDSelector::~AlidNdEtaAnalysisESDSelector()
42 {
43   //
44   // Destructor
45   //
46
47   // histograms are in the output list and deleted when the output
48   // list is deleted by the TSelector dtor
49 }
50
51 void AlidNdEtaAnalysisESDSelector::Begin(TTree* tree)
52 {
53   // Begin function
54
55   ReadUserObjects(tree);
56 }
57
58 void AlidNdEtaAnalysisESDSelector::ReadUserObjects(TTree* tree)
59 {
60   // read the user objects, called from slavebegin and begin
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 void AlidNdEtaAnalysisESDSelector::SlaveBegin(TTree* tree)
73 {
74   // The SlaveBegin() function is called after the Begin() function.
75   // When running with PROOF SlaveBegin() is called on each slave server.
76   // The tree argument is deprecated (on PROOF 0 is passed).
77
78   AliSelector::SlaveBegin(tree);
79
80   ReadUserObjects(tree);
81
82   fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
83   fMult = new TH1F("fMult", "fMult;Ntracks;Count", 201, -0.5, 200.5);
84 }
85
86 void AlidNdEtaAnalysisESDSelector::Init(TTree* tree)
87 {
88   // read the user objects
89
90   AliSelector::Init(tree);
91
92   // Enable only the needed branches
93   if (tree)
94   {
95     tree->SetBranchStatus("*", 0);
96     tree->SetBranchStatus("fTriggerMask", 1);
97     tree->SetBranchStatus("fSPDVertex*", 1);
98     tree->SetBranchStatus("fTracks.fLabel", 1);
99
100     AliESDtrackCuts::EnableNeededBranches(tree);
101   }
102 }
103
104 Bool_t AlidNdEtaAnalysisESDSelector::Process(Long64_t entry)
105 {
106   // loop over all events
107
108   if (AliSelector::Process(entry) == kFALSE)
109     return kFALSE;
110
111   // Check prerequisites
112   if (!fESD)
113   {
114     AliDebug(AliLog::kError, "ESD branch not available");
115     return kFALSE;
116   }
117
118   if (!fEsdTrackCuts)
119   {
120     AliDebug(AliLog::kError, "fESDTrackCuts not available");
121     return kFALSE;
122   }
123
124   if (AliPWG0Helper::IsEventTriggered(fESD) == kFALSE)
125   {
126     AliDebug(AliLog::kDebug+1, Form("Skipping event %d because it was not triggered", (Int_t) entry));
127     return kTRUE;
128   }
129
130   if (AliPWG0Helper::IsVertexReconstructed(fESD) == kFALSE)
131   {
132     AliDebug(AliLog::kDebug+1, Form("Skipping event %d because its vertex was not reconstructed", (Int_t) entry));
133     return kTRUE;
134   }
135
136   // get the ESD vertex
137   const AliESDVertex* vtxESD = fESD->GetVertex();
138   Double_t vtx[3];
139   vtxESD->GetXYZ(vtx);
140
141   // get number of "good" tracks
142   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
143   Int_t nGoodTracks = list->GetEntries();
144
145   // loop over esd tracks
146   for (Int_t t=0; t<nGoodTracks; t++)
147   {
148     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(t));
149     if (!esdTrack)
150     {
151       AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
152       continue;
153     }
154
155     Double_t p[3];
156     esdTrack->GetConstrainedPxPyPz(p); // ### TODO should be okay because we have a vertex, however GetInnerPxPyPy / GetOuterPxPyPy also exist
157     TVector3 vector(p);
158
159     Float_t theta = vector.Theta();
160     Float_t eta   = -TMath::Log(TMath::Tan(theta/2.));
161     Float_t pt = vector.Pt();
162
163     fdNdEtaAnalysis->FillTrack(vtx[2], eta, pt);
164   } // end of track loop
165
166   delete list;
167   list = 0;
168
169   // for event count per vertex
170   fdNdEtaAnalysis->FillEvent(vtx[2], nGoodTracks);
171   fMult->Fill(nGoodTracks);
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   // Add the objects to the output list and set them to 0, so that the destructor does not delete them.
192
193   fOutput->Add(fdNdEtaAnalysis);
194   fOutput->Add(fMult);
195
196   fdNdEtaAnalysis = 0;
197 }
198
199 void AlidNdEtaAnalysisESDSelector::Terminate()
200 {
201   // The Terminate() function is the last function to be called during
202   // a query. It always runs on the client, it can be used to present
203   // the results graphically or save the results to file.
204
205   AliSelector::Terminate();
206
207   fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
208   fMult = dynamic_cast<TH1F*> (fOutput->FindObject("fMult"));
209
210   if (!fdNdEtaAnalysis)
211   {
212     AliDebug(AliLog::kError, "ERROR: Histograms not available");
213     return;
214   }
215
216   TFile* fout = new TFile("analysis_esd_raw.root", "RECREATE");
217
218   if (fdNdEtaAnalysis)
219     fdNdEtaAnalysis->SaveHistograms();
220
221   if (fEsdTrackCuts)
222     fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");
223
224   if (fMult)
225     fMult->Write();
226
227   fout->Write();
228   fout->Close();
229 }