]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaAnalysisESDSelector.cxx
o) adding log tags to all files
[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
10 #include <AliLog.h>
11 #include <AliGenEventHeader.h>
12
13 #include "esdTrackCuts/AliESDtrackCuts.h"
14 #include "dNdEtaCorrection.h"
15 #include "dNdEtaAnalysis.h"
16
17 ClassImp(AlidNdEtaAnalysisESDSelector)
18
19 AlidNdEtaAnalysisESDSelector::AlidNdEtaAnalysisESDSelector() :
20   AlidNdEtaAnalysisSelector(),
21   fEsdTrackCuts(0),
22   fdNdEtaCorrection(0)
23 {
24   //
25   // Constructor. Initialization of pointers
26   //
27 }
28
29 AlidNdEtaAnalysisESDSelector::~AlidNdEtaAnalysisESDSelector()
30 {
31   //
32   // Destructor
33   //
34
35   // histograms are in the output list and deleted when the output
36   // list is deleted by the TSelector dtor
37 }
38
39 void AlidNdEtaAnalysisESDSelector::SlaveBegin(TTree * tree)
40 {
41   // The SlaveBegin() function is called after the Begin() function.
42   // When running with PROOF SlaveBegin() is called on each slave server.
43   // The tree argument is deprecated (on PROOF 0 is passed).
44
45   AlidNdEtaAnalysisSelector::SlaveBegin(tree);
46
47   if (fChain)
48   {
49     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fChain->GetUserInfo()->FindObject("AliESDtrackCuts"));
50     fdNdEtaCorrection = dynamic_cast<dNdEtaCorrection*> (fChain->GetUserInfo()->FindObject("dNdEtaCorrection"));
51   }
52
53   if (!fEsdTrackCuts)
54      AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from user info.");
55
56   if (!fEsdTrackCuts)
57      AliDebug(AliLog::kError, "ERROR: Could not read dNdEtaCorrection from user info.");
58 }
59
60 Bool_t AlidNdEtaAnalysisESDSelector::Process(Long64_t entry)
61 {
62   // The Process() function is called for each entry in the tree (or possibly
63   // keyed object in the case of PROOF) to be processed. The entry argument
64   // specifies which entry in the currently loaded tree is to be processed.
65   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
66   // to read either all or the required parts of the data. When processing
67   // keyed objects with PROOF, the object is already loaded and is available
68   // via the fObject pointer.
69   //
70   // This function should contain the "body" of the analysis. It can contain
71   // simple or elaborate selection criteria, run algorithms on the data
72   // of the event and typically fill histograms.
73
74   // WARNING when a selector is used with a TChain, you must use
75   //  the pointer to the current TTree to call GetEntry(entry).
76   //  The entry is always the local entry number in the current tree.
77   //  Assuming that fChain is the pointer to the TChain being processed,
78   //  use fChain->GetTree()->GetEntry(entry).
79
80   if (AlidNdEtaAnalysisSelector::Process(entry) == kFALSE)
81     return kFALSE;
82
83   // Check prerequisites
84   if (!fESD)
85   {
86     AliDebug(AliLog::kError, "ESD branch not available");
87     return kFALSE;
88   }
89
90   if (!fEsdTrackCuts)
91   {
92     AliDebug(AliLog::kError, "fESDTrackCuts not available");
93     return kFALSE;
94   }
95
96   // ########################################################
97   // get the EDS vertex
98   const AliESDVertex* vtxESD = fESD->GetVertex();
99
100   // the vertex should be reconstructed
101   if (strcmp(vtxESD->GetName(),"default")==0)
102     return kTRUE;
103
104   Double_t vtx_res[3];
105   vtx_res[0] = vtxESD->GetXRes();
106   vtx_res[1] = vtxESD->GetYRes();
107   vtx_res[2] = vtxESD->GetZRes();
108
109   // the resolution should be reasonable???
110   if (vtx_res[2]==0 || vtx_res[2]>0.1)
111     return kTRUE;
112
113   Double_t vtx[3];
114   vtxESD->GetXYZ(vtx);
115
116   // ########################################################
117   // loop over esd tracks
118   Int_t nTracks = fESD->GetNumberOfTracks();
119   for (Int_t t=0; t<nTracks; t++)
120   {
121     AliESDtrack* esdTrack = fESD->GetTrack(t);
122     if (!esdTrack)
123     {
124       AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
125       continue;
126     }
127
128     // cut the esd track?
129     if (!fEsdTrackCuts->AcceptTrack(esdTrack))
130       continue;
131
132     Double_t p[3];
133     esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
134     TVector3 vector(p);
135
136     Float_t theta = vector.Theta();
137     Float_t eta   = -TMath::Log(TMath::Tan(theta/2.));
138
139     Float_t correction = fdNdEtaCorrection->GetCorrection(vtx[2], eta);
140
141     fdNdEtaAnalysis->FillTrack(vtx[2], eta, correction);
142
143   } // end of track loop
144
145   // for event count per vertex
146   fdNdEtaAnalysis->FillEvent(vtx[2]);
147
148   return kTRUE;
149 }
150
151 void AlidNdEtaAnalysisESDSelector::WriteObjects()
152 {
153   AlidNdEtaAnalysisSelector::WriteObjects();
154
155   if (fEsdTrackCuts)
156     fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");
157
158   if (fdNdEtaCorrection)
159     fdNdEtaCorrection->SaveHistograms();
160 }