]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTask.cxx
Added new functionalities for projections and slices.
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFSingleTrackTask.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
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
21 // can be calculated
22 //-----------------------------------------------------------------------
23 // Author : R. Vernet, Consorzio Cometa - Catania (it)
24 //-----------------------------------------------------------------------
25
26
27 #ifndef ALICFSINGLETRACKTASK_CXX
28 #define ALICFSINGLETRACKTASK_CXX
29
30 #include "AliCFSingleTrackTask.h"
31 #include "TCanvas.h"
32 #include "AliStack.h"
33 #include "TParticle.h"
34 #include "TH1I.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"
42 #include "TChain.h"
43 #include "AliESDtrack.h"
44 #include "AliLog.h"
45
46 ClassImp(AliCFSingleTrackTask)
47
48 //__________________________________________________________________________
49 AliCFSingleTrackTask::AliCFSingleTrackTask() :
50   fReadTPCTracks(0),
51   fReadAODData(0),
52   fCFManager(0x0),
53   fQAHistList(0x0),
54   fHistEventsProcessed(0x0)
55 {
56   //
57   //Default ctor
58   //
59 }
60 //___________________________________________________________________________
61 AliCFSingleTrackTask::AliCFSingleTrackTask(const Char_t* name) :
62   AliAnalysisTaskSE(name),
63   fReadTPCTracks(0),
64   fReadAODData(0),
65   fCFManager(0x0),
66   fQAHistList(0x0),
67   fHistEventsProcessed(0x0)
68 {
69   //
70   // Constructor. Initialization of Inputs and Outputs
71   //
72   Info("AliCFSingleTrackTask","Calling Constructor");
73
74   /*
75     DefineInput(0) and DefineOutput(0)
76     are taken care of by AliAnalysisTaskSE constructor
77   */
78   DefineOutput(1,TH1I::Class());
79   DefineOutput(2,AliCFContainer::Class());
80   DefineOutput(3,TList::Class());
81 }
82
83 //___________________________________________________________________________
84 AliCFSingleTrackTask& AliCFSingleTrackTask::operator=(const AliCFSingleTrackTask& c) 
85 {
86   //
87   // Assignment operator
88   //
89   if (this!=&c) {
90     AliAnalysisTaskSE::operator=(c) ;
91     fReadTPCTracks = c.fReadTPCTracks ;
92     fReadAODData = c.fReadAODData ;
93     fCFManager  = c.fCFManager;
94     fQAHistList = c.fQAHistList ;
95     fHistEventsProcessed = c.fHistEventsProcessed;
96   }
97   return *this;
98 }
99
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)
108 {
109   //
110   // Copy Constructor
111   //
112 }
113
114 //___________________________________________________________________________
115 AliCFSingleTrackTask::~AliCFSingleTrackTask() {
116   //
117   //destructor
118   //
119   Info("~AliCFSingleTrackTask","Calling Destructor");
120   if (fCFManager)           delete fCFManager ;
121   if (fHistEventsProcessed) delete fHistEventsProcessed ;
122   if (fQAHistList) {fQAHistList->Clear(); delete fQAHistList;}
123 }
124
125 //_________________________________________________
126 void AliCFSingleTrackTask::UserExec(Option_t *)
127 {
128   //
129   // Main loop function
130   //
131   Info("UserExec","") ;
132
133   AliVEvent*    fEvent = fInputEvent ;
134   AliVParticle* track ;
135   
136   if (!fEvent) {
137     Error("UserExec","NO EVENT FOUND!");
138     return;
139   }
140
141   if (!fMCEvent) Error("UserExec","NO MC INFO FOUND");
142   
143   //pass the MC evt handler to the cuts that need it 
144   fCFManager->SetEventInfo(fMCEvent);
145
146   // MC-event selection
147   Double_t containerInput[2] ;
148         
149   //loop on the MC event
150   for (Int_t ipart=0; ipart<fMCEvent->GetNumberOfTracks(); ipart++) { 
151     AliMCParticle *mcPart  = fMCEvent->GetTrack(ipart);
152
153     //check the MC-level cuts
154     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
155
156     containerInput[0] = (Float_t)mcPart->Pt();
157     containerInput[1] = mcPart->Eta() ;
158     //fill the container for Gen-level selection
159     fCFManager->GetParticleContainer()->Fill(containerInput,kStepGenerated);
160
161     //check the Acceptance-level cuts
162     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
163     //fill the container for Acceptance-level selection
164     fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructible);
165   }    
166
167   //Now go to rec level
168   for (Int_t iTrack = 0; iTrack<fEvent->GetNumberOfTracks(); iTrack++) {
169     
170     track = fEvent->GetTrack(iTrack);
171     
172     if (fReadTPCTracks) {
173       if (fReadAODData) {
174         Error("UserExec","TPC-only tracks are not supported with AOD");
175         return ;
176       }
177       AliESDtrack* esdTrack    = (AliESDtrack*) track;
178       AliESDtrack* esdTrackTPC = new AliESDtrack();
179       if (!esdTrack->FillTPCOnlyTrack(*esdTrackTPC)) {
180         Error("UserExec","Could not retrieve TPC info");
181         continue;
182       }
183       track = esdTrackTPC ;
184     }
185
186     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,track)) continue;
187     
188     // is track associated to particle ?
189
190     Int_t label = track->GetLabel();
191
192     if (label<0) continue;
193     AliMCParticle *mcPart  = fMCEvent->GetTrack(label);
194     
195     // check if this track was part of the signal
196     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue; 
197     
198     //fill the container
199     Double_t mom[3];
200     track->PxPyPz(mom);
201     Double_t pt=TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]);
202     containerInput[0] = pt ;
203     containerInput[1] = track->Eta();
204     fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;   
205
206     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,track)) continue ;
207     fCFManager->GetParticleContainer()->Fill(containerInput,kStepSelected);
208
209     if (fReadTPCTracks) delete track;
210   }
211   
212   fHistEventsProcessed->Fill(0);
213
214   /* PostData(0) is taken care of by AliAnalysisTaskSE */
215   PostData(1,fHistEventsProcessed) ;
216   PostData(2,fCFManager->GetParticleContainer()) ;
217   PostData(3,fQAHistList) ;
218 }
219
220
221 //___________________________________________________________________________
222 void AliCFSingleTrackTask::Terminate(Option_t*)
223 {
224   // The Terminate() function is the last function to be called during
225   // a query. It always runs on the client, it can be used to present
226   // the results graphically or save the results to file.
227
228   Info("Terminate","");
229   AliAnalysisTaskSE::Terminate();
230
231   //draw some example plots....
232
233   AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(2));
234
235   TH1D* h00 =   cont->ShowProjection(0,0) ;
236   TH1D* h01 =   cont->ShowProjection(0,1) ;
237   TH1D* h02 =   cont->ShowProjection(0,2) ;
238   TH1D* h03 =   cont->ShowProjection(0,3) ;
239
240   TH1D* h10 =   cont->ShowProjection(1,0) ;
241   TH1D* h11 =   cont->ShowProjection(1,1) ;
242   TH1D* h12 =   cont->ShowProjection(1,2) ;
243   TH1D* h13 =   cont->ShowProjection(1,3) ;
244
245   Double_t max1 = h00->GetMaximum();
246   Double_t max2 = h10->GetMaximum();
247
248   h00->GetYaxis()->SetRangeUser(0,max1*1.2);
249   h01->GetYaxis()->SetRangeUser(0,max1*1.2);
250   h02->GetYaxis()->SetRangeUser(0,max1*1.2);
251   h03->GetYaxis()->SetRangeUser(0,max1*1.2);
252
253   h10->GetYaxis()->SetRangeUser(0,max2*1.2);
254   h11->GetYaxis()->SetRangeUser(0,max2*1.2);
255   h12->GetYaxis()->SetRangeUser(0,max2*1.2);
256   h13->GetYaxis()->SetRangeUser(0,max2*1.2);
257
258   h00->SetMarkerStyle(23) ;
259   h01->SetMarkerStyle(24) ;
260   h02->SetMarkerStyle(25) ;
261   h03->SetMarkerStyle(26) ;
262
263   h10->SetMarkerStyle(23) ;
264   h11->SetMarkerStyle(24) ;
265   h12->SetMarkerStyle(25) ;
266   h13->SetMarkerStyle(26) ;
267
268   TCanvas * c =new TCanvas("c","",1400,800);
269   c->Divide(4,2);
270
271   c->cd(1);
272   h00->Draw("p");
273   c->cd(2);
274   h01->Draw("p");
275   c->cd(3);
276   h02->Draw("p");
277   c->cd(4);
278   h03->Draw("p");
279   c->cd(5);
280   h10->Draw("p");
281   c->cd(6);
282   h11->Draw("p");
283   c->cd(7);
284   h12->Draw("p");
285   c->cd(8);
286   h13->Draw("p");
287
288   c->SaveAs("plots.eps");
289 }
290
291
292 //___________________________________________________________________________
293 void AliCFSingleTrackTask::UserCreateOutputObjects() {
294   //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
295   //TO BE SET BEFORE THE EXECUTION OF THE TASK
296   //
297   Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
298
299   //slot #1
300   OpenFile(1);
301   fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
302
303 //   OpenFile(2);
304 //   OpenFile(3);
305 }
306
307 #endif