]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTask.cxx
895d12e19fc0271a1cb2740c6509ad29dc9b8746
[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 #include <TROOT.h>
30 #include <TInterpreter.h>
31
32 #include "AliCFSingleTrackTask.h"
33 #include "TCanvas.h"
34 #include "AliStack.h"
35 #include "TParticle.h"
36 #include "TH1I.h"
37 #include "TChain.h"
38 #include "AliMCEventHandler.h"
39 #include "AliMCEvent.h"
40 #include "AliAnalysisManager.h"
41 #include "AliESDEvent.h"
42 #include "AliCFManager.h"
43 #include "AliCFCutBase.h"
44 #include "AliCFContainer.h"
45 #include "TChain.h"
46 #include "AliESDtrack.h"
47 #include "AliLog.h"
48 ClassImp(AliCFSingleTrackTask)
49 //__________________________________________________________________________
50 AliCFSingleTrackTask::AliCFSingleTrackTask() :
51   fChain(0x0),
52   fESD(0x0),
53   fCFManager(0x0),
54   fQAHistList(0x0),
55   fHistEventsProcessed(0x0)
56 {
57 //Default ctor
58 }
59 //___________________________________________________________________________
60 AliCFSingleTrackTask::AliCFSingleTrackTask(const Char_t* name) :
61   AliAnalysisTask(name,"AliCFSingleTrackTask"),
62   fChain(0x0),
63   fESD(0x0),
64   fCFManager(0x0),
65   fQAHistList(0x0),
66   fHistEventsProcessed(0x0)
67 {
68   //
69   // Constructor. Initialization of Inputs and Outputs
70   //
71   Info("AliCFSingleTrackTask","Calling Constructor");
72   DefineInput (0,TChain::Class());
73   DefineOutput(0,TH1I::Class());
74   DefineOutput(1,AliCFContainer::Class());
75   DefineOutput(2,TList::Class());
76 }
77
78 //___________________________________________________________________________
79 AliCFSingleTrackTask& AliCFSingleTrackTask::operator=(const AliCFSingleTrackTask& c) 
80 {
81   //
82   // Assignment operator
83   //
84   if (this!=&c) {
85     AliAnalysisTask::operator=(c) ;
86     fChain      = c.fChain;
87     fESD        = c.fESD;
88     fCFManager  = c.fCFManager;
89     fQAHistList = c.fQAHistList ;
90     fHistEventsProcessed = c.fHistEventsProcessed;
91   }
92   return *this;
93 }
94
95 //___________________________________________________________________________
96 AliCFSingleTrackTask::AliCFSingleTrackTask(const AliCFSingleTrackTask& c) :
97   AliAnalysisTask(c),
98   fChain(c.fChain),
99   fESD(c.fESD),
100   fCFManager(c.fCFManager),
101   fQAHistList(c.fQAHistList),
102   fHistEventsProcessed(c.fHistEventsProcessed)
103 {
104   //
105   // Copy Constructor
106   //
107 }
108
109 //___________________________________________________________________________
110 AliCFSingleTrackTask::~AliCFSingleTrackTask() {
111   //
112   //destructor
113   //
114   Info("~AliCFSingleTrackTask","Calling Destructor");
115   if (fChain)               delete fChain ;
116   if (fESD)                 delete fESD ;
117   if (fCFManager)           delete fCFManager ;
118   if (fHistEventsProcessed) delete fHistEventsProcessed ;
119   if (fQAHistList) {fQAHistList->Clear(); delete fQAHistList;}
120 }
121
122 //___________________________________________________________________________
123
124 void AliCFSingleTrackTask::Init()
125 {
126
127 }
128 //_________________________________________________
129 void AliCFSingleTrackTask::Exec(Option_t *)
130 {
131   //
132   // Main loop function
133   //
134   Info("Exec","") ;
135   // Get the mc truth
136   AliMCEventHandler* mcTruth = (AliMCEventHandler*)((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
137
138   if (!mcTruth) Error("Exec","NO MC INFO FOUND... EXITING\n");
139
140   // transform possible old AliESD into AliESDEvent
141
142   if (fESD->GetAliESDOld()) fESD->CopyFromOldESD(); //transition to new ESD format
143
144   //pass the MC evt handler to the cuts that need it 
145
146   fCFManager->SetEventInfo(mcTruth);
147
148   // Get the MC event 
149   AliMCEvent* mcEvent = mcTruth->MCEvent();
150
151   // MC-event selection
152   Double_t containerInput[2] ;
153         
154   //loop on the MC event
155   for (Int_t ipart=0; ipart<mcEvent->GetNumberOfTracks(); ipart++) { 
156     AliMCParticle *mcPart  = mcEvent->GetTrack(ipart);
157
158     //check the MC-level cuts
159     fCFManager->FillQABeforeParticleCuts(AliCFManager::kPartGenCuts,mcPart);  
160     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
161     fCFManager->FillQAAfterParticleCuts(AliCFManager::kPartGenCuts,mcPart);  
162
163     containerInput[0] = (Float_t)mcPart->Pt();
164     containerInput[1] = mcPart->Eta() ;
165     //fill the container for Gen-level selection
166     fCFManager->GetParticleContainer()->Fill(containerInput,kStepGenerated);
167     
168     //check the Acceptance-level cuts
169     fCFManager->FillQABeforeParticleCuts(AliCFManager::kPartAccCuts,mcPart);  
170     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
171     fCFManager->FillQAAfterParticleCuts(AliCFManager::kPartAccCuts,mcPart);  
172     //fill the container for Acceptance-level selection
173     fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructible);
174   }    
175
176   //Now go to rec level
177   for (Int_t iTrack = 0; iTrack<fESD->GetNumberOfTracks(); iTrack++) {
178
179     AliESDtrack* track = fESD->GetTrack(iTrack);
180     
181     fCFManager->FillQABeforeParticleCuts(AliCFManager::kPartRecCuts,track);  
182     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,track)) continue;
183     fCFManager->FillQAAfterParticleCuts(AliCFManager::kPartRecCuts,track);
184
185     // is track associated to particle ?
186     if (track->GetLabel()<0) continue;
187     AliMCParticle *mcPart  = mcEvent->GetTrack(track->GetLabel());
188     
189     // check if this track was part of the signal
190     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue; 
191     
192     //fill the container
193     Double_t mom[3];
194     track->GetPxPyPz(mom);
195     Double_t pt=TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]);
196     containerInput[0] = pt ;
197     containerInput[1] = track->Eta();
198     fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;   
199
200
201     fCFManager->FillQABeforeParticleCuts(AliCFManager::kPartSelCuts,track);  
202     if (!fCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,track)) continue ;
203     fCFManager->FillQAAfterParticleCuts(AliCFManager::kPartSelCuts,track);  
204
205     fCFManager->GetParticleContainer()->Fill(containerInput,kStepSelected);
206   }
207   
208   fHistEventsProcessed->Fill(0);
209   PostData(0,fHistEventsProcessed) ;
210   PostData(1,fCFManager->GetParticleContainer()) ;
211   PostData(2,fQAHistList) ;
212 }
213
214
215 //___________________________________________________________________________
216 void AliCFSingleTrackTask::Terminate(Option_t*)
217 {
218   // The Terminate() function is the last function to be called during
219   // a query. It always runs on the client, it can be used to present
220   // the results graphically or save the results to file.
221
222   Info("Terminate","");
223   AliAnalysisTask::Terminate();
224
225   //draw some example plots....
226
227   AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(1));
228
229   Double_t max1 = cont->ShowProjection(0,0)->GetMaximum();
230   Double_t max2 = cont->ShowProjection(1,0)->GetMaximum();
231
232   cont->ShowProjection(0,0)->GetYaxis()->SetRangeUser(0,max1*1.2);
233   cont->ShowProjection(0,1)->GetYaxis()->SetRangeUser(0,max1*1.2);
234   cont->ShowProjection(0,2)->GetYaxis()->SetRangeUser(0,max1*1.2);
235   cont->ShowProjection(0,3)->GetYaxis()->SetRangeUser(0,max1*1.2);
236
237   cont->ShowProjection(1,0)->GetYaxis()->SetRangeUser(0,max2*1.2);
238   cont->ShowProjection(1,1)->GetYaxis()->SetRangeUser(0,max2*1.2);
239   cont->ShowProjection(1,2)->GetYaxis()->SetRangeUser(0,max2*1.2);
240   cont->ShowProjection(1,3)->GetYaxis()->SetRangeUser(0,max2*1.2);
241
242   cont->ShowProjection(0,0)->SetMarkerStyle(23) ;
243   cont->ShowProjection(0,1)->SetMarkerStyle(24) ;
244   cont->ShowProjection(0,2)->SetMarkerStyle(25) ;
245   cont->ShowProjection(0,3)->SetMarkerStyle(26) ;
246
247   cont->ShowProjection(1,0)->SetMarkerStyle(23) ;
248   cont->ShowProjection(1,1)->SetMarkerStyle(24) ;
249   cont->ShowProjection(1,2)->SetMarkerStyle(25) ;
250   cont->ShowProjection(1,3)->SetMarkerStyle(26) ;
251
252   TCanvas * c =new TCanvas("c","",1400,800);
253   c->Divide(4,2);
254
255   c->cd(1);
256   cont->ShowProjection(0,0)->Draw("p");
257   c->cd(2);
258   cont->ShowProjection(0,1)->Draw("p");
259   c->cd(3);
260   cont->ShowProjection(0,2)->Draw("p");
261   c->cd(4);
262   cont->ShowProjection(0,3)->Draw("p");
263   c->cd(5);
264   cont->ShowProjection(1,0)->Draw("p");
265   c->cd(6);
266   cont->ShowProjection(1,1)->Draw("p");
267   c->cd(7);
268   cont->ShowProjection(1,2)->Draw("p");
269   c->cd(8);
270   cont->ShowProjection(1,3)->Draw("p");
271
272   c->SaveAs("plots.eps");
273
274   delete fHistEventsProcessed ;
275 }
276
277 //___________________________________________________________________________
278 void AliCFSingleTrackTask::ConnectInputData(Option_t *) {
279   //
280   // Initialize branches.
281   //
282   Info("ConnectInputData","ConnectInputData of task %s\n",GetName());
283
284   fChain = (TChain*)GetInputData(0);
285   fChain->SetBranchStatus("*FMD*",0);
286   fChain->SetBranchStatus("*CaloClusters*",0);
287   fESD = new AliESDEvent();
288   fESD->ReadFromTree(fChain);
289 }
290
291 //___________________________________________________________________________
292 void AliCFSingleTrackTask::CreateOutputObjects() {
293   //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
294   //TO BE SET BEFORE THE EXECUTION OF THE TASK
295   //
296   Info("CreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
297
298   //slot #0
299   fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
300
301   //slot #2
302   fQAHistList = new TList();
303   fCFManager->InitQAHistos();
304   fCFManager->AddQAHistosToList(fQAHistList);
305 }
306
307 #endif