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 /* $Id: AliAnalysisTaskEx02.cxx 46301 2011-01-06 14:25:27Z agheata $ */
18 /* AliAnalysisTaskEx02.cxx
20 * Template task producing a P_t spectrum and pseudorapidity distribution.
21 * Includes explanations of physics and primary track selections
23 * Instructions for adding histograms can be found below, starting with NEW HISTO
25 * Based on tutorial example from offline pages using AliMultiInputHandler
26 * Edited by Martin Vala
33 #include "AliAnalysisManager.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
37 #include "AliAODTrack.h"
39 #include "AliAODEvent.h"
41 #include "AliAnalysisTaskEx02.h"
42 #include <AliMultiInputEventHandler.h>
43 #include <AliMixInputEventHandler.h>
45 ClassImp(AliAnalysisTaskEx02)
47 //________________________________________________________________________
48 AliAnalysisTaskEx02::AliAnalysisTaskEx02() // All data members should be initialised here
49 : AliAnalysisTaskSE(),
55 fUseLoopInUserExecMix(kFALSE),
56 fUseLoopMixedEvent(kFALSE),
59 fMixingInputHandler(0) // The last in the above list should not have a comma after it
61 // Dummy constructor ALWAYS needed for I/O.
64 //________________________________________________________________________
65 AliAnalysisTaskEx02::AliAnalysisTaskEx02(const char *name) // All data members should be initialised here
66 : AliAnalysisTaskSE(name),
72 fUseLoopInUserExecMix(kFALSE),
73 fUseLoopMixedEvent(kFALSE),
76 fMixingInputHandler(0) // The last in the above list should not have a comma after it
79 // Define input and output slots here (never in the dummy constructor)
80 // Input slot #0 works with a TChain - it is connected to the default input container
81 // Output slot #1 writes into a TH1 container
82 DefineOutput(1, TList::Class()); // for output list
85 //________________________________________________________________________
86 AliAnalysisTaskEx02::~AliAnalysisTaskEx02()
88 // Destructor. Clean-up the output list, but not the histograms that are put inside
89 // (the list is owner and will clean-up these histograms). Protect in PROOF case.
90 if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
95 //________________________________________________________________________
96 void AliAnalysisTaskEx02::UserCreateOutputObjects()
99 // Called once (on the worker node)
101 fOutput = new TList();
102 fOutput->SetOwner(); // IMPORTANT!
106 Float_t ptlow = 0.1, ptup = 3.1;
107 fHistPt = new TH1F("fHistPt", "P_{T} distribution for reconstructed", ptbins, ptlow, ptup);
108 fHistPt->GetXaxis()->SetTitle("P_{T} (GeV/c)");
109 fHistPt->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");
110 fHistPt->SetMarkerStyle(kFullCircle);
113 Float_t etalow = -2.0, etaup = 2.0;
114 fHistEta = new TH1F("fHistEta", "#eta distribution for reconstructed", etabins, etalow, etaup);
115 fHistEta->GetXaxis()->SetTitle("#eta");
116 fHistEta->GetYaxis()->SetTitle("counts");
118 fHistMultiDiff = new TH1F("fHistMultiDiff", "Multiplicity Difference", 100, 0, 100);
119 fHistMultiDiff->GetXaxis()->SetTitle("MultiDiff");
120 fHistMultiDiff->GetYaxis()->SetTitle("counts");
122 fHistZVertexDiff = new TH1F("fHistZVertexDiff", "Multiplicity Difference", 100, 0, 10);
123 fHistZVertexDiff->GetXaxis()->SetTitle("ZVertexDiff");
124 fHistZVertexDiff->GetYaxis()->SetTitle("counts");
126 // NEW HISTO should be defined here, with a sensible name,
128 fOutput->Add(fHistPt);
129 fOutput->Add(fHistEta);
130 fOutput->Add(fHistMultiDiff);
131 fOutput->Add(fHistZVertexDiff);
132 // NEW HISTO added to fOutput here
135 // sets helper pointers for Mixing
136 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
137 SetMainInputHandler(mgr);
138 if (fMainInputHandler) SetMixingInputHandler(fMainInputHandler);
140 PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
143 //________________________________________________________________________
144 void AliAnalysisTaskEx02::UserExec(Option_t *)
147 // Called for each eventn (Simply standard UserExec() like you would do without mixing)
148 // If you want to process events which are used for mixing,
149 // then you have to check for (multiplicity, Vz, ... ) ranges set in event pool.
150 // This action is up to user
154 AliESDEvent *esd = dynamic_cast<AliESDEvent *>(GetMainEvent());
156 if (!fUseLoopInUserExecMix) {
163 AliAODEvent *aod = dynamic_cast<AliAODEvent *>(GetMainEvent());
165 if (!fUseLoopInUserExecMix) {
174 PostData(1, fOutput);
177 //________________________________________________________________________
178 void AliAnalysisTaskEx02::UserExecMix(Option_t *)
181 // Running mixing function
184 if (!fMixingInputHandler) return;
186 Int_t bufferSize = fMixingInputHandler->BufferSize();
187 Int_t numberMixed = fMixingInputHandler->NumberMixed();
188 if (numberMixed == 1) {
189 // you may process Main Event here instead of UserExec()
190 // Here you will have events which are in range of bins in event pool
191 // but only in case when other N (bufferSize) events are found for mix
194 AliESDEvent *esdEvent = dynamic_cast<AliESDEvent *>(GetMainEvent());
197 for(Int_t iBuffer=0; iBuffer<bufferSize; iBuffer++) {
198 AliESDEvent *esdEventMix = dynamic_cast<AliESDEvent *>(GetMixedEvent(iBuffer));
199 AliDebug(AliLog::kDebug, Form("Multi=%d MultiMix=%d", esdEvent->GetNumberOfTracks(), esdEventMix->GetNumberOfTracks()));
200 // AliDebug(AliLog::kDebug, Form("Mask=%lld MaskMix=%lld", esdEvent->GetTriggerMask(), esdEvent->GetTriggerMask()));
201 fHistMultiDiff->Fill(TMath::Abs(esdEvent->GetNumberOfTracks() - esdEventMix->GetNumberOfTracks()));
203 // For tesing purpose only (no physics)
204 if (fUseLoopInUserExecMix) {
206 if (fUseLoopMixedEvent) LoopV0(esdEventMix);
207 else LoopV0(esdEvent);
210 if (fUseLoopMixedEvent) Loop(esdEventMix);
217 AliAODEvent *aodEvent = dynamic_cast<AliAODEvent *>(GetMainEvent());
219 for(Int_t iBuffer=0; iBuffer<bufferSize; iBuffer++) {
220 AliAODEvent *aodEventMix = dynamic_cast<AliAODEvent *>(GetMixedEvent(iBuffer));
221 AliDebug(AliLog::kDebug, Form("Multi=%d MultiMix=%d", aodEvent->GetNumberOfTracks(), aodEventMix->GetNumberOfTracks()));
222 // AliDebug(AliLog::kDebug, Form("Mask=%lld MaskMix=%lld", aodEvent->GetTriggerMask(), aodEventMix->GetTriggerMask()));
223 fHistMultiDiff->Fill(TMath::Abs(aodEvent->GetNumberOfTracks() - aodEventMix->GetNumberOfTracks()));
225 // For tesing purpose only (no physics)
226 if (fUseLoopInUserExecMix) {
228 if (fUseLoopMixedEvent) LoopV0(aodEventMix);
229 else LoopV0(aodEvent);
232 if (fUseLoopMixedEvent) Loop(aodEventMix);
236 AliAODVertex *aodVertex = aodEvent->GetPrimaryVertex();
237 AliAODVertex *aodVertexMix = aodEventMix->GetPrimaryVertex();
238 if ( aodVertex && aodVertexMix) {
239 AliDebug(AliLog::kDebug, Form("Vz=%f VzMix=%f", aodVertex->GetZ(), aodVertexMix->GetZ()));
240 fHistZVertexDiff->Fill(TMath::Abs( aodVertex->GetZ()-aodVertexMix->GetZ()));
248 PostData(1, fOutput);
251 //________________________________________________________________________
252 void AliAnalysisTaskEx02::Terminate(Option_t *)
254 // Draw result to screen, or perform fitting, normalizations
255 // Called once at the end of the query
257 fOutput = dynamic_cast<TList *>(GetOutputData(1));
258 if (!fOutput) { AliError("Could not retrieve TList fOutput"); return; }
260 fHistPt = dynamic_cast<TH1F *>(fOutput->FindObject("fHistPt"));
261 if (!fHistPt) { AliError("Could not retrieve fHistPt"); return;}
263 fHistEta = dynamic_cast<TH1F *>(fOutput->FindObject("fHistEta"));
264 if (!fHistEta) { AliError("Could not retrieve fHistEta"); return;}
266 fHistMultiDiff = dynamic_cast<TH1F *>(fOutput->FindObject("fHistMultiDiff"));
267 if (!fHistMultiDiff) { AliError("Could not retrieve fHistMultiDiff"); return;}
269 fHistZVertexDiff = dynamic_cast<TH1F *>(fOutput->FindObject("fHistZVertexDiff"));
270 if (!fHistZVertexDiff) { AliError("Could not retrieve fHistZVertexDiff"); return;}
273 // NEW HISTO should be retrieved from the TList container in the above way,
274 // so it is available to draw on a canvas such as below
276 TCanvas *c = new TCanvas("AliAnalysisTaskEx02", "P_{T} & #eta", 10, 10, 820, 410);
279 fHistPt->DrawCopy("E");
281 fHistEta->DrawCopy("E");
283 fHistMultiDiff->DrawCopy();
285 fHistZVertexDiff->DrawCopy();
290 //________________________________________________________________________
291 void AliAnalysisTaskEx02::Loop(AliESDEvent *esd)
293 // Track loop for reconstructed event
294 Int_t ntracks = esd->GetNumberOfTracks();
295 for (Int_t i = 0; i < ntracks; i++) {
296 AliESDtrack *esdtrack = esd->GetTrack(i); // pointer to reconstructed to track
298 AliError(Form("ERROR: Could not retrieve esdtrack %d", i));
301 fHistPt->Fill(esdtrack->Pt());
302 fHistEta->Fill(esdtrack->Eta());
306 //________________________________________________________________________
307 void AliAnalysisTaskEx02::LoopV0(AliESDEvent *esd)
310 // V0 loop for reconstructed event (ESD)
313 Int_t nv0 = esd->GetNumberOfV0s();
314 Int_t lIndexTrackPos;
315 AliESDtrack *myTrackPosTest;
316 for (Int_t i = 0; i < nv0; i++)
318 AliESDv0 *esdv0 = esd->GetV0(i);
320 AliError(Form("ERROR: Could not retrieve esdv0 %d", i));
323 lIndexTrackPos = TMath::Abs(esdv0->GetPindex());
324 myTrackPosTest = esd->GetTrack(lIndexTrackPos);
325 fHistPt->Fill(myTrackPosTest->Pt());
326 fHistEta->Fill(myTrackPosTest->Eta());
331 //________________________________________________________________________
332 void AliAnalysisTaskEx02::LoopESDMC()
339 //________________________________________________________________________
340 void AliAnalysisTaskEx02::Loop(AliAODEvent *aod)
343 // Loops over AOD event
346 // Track loop for reconstructed event
347 Int_t ntracks = aod->GetNumberOfTracks();
348 for (Int_t i = 0; i < ntracks; i++) {
349 AliAODTrack *aodTrack = aod->GetTrack(i); // pointer to reconstructed to track
351 AliError(Form("ERROR: Could not retrieve esdtrack %d", i));
355 fHistPt->Fill(aodTrack->Pt());
356 fHistEta->Fill(aodTrack->Eta());
362 //________________________________________________________________________
363 void AliAnalysisTaskEx02::LoopV0(AliAODEvent *aod)
366 // V0 loop for reconstructed event (AOD)
369 Int_t nv0 = aod->GetNumberOfV0s();
370 AliAODTrack *postrackmix;
371 for (Int_t i = 0; i < nv0; i++)
373 AliAODv0 *aodv0 = dynamic_cast<AliAODv0 *>(aod->GetV0(i));
375 AliError(Form("ERROR: Could not retrieve aodv0 %d", i));
378 postrackmix = (AliAODTrack *)aodv0->GetDaughter(0);
379 fHistPt->Fill(postrackmix->Pt());
380 fHistEta->Fill(postrackmix->Eta());
385 //________________________________________________________________________
386 void AliAnalysisTaskEx02::LoopAODMC()
393 //________________________________________________________________________
394 AliVEvent *AliAnalysisTaskEx02::GetMainEvent()
397 // Access to MainEvent
400 AliMultiInputEventHandler *inEvHMainMulti = fMainInputHandler;
401 if (inEvHMainMulti) {
402 AliInputEventHandler *inEvMain = dynamic_cast<AliInputEventHandler *>(inEvHMainMulti->GetFirstInputEventHandler());
403 if (inEvMain) return inEvMain->GetEvent();
409 //________________________________________________________________________
410 AliVEvent *AliAnalysisTaskEx02::GetMixedEvent(Int_t buffId)
413 // Access to Mixed event with buffer id
416 AliMultiInputEventHandler *inEvHMain = fMainInputHandler;
419 AliMixInputEventHandler *mixIH = fMixingInputHandler;
420 if (!mixIH) return 0;
421 if (mixIH->CurrentBinIndex() < 0) {
422 AliDebug(AliLog::kDebug + 1, "Current event mixEH->CurrentEntry() == -1");
425 // AliMultiInputEventHandler *inEvHMixedCurrent = mixEH->GetFirstMultiInputHandler();
426 AliMultiInputEventHandler *inEvHMixedCurrent = dynamic_cast<AliMultiInputEventHandler *>(mixIH->InputEventHandler(buffId));
427 if (!inEvHMixedCurrent) return 0;
428 AliInputEventHandler *ihMixedCurrent = inEvHMixedCurrent->GetFirstInputEventHandler();
429 if (ihMixedCurrent) return ihMixedCurrent->GetEvent();
435 //________________________________________________________________________
436 AliMultiInputEventHandler *AliAnalysisTaskEx02::SetMainInputHandler(AliAnalysisManager *mgr)
439 // Sets main input handler
442 if (!fMainInputHandler) fMainInputHandler = dynamic_cast<AliMultiInputEventHandler *>(mgr->GetInputEventHandler());
444 return fMainInputHandler;
447 //________________________________________________________________________
448 AliMixInputEventHandler *AliAnalysisTaskEx02::SetMixingInputHandler(AliMultiInputEventHandler *mainIH)
451 // Sets mixing input handler
454 if (!fMixingInputHandler) fMixingInputHandler = dynamic_cast<AliMixInputEventHandler *>(mainIH->GetFirstMultiInputHandler());
456 return fMixingInputHandler;