]> git.uio.no Git - u/mrichter/AliRoot.git/blame - CORRFW/test/AliCFSingleTrackTask.cxx
in Terminate(), to correctly display the example plots on the container, retrieve...
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFSingleTrackTask.cxx
CommitLineData
563113d0 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//-----------------------------------------------------------------------
3170dff4 17// Example of task (running locally, on AliEn and CAF),
563113d0 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"
48ClassImp(AliCFSingleTrackTask)
49//__________________________________________________________________________
50AliCFSingleTrackTask::AliCFSingleTrackTask() :
51 fChain(0x0),
52 fESD(0x0),
53 fCFManager(0x0),
54 fQAHistList(0x0),
55 fHistEventsProcessed(0x0)
56{
3170dff4 57//Default ctor
563113d0 58}
59//___________________________________________________________________________
60AliCFSingleTrackTask::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//___________________________________________________________________________
79AliCFSingleTrackTask& 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//___________________________________________________________________________
96AliCFSingleTrackTask::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//___________________________________________________________________________
110AliCFSingleTrackTask::~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
124void AliCFSingleTrackTask::Init()
125{
126
127}
128//_________________________________________________
129void 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
1e9dad92 152 Double_t containerInput[2] ;
563113d0 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 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
160
161 containerInput[0] = (Float_t)mcPart->Pt();
162 containerInput[1] = mcPart->Eta() ;
163 //fill the container for Gen-level selection
164 fCFManager->GetParticleContainer()->Fill(containerInput,kStepGenerated);
165
166 //check the Acceptance-level cuts
167 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
168 //fill the container for Acceptance-level selection
169 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructible);
170 }
171
172 //Now go to rec level
173 for (Int_t iTrack = 0; iTrack<fESD->GetNumberOfTracks(); iTrack++) {
174
175 AliESDtrack* track = fESD->GetTrack(iTrack);
176
177 fCFManager->FillQABeforeParticleCuts(AliCFManager::kPartRecCuts,track);
178 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,track)) continue;
179 fCFManager->FillQAAfterParticleCuts(AliCFManager::kPartRecCuts,track);
180
181 // is track associated to particle ?
182 if (track->GetLabel()<0) continue;
183 AliMCParticle *mcPart = mcEvent->GetTrack(track->GetLabel());
184
185 // check if this track was part of the signal
186 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
187
188 //fill the container
189 Double_t mom[3];
190 track->GetPxPyPz(mom);
191 Double_t pt=TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]);
192 containerInput[0] = pt ;
193 containerInput[1] = track->Eta();
194 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;
195 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,track)) continue ;
196 fCFManager->GetParticleContainer()->Fill(containerInput,kStepSelected);
197 }
198
199 fHistEventsProcessed->Fill(0);
200 PostData(0,fHistEventsProcessed) ;
201 PostData(1,fCFManager->GetParticleContainer()) ;
202 PostData(2,fQAHistList) ;
203}
204
205
206//___________________________________________________________________________
207void AliCFSingleTrackTask::Terminate(Option_t*)
208{
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.
212
213 Info("Terminate","");
214 AliAnalysisTask::Terminate();
215
3170dff4 216 //draw some example plots....
563113d0 217
3170dff4 218 AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(1));
563113d0 219
3170dff4 220 Double_t max1 = cont->ShowProjection(0,0)->GetMaximum();
221 Double_t max2 = cont->ShowProjection(1,0)->GetMaximum();
563113d0 222
3170dff4 223 cont->ShowProjection(0,0)->GetYaxis()->SetRangeUser(0,max1*1.2);
224 cont->ShowProjection(0,1)->GetYaxis()->SetRangeUser(0,max1*1.2);
225 cont->ShowProjection(0,2)->GetYaxis()->SetRangeUser(0,max1*1.2);
226 cont->ShowProjection(0,3)->GetYaxis()->SetRangeUser(0,max1*1.2);
563113d0 227
3170dff4 228 cont->ShowProjection(1,0)->GetYaxis()->SetRangeUser(0,max2*1.2);
229 cont->ShowProjection(1,1)->GetYaxis()->SetRangeUser(0,max2*1.2);
230 cont->ShowProjection(1,2)->GetYaxis()->SetRangeUser(0,max2*1.2);
231 cont->ShowProjection(1,3)->GetYaxis()->SetRangeUser(0,max2*1.2);
563113d0 232
3170dff4 233 cont->ShowProjection(0,0)->SetMarkerStyle(23) ;
234 cont->ShowProjection(0,1)->SetMarkerStyle(24) ;
235 cont->ShowProjection(0,2)->SetMarkerStyle(25) ;
236 cont->ShowProjection(0,3)->SetMarkerStyle(26) ;
237
238 cont->ShowProjection(1,0)->SetMarkerStyle(23) ;
239 cont->ShowProjection(1,1)->SetMarkerStyle(24) ;
240 cont->ShowProjection(1,2)->SetMarkerStyle(25) ;
241 cont->ShowProjection(1,3)->SetMarkerStyle(26) ;
563113d0 242
243 TCanvas * c =new TCanvas("c","",1400,800);
244 c->Divide(4,2);
245
563113d0 246 c->cd(1);
3170dff4 247 cont->ShowProjection(0,0)->Draw("p");
563113d0 248 c->cd(2);
3170dff4 249 cont->ShowProjection(0,1)->Draw("p");
563113d0 250 c->cd(3);
3170dff4 251 cont->ShowProjection(0,2)->Draw("p");
563113d0 252 c->cd(4);
3170dff4 253 cont->ShowProjection(0,3)->Draw("p");
563113d0 254 c->cd(5);
3170dff4 255 cont->ShowProjection(1,0)->Draw("p");
563113d0 256 c->cd(6);
3170dff4 257 cont->ShowProjection(1,1)->Draw("p");
563113d0 258 c->cd(7);
3170dff4 259 cont->ShowProjection(1,2)->Draw("p");
563113d0 260 c->cd(8);
3170dff4 261 cont->ShowProjection(1,3)->Draw("p");
563113d0 262
263 c->SaveAs("plots.eps");
264
265 delete fHistEventsProcessed ;
266}
267
268//___________________________________________________________________________
269void AliCFSingleTrackTask::ConnectInputData(Option_t *) {
270 //
271 // Initialize branches.
272 //
273 Info("ConnectInputData","ConnectInputData of task %s\n",GetName());
274
275 fChain = (TChain*)GetInputData(0);
276 fChain->SetBranchStatus("*FMD*",0);
277 fChain->SetBranchStatus("*CaloClusters*",0);
278 fESD = new AliESDEvent();
279 fESD->ReadFromTree(fChain);
280}
281
282//___________________________________________________________________________
283void AliCFSingleTrackTask::CreateOutputObjects() {
284 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
285 //TO BE SET BEFORE THE EXECUTION OF THE TASK
286 //
287 Info("CreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
288
289 //slot #0
290 fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
291
292 //slot #2
293 fQAHistList = new TList();
294 fCFManager->InitQAHistos();
295 fCFManager->AddQAHistosToList(fQAHistList);
296}
297
298#endif