3 #include "AlidNdEtaAnalysisESDSelector.h"
13 #include <AliESDVertex.h>
16 #include "esdTrackCuts/AliESDtrackCuts.h"
17 #include "dNdEta/dNdEtaAnalysis.h"
18 #include "dNdEta/AlidNdEtaCorrection.h"
19 #include "AliPWG0Helper.h"
21 ClassImp(AlidNdEtaAnalysisESDSelector)
23 AlidNdEtaAnalysisESDSelector::AlidNdEtaAnalysisESDSelector() :
29 // Constructor. Initialization of pointers
32 //AliLog::SetClassDebugLevel("AlidNdEtaAnalysisESDSelector", AliLog::kDebug);
35 AlidNdEtaAnalysisESDSelector::~AlidNdEtaAnalysisESDSelector()
41 // histograms are in the output list and deleted when the output
42 // list is deleted by the TSelector dtor
45 void AlidNdEtaAnalysisESDSelector::SlaveBegin(TTree* tree)
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).
51 AliSelector::SlaveBegin(tree);
55 printf("Printing input list:\n");
59 if (!fEsdTrackCuts && fInput)
60 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
63 AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
65 fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
68 void AlidNdEtaAnalysisESDSelector::Init(TTree* tree)
70 // read the user objects
72 AliSelector::Init(tree);
74 if (!fEsdTrackCuts && fTree)
75 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fTree->GetUserInfo()->FindObject("AliESDtrackCuts"));
78 AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from user info.");
80 if (!fdNdEtaCorrection && fTree)
81 fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fTree->GetUserInfo()->FindObject("dndeta_correction"));
84 Bool_t AlidNdEtaAnalysisESDSelector::Process(Long64_t entry)
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.
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.
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).
104 if (AliSelector::Process(entry) == kFALSE)
107 // Check prerequisites
110 AliDebug(AliLog::kError, "ESD branch not available");
116 AliDebug(AliLog::kError, "fESDTrackCuts not available");
120 if (!fdNdEtaCorrection)
122 AliDebug(AliLog::kError, "fdNdEtaCorrection not available");
126 if (AliPWG0Helper::IsEventTriggered(fESD) == kFALSE)
129 if (AliPWG0Helper::IsVertexReconstructed(fESD) == kFALSE)
132 // ########################################################
133 // get the EDS vertex
134 const AliESDVertex* vtxESD = fESD->GetVertex();
138 // get number of "good" tracks
139 TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
140 Int_t nGoodTracks = list->GetEntries();
142 Float_t vertexRecoCorr = fdNdEtaCorrection->GetVertexRecoCorrection(vtx[2], nGoodTracks);
143 if (vertexRecoCorr <= 0)
145 AliDebug(AliLog::kError, Form("INFO: Skipping event because vertexRecoCorr is <= 0 (%f)", vertexRecoCorr));
150 // loop over esd tracks
151 for (Int_t t=0; t<nGoodTracks; t++)
153 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(t));
156 AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
161 esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
164 Float_t theta = vector.Theta();
165 Float_t eta = -TMath::Log(TMath::Tan(theta/2.));
166 Float_t pt = vector.Pt();
168 Float_t track2particleCorr = fdNdEtaCorrection->GetTrack2ParticleCorrection(vtx[2], eta, pt);
170 Float_t weight = vertexRecoCorr * track2particleCorr;
173 AliDebug(AliLog::kError, Form("INFO: Skipping track because weight is <= 0 (track %d, weight %f)", t, weight));
177 fdNdEtaAnalysis->FillTrack(vtx[2], eta, pt, weight);
178 } // end of track loop
183 // for event count per vertex
184 fdNdEtaAnalysis->FillEvent(vtx[2], vertexRecoCorr);
189 void AlidNdEtaAnalysisESDSelector::SlaveTerminate()
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.
195 AliSelector::SlaveTerminate();
197 // Add the histograms to the output on each slave server
200 AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
204 fOutput->Add(fdNdEtaAnalysis);
207 void AlidNdEtaAnalysisESDSelector::Terminate()
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.
213 AliSelector::Terminate();
215 fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
217 if (!fdNdEtaAnalysis)
219 AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fdNdEtaAnalysis));
223 TFile* fout = new TFile("analysis_esd.root","RECREATE");
226 fdNdEtaAnalysis->SaveHistograms();
229 fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");