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