]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/examples/EventMixing/AliAnalysisTaskEx02.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / ANALYSIS / examples / EventMixing / AliAnalysisTaskEx02.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 /* $Id: AliAnalysisTaskEx02.cxx 46301 2011-01-06 14:25:27Z agheata $ */
17
18 /* AliAnalysisTaskEx02.cxx
19  *
20  * Template task producing a P_t spectrum and pseudorapidity distribution.
21  * Includes explanations of physics and primary track selections
22  *
23  * Instructions for adding histograms can be found below, starting with NEW HISTO
24  *
25  * Based on tutorial example from offline pages using AliMultiInputHandler
26  * Edited by Martin Vala
27  */
28 #include <TCanvas.h>
29 #include <TList.h>
30 #include <TH1.h>
31
32 #include "AliLog.h"
33 #include "AliAnalysisManager.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliESDv0.h"
37 #include "AliAODTrack.h"
38 #include "AliAODv0.h"
39 #include "AliAODEvent.h"
40
41 #include "AliAnalysisTaskEx02.h"
42 #include <AliMultiInputEventHandler.h>
43 #include <AliMixInputEventHandler.h>
44
45 ClassImp(AliAnalysisTaskEx02)
46
47 //________________________________________________________________________
48 AliAnalysisTaskEx02::AliAnalysisTaskEx02() // All data members should be initialised here
49    : AliAnalysisTaskSE(),
50      fOutput(0),
51      fHistPt(0),
52      fHistEta(0),
53      fHistMultiDiff(0),
54      fHistZVertexDiff(0),
55      fUseLoopInUserExecMix(kFALSE),
56      fUseLoopMixedEvent(kFALSE),
57      fUseLoopV0(kFALSE),
58      fMainInputHandler(0),
59      fMixingInputHandler(0) // The last in the above list should not have a comma after it
60 {
61    // Dummy constructor ALWAYS needed for I/O.
62 }
63
64 //________________________________________________________________________
65 AliAnalysisTaskEx02::AliAnalysisTaskEx02(const char *name) // All data members should be initialised here
66    : AliAnalysisTaskSE(name),
67      fOutput(0),
68      fHistPt(0),
69      fHistEta(0),
70      fHistMultiDiff(0),
71      fHistZVertexDiff(0),
72      fUseLoopInUserExecMix(kFALSE),
73      fUseLoopMixedEvent(kFALSE),
74      fUseLoopV0(kFALSE),
75      fMainInputHandler(0),
76      fMixingInputHandler(0) // The last in the above list should not have a comma after it
77 {
78    // Constructor
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
83 }
84
85 //________________________________________________________________________
86 AliAnalysisTaskEx02::~AliAnalysisTaskEx02()
87 {
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()) {
91       delete fOutput;
92    }
93 }
94
95 //________________________________________________________________________
96 void AliAnalysisTaskEx02::UserCreateOutputObjects()
97 {
98    // Create histograms
99    // Called once (on the worker node)
100
101    fOutput = new TList();
102    fOutput->SetOwner();  // IMPORTANT!
103
104    // Create histograms
105    Int_t ptbins = 15;
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);
111
112    Int_t etabins = 40;
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");
117
118    fHistMultiDiff = new TH1F("fHistMultiDiff", "Multiplicity Difference", 100, 0, 100);
119    fHistMultiDiff->GetXaxis()->SetTitle("MultiDiff");
120    fHistMultiDiff->GetYaxis()->SetTitle("counts");
121
122    fHistZVertexDiff = new TH1F("fHistZVertexDiff", "Multiplicity Difference", 100, 0, 10);
123    fHistZVertexDiff->GetXaxis()->SetTitle("ZVertexDiff");
124    fHistZVertexDiff->GetYaxis()->SetTitle("counts");
125
126    // NEW HISTO should be defined here, with a sensible name,
127
128    fOutput->Add(fHistPt);
129    fOutput->Add(fHistEta);
130    fOutput->Add(fHistMultiDiff);
131    fOutput->Add(fHistZVertexDiff);
132    // NEW HISTO added to fOutput here
133
134
135    // sets helper pointers for Mixing
136    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
137    SetMainInputHandler(mgr);
138    if (fMainInputHandler) SetMixingInputHandler(fMainInputHandler);
139
140    PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
141 }
142
143 //________________________________________________________________________
144 void AliAnalysisTaskEx02::UserExec(Option_t *)
145 {
146    // Main loop
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
151    //
152
153
154    AliESDEvent *esd = dynamic_cast<AliESDEvent *>(GetMainEvent());
155    if (esd) {
156       if (!fUseLoopInUserExecMix) {
157          if (fUseLoopV0)
158             LoopV0(esd);
159          else
160             Loop(esd);
161       }
162    } else {
163       AliAODEvent *aod = dynamic_cast<AliAODEvent *>(GetMainEvent());
164       if (aod) {
165          if (!fUseLoopInUserExecMix) {
166             if (fUseLoopV0)
167                LoopV0(aod);
168             else
169                Loop(aod);
170          }
171       }
172    }
173
174    PostData(1, fOutput);
175 }
176
177 //________________________________________________________________________
178 void AliAnalysisTaskEx02::UserExecMix(Option_t *)
179 {
180    //
181    // Running mixing function
182    //
183
184    if (!fMixingInputHandler) return;
185
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
192    }
193
194    AliESDEvent *esdEvent = dynamic_cast<AliESDEvent *>(GetMainEvent());
195    if (esdEvent) {
196
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()));
202
203          // For tesing purpose only (no physics)
204          if (fUseLoopInUserExecMix) {
205             if (fUseLoopV0) {
206                if (fUseLoopMixedEvent) LoopV0(esdEventMix);
207                else LoopV0(esdEvent);
208             }
209             else {
210                if (fUseLoopMixedEvent) Loop(esdEventMix);
211                else Loop(esdEvent);
212             }
213          }
214
215       }
216    } else {
217       AliAODEvent *aodEvent = dynamic_cast<AliAODEvent *>(GetMainEvent());
218       if (aodEvent) {
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()));
224
225             // For tesing purpose only (no physics)
226             if (fUseLoopInUserExecMix) {
227                if (fUseLoopV0) {
228                   if (fUseLoopMixedEvent) LoopV0(aodEventMix);
229                   else LoopV0(aodEvent);
230                }
231                else {
232                   if (fUseLoopMixedEvent) Loop(aodEventMix);
233                   else Loop(aodEvent);
234                }
235             }
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()));
241             }
242
243          }
244
245       }
246    }
247
248    PostData(1, fOutput);
249 }
250
251 //________________________________________________________________________
252 void AliAnalysisTaskEx02::Terminate(Option_t *)
253 {
254    // Draw result to screen, or perform fitting, normalizations
255    // Called once at the end of the query
256
257    fOutput = dynamic_cast<TList *>(GetOutputData(1));
258    if (!fOutput) { AliError("Could not retrieve TList fOutput"); return; }
259
260    fHistPt = dynamic_cast<TH1F *>(fOutput->FindObject("fHistPt"));
261    if (!fHistPt) { AliError("Could not retrieve fHistPt"); return;}
262
263    fHistEta = dynamic_cast<TH1F *>(fOutput->FindObject("fHistEta"));
264    if (!fHistEta) { AliError("Could not retrieve fHistEta"); return;}
265
266    fHistMultiDiff = dynamic_cast<TH1F *>(fOutput->FindObject("fHistMultiDiff"));
267    if (!fHistMultiDiff) { AliError("Could not retrieve fHistMultiDiff"); return;}
268
269    fHistZVertexDiff = dynamic_cast<TH1F *>(fOutput->FindObject("fHistZVertexDiff"));
270    if (!fHistZVertexDiff) { AliError("Could not retrieve fHistZVertexDiff"); return;}
271
272
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
275
276    TCanvas *c = new TCanvas("AliAnalysisTaskEx02", "P_{T} & #eta", 10, 10, 820, 410);
277    c->Divide(2, 2);
278    c->cd(1)->SetLogy();
279    fHistPt->DrawCopy("E");
280    c->cd(2);
281    fHistEta->DrawCopy("E");
282    c->cd(3)->SetLogy();
283    fHistMultiDiff->DrawCopy();
284    c->cd(4)->SetLogy();
285    fHistZVertexDiff->DrawCopy();
286
287
288 }
289
290 //________________________________________________________________________
291 void AliAnalysisTaskEx02::Loop(AliESDEvent *esd)
292 {
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
297       if (!esdtrack) {
298          AliError(Form("ERROR: Could not retrieve esdtrack %d", i));
299          continue;
300       }
301       fHistPt->Fill(esdtrack->Pt());
302       fHistEta->Fill(esdtrack->Eta());
303    }
304 }
305
306 //________________________________________________________________________
307 void AliAnalysisTaskEx02::LoopV0(AliESDEvent *esd)
308 {
309    //
310    // V0 loop for reconstructed event (ESD)
311    //
312
313    Int_t nv0 = esd->GetNumberOfV0s();
314    Int_t lIndexTrackPos;
315    AliESDtrack *myTrackPosTest;
316    for (Int_t i = 0; i < nv0; i++)
317    {
318       AliESDv0 *esdv0 = esd->GetV0(i);
319       if (!esdv0) {
320          AliError(Form("ERROR: Could not retrieve esdv0 %d", i));
321          continue;
322       }
323       lIndexTrackPos = TMath::Abs(esdv0->GetPindex());
324       myTrackPosTest = esd->GetTrack(lIndexTrackPos);
325       fHistPt->Fill(myTrackPosTest->Pt());
326       fHistEta->Fill(myTrackPosTest->Eta());
327    }
328
329 }
330
331 //________________________________________________________________________
332 void AliAnalysisTaskEx02::LoopESDMC()
333 {
334    //
335    // TODO
336    //
337 }
338
339 //________________________________________________________________________
340 void AliAnalysisTaskEx02::Loop(AliAODEvent *aod)
341 {
342    //
343    // Loops over AOD event
344    //
345
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
350       if (!aodTrack) {
351          AliError(Form("ERROR: Could not retrieve esdtrack %d", i));
352          continue;
353       }
354
355       fHistPt->Fill(aodTrack->Pt());
356       fHistEta->Fill(aodTrack->Eta());
357    }
358 }
359
360
361
362 //________________________________________________________________________
363 void AliAnalysisTaskEx02::LoopV0(AliAODEvent *aod)
364 {
365    //
366    // V0 loop for reconstructed event (AOD)
367    //
368
369    Int_t nv0 = aod->GetNumberOfV0s();
370    AliAODTrack *postrackmix;
371    for (Int_t i = 0; i < nv0; i++)
372    {
373       AliAODv0 *aodv0 = dynamic_cast<AliAODv0 *>(aod->GetV0(i));
374       if (!aodv0) {
375          AliError(Form("ERROR: Could not retrieve aodv0 %d", i));
376          continue;
377       }
378       postrackmix = (AliAODTrack *)aodv0->GetDaughter(0);
379       fHistPt->Fill(postrackmix->Pt());
380       fHistEta->Fill(postrackmix->Eta());
381    }
382
383 }
384
385 //________________________________________________________________________
386 void AliAnalysisTaskEx02::LoopAODMC()
387 {
388    //
389    // TODO
390    //
391 }
392
393 //________________________________________________________________________
394 AliVEvent *AliAnalysisTaskEx02::GetMainEvent()
395 {
396    //
397    // Access to MainEvent
398    //
399
400    AliMultiInputEventHandler *inEvHMainMulti = fMainInputHandler;
401    if (inEvHMainMulti) {
402       AliInputEventHandler *inEvMain = dynamic_cast<AliInputEventHandler *>(inEvHMainMulti->GetFirstInputEventHandler());
403       if (inEvMain) return inEvMain->GetEvent();
404    }
405
406    return 0;
407 }
408
409 //________________________________________________________________________
410 AliVEvent *AliAnalysisTaskEx02::GetMixedEvent(Int_t buffId)
411 {
412    //
413    // Access to Mixed event with buffer id
414    //
415
416    AliMultiInputEventHandler *inEvHMain = fMainInputHandler;
417    if (inEvHMain) {
418
419       AliMixInputEventHandler *mixIH = fMixingInputHandler;
420       if (!mixIH) return 0;
421       if (mixIH->CurrentBinIndex() < 0) {
422          AliDebug(AliLog::kDebug + 1, "Current event mixEH->CurrentEntry() == -1");
423          return 0;
424       }
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();
430    }
431
432    return 0;
433 }
434
435 //________________________________________________________________________
436 AliMultiInputEventHandler *AliAnalysisTaskEx02::SetMainInputHandler(AliAnalysisManager *mgr)
437 {
438    //
439    // Sets main input handler
440    //
441
442    if (!fMainInputHandler) fMainInputHandler = dynamic_cast<AliMultiInputEventHandler *>(mgr->GetInputEventHandler());
443
444    return fMainInputHandler;
445 }
446
447 //________________________________________________________________________
448 AliMixInputEventHandler *AliAnalysisTaskEx02::SetMixingInputHandler(AliMultiInputEventHandler *mainIH)
449 {
450    //
451    // Sets mixing input handler
452    //
453
454    if (!fMixingInputHandler) fMixingInputHandler = dynamic_cast<AliMixInputEventHandler *>(mainIH->GetFirstMultiInputHandler());
455
456    return fMixingInputHandler;
457 }
458
459