1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 //-----------------------------------------------------------------------
17 // Example of task (running locally, on AliEn and CAF),
18 // which provides standard way of calculating acceptance and efficiency
19 // between different steps of the procedure.
20 // The ouptut of the task is a AliCFContainer from which the efficiencies
22 //-----------------------------------------------------------------------
23 // Author : R. Vernet, Consorzio Cometa - Catania (it)
24 //-----------------------------------------------------------------------
27 #ifndef ALICFSINGLETRACKTASK_CXX
28 #define ALICFSINGLETRACKTASK_CXX
30 #include "AliCFSingleTrackTask.h"
33 #include "TParticle.h"
35 #include "AliMCEvent.h"
36 #include "AliAnalysisManager.h"
37 #include "AliESDEvent.h"
38 #include "AliAODEvent.h"
39 #include "AliCFManager.h"
40 #include "AliCFCutBase.h"
41 #include "AliCFContainer.h"
43 #include "AliESDtrack.h"
46 ClassImp(AliCFSingleTrackTask)
48 //__________________________________________________________________________
49 AliCFSingleTrackTask::AliCFSingleTrackTask() :
54 fHistEventsProcessed(0x0)
60 //___________________________________________________________________________
61 AliCFSingleTrackTask::AliCFSingleTrackTask(const Char_t* name) :
62 AliAnalysisTaskSE(name),
67 fHistEventsProcessed(0x0)
70 // Constructor. Initialization of Inputs and Outputs
72 Info("AliCFSingleTrackTask","Calling Constructor");
75 DefineInput(0) and DefineOutput(0)
76 are taken care of by AliAnalysisTaskSE constructor
78 DefineOutput(1,TH1I::Class());
79 DefineOutput(2,AliCFContainer::Class());
80 DefineOutput(3,TList::Class());
83 //___________________________________________________________________________
84 AliCFSingleTrackTask& AliCFSingleTrackTask::operator=(const AliCFSingleTrackTask& c)
87 // Assignment operator
90 AliAnalysisTaskSE::operator=(c) ;
91 fReadTPCTracks = c.fReadTPCTracks ;
92 fReadAODData = c.fReadAODData ;
93 fCFManager = c.fCFManager;
94 fQAHistList = c.fQAHistList ;
95 fHistEventsProcessed = c.fHistEventsProcessed;
100 //___________________________________________________________________________
101 AliCFSingleTrackTask::AliCFSingleTrackTask(const AliCFSingleTrackTask& c) :
102 AliAnalysisTaskSE(c),
103 fReadTPCTracks(c.fReadTPCTracks),
104 fReadAODData(c.fReadAODData),
105 fCFManager(c.fCFManager),
106 fQAHistList(c.fQAHistList),
107 fHistEventsProcessed(c.fHistEventsProcessed)
114 //___________________________________________________________________________
115 AliCFSingleTrackTask::~AliCFSingleTrackTask() {
119 Info("~AliCFSingleTrackTask","Calling Destructor");
120 if (fCFManager) delete fCFManager ;
121 if (fHistEventsProcessed) delete fHistEventsProcessed ;
122 if (fQAHistList) {fQAHistList->Clear(); delete fQAHistList;}
125 //_________________________________________________
126 void AliCFSingleTrackTask::UserExec(Option_t *)
129 // Main loop function
131 Info("UserExec","") ;
133 AliVEvent* fEvent = fInputEvent ;
134 AliVParticle* track ;
137 Error("UserExec","NO EVENT FOUND!");
142 Error("UserExec","NO MC INFO FOUND");
146 //pass the evt info to the cuts that need it
147 fCFManager->SetMCEventInfo (fMCEvent);
148 fCFManager->SetRecEventInfo(fEvent);
150 // MC-event selection
151 Double_t containerInput[2] ;
153 //loop on the MC event
154 for (Int_t ipart=0; ipart<fMCEvent->GetNumberOfTracks(); ipart++) {
155 AliMCParticle *mcPart = (AliMCParticle*)fMCEvent->GetTrack(ipart);
157 //check the MC-level cuts
158 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
160 containerInput[0] = (Float_t)mcPart->Pt();
161 containerInput[1] = mcPart->Eta() ;
162 //fill the container for Gen-level selection
163 fCFManager->GetParticleContainer()->Fill(containerInput,kStepGenerated);
165 //check the Acceptance-level cuts
166 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
167 //fill the container for Acceptance-level selection
168 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructible);
171 //Now go to rec level
172 for (Int_t iTrack = 0; iTrack<fEvent->GetNumberOfTracks(); iTrack++) {
174 track = fEvent->GetTrack(iTrack);
176 if (fReadTPCTracks) {
178 Error("UserExec","TPC-only tracks are not supported with AOD");
181 AliESDtrack* esdTrack = (AliESDtrack*) track;
182 AliESDtrack* esdTrackTPC = new AliESDtrack();
183 if (!esdTrack->FillTPCOnlyTrack(*esdTrackTPC)) {
184 Error("UserExec","Could not retrieve TPC info");
187 track = esdTrackTPC ;
190 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,track)) continue;
192 // is track associated to particle ?
194 Int_t label = track->GetLabel();
196 if (label<0) continue;
197 AliMCParticle *mcPart = (AliMCParticle*)fMCEvent->GetTrack(label);
199 // check if this track was part of the signal
200 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
205 Double_t pt=TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]);
206 containerInput[0] = pt ;
207 containerInput[1] = track->Eta();
208 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;
210 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,track)) continue ;
211 fCFManager->GetParticleContainer()->Fill(containerInput,kStepSelected);
213 if (fReadTPCTracks) delete track;
216 fHistEventsProcessed->Fill(0);
218 /* PostData(0) is taken care of by AliAnalysisTaskSE */
219 PostData(1,fHistEventsProcessed) ;
220 PostData(2,fCFManager->GetParticleContainer()) ;
221 PostData(3,fQAHistList) ;
225 //___________________________________________________________________________
226 void AliCFSingleTrackTask::Terminate(Option_t*)
228 // The Terminate() function is the last function to be called during
229 // a query. It always runs on the client, it can be used to present
230 // the results graphically or save the results to file.
232 Info("Terminate","");
233 AliAnalysisTaskSE::Terminate();
235 //draw some example plots....
237 AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(2));
239 TH1D* h00 = cont->ShowProjection(0,0) ;
240 TH1D* h01 = cont->ShowProjection(0,1) ;
241 TH1D* h02 = cont->ShowProjection(0,2) ;
242 TH1D* h03 = cont->ShowProjection(0,3) ;
244 TH1D* h10 = cont->ShowProjection(1,0) ;
245 TH1D* h11 = cont->ShowProjection(1,1) ;
246 TH1D* h12 = cont->ShowProjection(1,2) ;
247 TH1D* h13 = cont->ShowProjection(1,3) ;
249 Double_t max1 = h00->GetMaximum();
250 Double_t max2 = h10->GetMaximum();
252 h00->GetYaxis()->SetRangeUser(0,max1*1.2);
253 h01->GetYaxis()->SetRangeUser(0,max1*1.2);
254 h02->GetYaxis()->SetRangeUser(0,max1*1.2);
255 h03->GetYaxis()->SetRangeUser(0,max1*1.2);
257 h10->GetYaxis()->SetRangeUser(0,max2*1.2);
258 h11->GetYaxis()->SetRangeUser(0,max2*1.2);
259 h12->GetYaxis()->SetRangeUser(0,max2*1.2);
260 h13->GetYaxis()->SetRangeUser(0,max2*1.2);
262 h00->SetMarkerStyle(23) ;
263 h01->SetMarkerStyle(24) ;
264 h02->SetMarkerStyle(25) ;
265 h03->SetMarkerStyle(26) ;
267 h10->SetMarkerStyle(23) ;
268 h11->SetMarkerStyle(24) ;
269 h12->SetMarkerStyle(25) ;
270 h13->SetMarkerStyle(26) ;
272 TCanvas * c =new TCanvas("c","",1400,800);
292 c->SaveAs("plots.eps");
296 //___________________________________________________________________________
297 void AliCFSingleTrackTask::UserCreateOutputObjects() {
298 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
299 //TO BE SET BEFORE THE EXECUTION OF THE TASK
301 Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
305 fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;