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() :
25 fdNdEtaAnalysisMBVtx(0),
32 // Constructor. Initialization of pointers
35 //AliLog::SetClassDebugLevel("AlidNdEtaAnalysisESDSelector", AliLog::kDebug);
38 AlidNdEtaAnalysisESDSelector::~AlidNdEtaAnalysisESDSelector()
44 // histograms are in the output list and deleted when the output
45 // list is deleted by the TSelector dtor
48 void AlidNdEtaAnalysisESDSelector::SlaveBegin(TTree* tree)
50 // The SlaveBegin() function is called after the Begin() function.
51 // When running with PROOF SlaveBegin() is called on each slave server.
52 // The tree argument is deprecated (on PROOF 0 is passed).
54 AliSelector::SlaveBegin(tree);
58 printf("Printing input list:\n");
62 if (!fEsdTrackCuts && fInput)
63 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
65 if (!fEsdTrackCuts && tree)
66 fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
69 AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
72 if (!fdNdEtaCorrection && fInput)
73 fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fInput->FindObject("dndeta_correction"));
75 if (!fdNdEtaCorrection && fTree)
76 fdNdEtaCorrection = dynamic_cast<AlidNdEtaCorrection*> (fTree->GetUserInfo()->FindObject("dndeta_correction"));
78 if (!fdNdEtaCorrection)
79 AliDebug(AliLog::kError, "ERROR: Could not read dndeta_correction from input list.");
82 fdNdEtaAnalysisMBVtx = new dNdEtaAnalysis("dndeta_mbvtx", "dndeta_mbvtx");
83 fdNdEtaAnalysisMB = new dNdEtaAnalysis("dndeta_mb", "dndeta_mb");
84 fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
87 void AlidNdEtaAnalysisESDSelector::Init(TTree* tree)
89 // read the user objects
91 AliSelector::Init(tree);
94 Bool_t AlidNdEtaAnalysisESDSelector::Process(Long64_t entry)
96 // The Process() function is called for each entry in the tree (or possibly
97 // keyed object in the case of PROOF) to be processed. The entry argument
98 // specifies which entry in the currently loaded tree is to be processed.
99 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
100 // to read either all or the required parts of the data. When processing
101 // keyed objects with PROOF, the object is already loaded and is available
102 // via the fObject pointer.
104 // This function should contain the "body" of the analysis. It can contain
105 // simple or elaborate selection criteria, run algorithms on the data
106 // of the event and typically fill histograms.
108 // WARNING when a selector is used with a TChain, you must use
109 // the pointer to the current TTree to call GetEntry(entry).
110 // The entry is always the local entry number in the current tree.
111 // Assuming that fTree is the pointer to the TChain being processed,
112 // use fTree->GetTree()->GetEntry(entry).
114 if (AliSelector::Process(entry) == kFALSE)
117 // Check prerequisites
120 AliDebug(AliLog::kError, "ESD branch not available");
126 AliDebug(AliLog::kError, "fESDTrackCuts not available");
130 if (!fdNdEtaCorrection)
132 AliDebug(AliLog::kError, "fdNdEtaCorrection not available");
136 if (AliPWG0Helper::IsEventTriggered(fESD) == kFALSE)
139 if (AliPWG0Helper::IsVertexReconstructed(fESD) == kFALSE)
142 // ########################################################
143 // get the EDS vertex
144 const AliESDVertex* vtxESD = fESD->GetVertex();
148 // get number of "good" tracks
149 TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
150 Int_t nGoodTracks = list->GetEntries();
152 Float_t vertexRecoCorr = fdNdEtaCorrection->GetVertexRecoCorrection(vtx[2], nGoodTracks);
153 if (vertexRecoCorr <= 0)
155 AliDebug(AliLog::kError, Form("INFO: Skipping event because vertexRecoCorr is <= 0 (%f)", vertexRecoCorr));
160 Float_t triggerCorr = fdNdEtaCorrection->GetTriggerCorrection(vtx[2], nGoodTracks);
161 if (triggerCorr <= 0)
163 AliDebug(AliLog::kError, Form("INFO: Skipping event because triggerCorr is <= 0 (%f)", triggerCorr));
168 // loop over esd tracks
169 for (Int_t t=0; t<nGoodTracks; t++)
171 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(t));
174 AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", t));
179 esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
182 Float_t theta = vector.Theta();
183 Float_t eta = -TMath::Log(TMath::Tan(theta/2.));
184 Float_t pt = vector.Pt();
186 Float_t track2particleCorr = fdNdEtaCorrection->GetTrack2ParticleCorrection(vtx[2], eta, pt);
188 Float_t weight = track2particleCorr * vertexRecoCorr * triggerCorr;
191 AliDebug(AliLog::kError, Form("INFO: Skipping track because weight is <= 0 (track %d, weight %f) (vtx %f, eta %f, pt %f)", t, weight, vtx[2], eta, pt));
195 fdNdEtaAnalysisMBVtx->FillTrack(vtx[2], eta, pt, track2particleCorr);
196 fdNdEtaAnalysisMB->FillTrack(vtx[2], eta, pt, track2particleCorr * vertexRecoCorr);
197 fdNdEtaAnalysis->FillTrack(vtx[2], eta, pt, weight);
198 } // end of track loop
203 // for event count per vertex
204 fdNdEtaAnalysisMBVtx->FillEvent(vtx[2], 1);
205 fdNdEtaAnalysisMB->FillEvent(vtx[2], vertexRecoCorr);
206 fdNdEtaAnalysis->FillEvent(vtx[2], vertexRecoCorr * triggerCorr);
211 void AlidNdEtaAnalysisESDSelector::SlaveTerminate()
213 // The SlaveTerminate() function is called after all entries or objects
214 // have been processed. When running with PROOF SlaveTerminate() is called
215 // on each slave server.
217 AliSelector::SlaveTerminate();
219 // Add the histograms to the output on each slave server
222 AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
226 fOutput->Add(fdNdEtaAnalysis);
229 void AlidNdEtaAnalysisESDSelector::Terminate()
231 // The Terminate() function is the last function to be called during
232 // a query. It always runs on the client, it can be used to present
233 // the results graphically or save the results to file.
235 AliSelector::Terminate();
237 fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
239 if (!fdNdEtaAnalysis)
241 AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fdNdEtaAnalysis));
245 fdNdEtaAnalysis->Finish(fdNdEtaCorrection, 0.3);
246 fdNdEtaAnalysisMB->Finish(fdNdEtaCorrection, 0.3);
247 fdNdEtaAnalysisMBVtx->Finish(fdNdEtaCorrection, 0.3);
249 TFile* fout = new TFile("analysis_esd.root","RECREATE");
252 fdNdEtaAnalysis->SaveHistograms();
254 if (fdNdEtaAnalysisMB)
255 fdNdEtaAnalysisMB->SaveHistograms();
257 if (fdNdEtaAnalysisMBVtx)
258 fdNdEtaAnalysisMBVtx->SaveHistograms();
261 fEsdTrackCuts->SaveHistograms("esd_tracks_cuts");