]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalysisTaskPhiCorrelations.cxx
samll fixes (Marta)
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalysisTaskPhiCorrelations.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:$ */
17
18 #include <TROOT.h>
19 #include <TChain.h>
20 #include <TFile.h>
21 #include <TList.h>
22 #include <TMath.h>
23 #include <TTree.h>
24 #include <TH2F.h>
25 #include <TRandom.h>
26 #include <TMCProcess.h>
27
28 #include "AliAnalysisTaskPhiCorrelations.h"
29 #include "AliAnalyseLeadingTrackUE.h"
30 #include "AliUEHistograms.h"
31 #include "AliUEHist.h"
32
33 #include "AliAnalysisHelperJetTasks.h"
34 #include "AliAnalysisManager.h"
35 #include "AliAODHandler.h"
36 #include "AliAODInputHandler.h"
37 #include "AliAODMCParticle.h"
38 #include "AliGenPythiaEventHeader.h"
39 #include "AliInputEventHandler.h"
40 #include "AliLog.h"
41 #include "AliMCEventHandler.h"
42 #include "AliVParticle.h"
43 #include "AliCFContainer.h"
44
45 #include "AliESDEvent.h"
46 #include "AliESDInputHandler.h"
47 #include "AliMultiplicity.h"
48 #include "AliCentrality.h"
49 #include "AliStack.h"
50
51 #include "AliEventPoolManager.h"
52
53
54 ////////////////////////////////////////////////////////////////////////
55 //
56 // Analysis class for azimuthal correlation studies
57 // Based on the UE task from Sara Vallero and Jan Fiete
58 //
59 // This class needs input AODs.
60 // The output is a list of analysis-specific containers.
61 //
62 // The AOD can be either connected to the InputEventHandler  
63 // for a chain of AOD files 
64 // or 
65 // to the OutputEventHandler
66 // for a chain of ESD files,
67 // in this case the class should be in the train after the jet-finder
68 //
69 //    Authors:
70 //    Jan Fiete Grosse-Oetringhaus
71 // 
72 ////////////////////////////////////////////////////////////////////////
73
74
75 ClassImp( AliAnalysisTaskPhiCorrelations )
76
77 //____________________________________________________________________
78 AliAnalysisTaskPhiCorrelations:: AliAnalysisTaskPhiCorrelations(const char* name):
79 AliAnalysisTask(name,""),
80 // general configuration
81 fDebug(0),
82 fMode(0),
83 fReduceMemoryFootprint(kFALSE),
84 fFillMixed(kTRUE),
85 fCompareCentralities(kFALSE),
86 // pointers to UE classes
87 fAnalyseUE(0x0),
88 fHistos(0x0),
89 fHistosMixed(0),
90 fkTrackingEfficiency(0x0),
91 // handlers and events
92 fAOD(0x0),
93 fESD(0x0),
94 fArrayMC(0x0),
95 fInputHandler(0x0),
96 fMcEvent(0x0),
97 fMcHandler(0x0),
98 fPoolMgr(0x0),
99 // histogram settings
100 fListOfHistos(0x0), 
101 // event QA
102 fnTracksVertex(1),  // QA tracks pointing to principal vertex (= 3 default) 
103 fZVertex(7.),
104 fCentralityMethod("V0M"),
105 // track cuts
106 fTrackEtaCut(0.8),
107 fPtMin(0.5),
108 fFilterBit(0xFF),
109 fSelectBit(0),
110 fUseChargeHadrons(kFALSE),
111 fSelectCharge(0)
112 {
113   // Default constructor
114   // Define input and output slots here
115   // Input slot #0 works with a TChain
116   DefineInput(0, TChain::Class());
117   // Output slot #0 writes into a TList container
118   DefineOutput(0, TList::Class());
119 }
120
121 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations() 
122
123   // destructor
124   
125   if (fListOfHistos  && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) 
126     delete fListOfHistos;
127 }
128
129 //____________________________________________________________________
130 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
131 {
132   
133   // Connect the input data  
134   if (fDebug > 1) AliInfo("ConnectInputData() ");
135   
136   // Since AODs can either be connected to the InputEventHandler
137   // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
138   // we need to get the pointer to the AODEvent correctly.
139   
140   // Delta AODs are also accepted.
141   
142   TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
143   
144   if( handler && handler->InheritsFrom("AliAODInputHandler") ) 
145   { // input AOD
146     fAOD = ((AliAODInputHandler*)handler)->GetEvent();
147     if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
148   } 
149   else 
150   {  //output AOD
151     handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
152     if (handler && handler->InheritsFrom("AliAODHandler") ) 
153     {
154       fAOD = ((AliAODHandler*)handler)->GetAOD();
155       if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
156     } 
157     else 
158     {  // no AOD
159       AliWarning("I can't get any AOD Event Handler");
160     }
161   }
162   
163   if (handler && handler->InheritsFrom("AliESDInputHandler") ) 
164   { // input ESD
165     // pointer received per event in ::Exec
166     if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
167   } 
168   
169   // Initialize common pointers
170   Initialize();
171 }
172
173 //____________________________________________________________________
174 void  AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
175 {
176   // Create the output container
177   
178   if (fDebug > 1) AliInfo("CreateOutputObjects()");
179    
180   // Initialize class with main algorithms, event and track selection. 
181   fAnalyseUE = new AliAnalyseLeadingTrackUE();
182   fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
183   fAnalyseUE->SetDebug(fDebug); 
184   fAnalyseUE->DefineESDCuts(fFilterBit);
185
186   // Initialize output list of containers
187   if (fListOfHistos != NULL){
188         delete fListOfHistos;
189         fListOfHistos = NULL;
190         }
191   if (!fListOfHistos){
192         fListOfHistos = new TList();
193         fListOfHistos->SetOwner(kTRUE); 
194         }
195
196   // Initialize class to handle histograms 
197   fHistos = new AliUEHistograms("AliUEHistogramsSame", "4");
198   fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", "4");
199   
200   fHistos->SetSelectCharge(fSelectCharge);
201   fHistosMixed->SetSelectCharge(fSelectCharge);
202   
203   // add histograms to list
204   fListOfHistos->Add(fHistos);
205   fListOfHistos->Add(fHistosMixed);
206   
207   fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
208   fListOfHistos->Add(new TH2F("processIDs", ";#Delta#phi;process id", 100, -0.5 * TMath::Pi(), 1.5 * TMath::Pi(), kPNoProcess + 1, -0.5, kPNoProcess + 0.5));
209   
210   PostData(0,fListOfHistos);
211   
212   // Add task configuration to output list 
213   AddSettingsTree();
214
215   // event mixing
216   //  Int_t trackDepth = 100; // Require e.g. 20 5-track events, or 2 50-track events
217   Int_t trackDepth = 5000; 
218   Int_t poolsize   = 100;  // Maximum number of events
219   
220   Int_t nCentralityBins  = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
221   Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
222   
223   Int_t nZvtxBins  = 7;
224   Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7 };
225   Double_t* zvtxbin = vertexBins;
226   
227   if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
228   {
229     nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
230     zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
231   }
232
233   fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
234 }
235
236 //____________________________________________________________________
237 void  AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
238 {
239   // array of MC particles
240   if (fMcHandler) {
241     if (fAOD)
242     {
243       fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
244       if (!fArrayMC)
245         AliFatal("No array of MC particles found !!!");
246     }
247     fMcEvent = fMcHandler->MCEvent();
248   }
249
250   // receive ESD pointer if we are not running AOD analysis
251   if (!fAOD)
252   {
253     AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
254     if (handler && handler->InheritsFrom("AliESDInputHandler"))
255       fESD = ((AliESDInputHandler*)handler)->GetEvent();
256   }
257
258   // Analyse the event
259   if (fMode) AnalyseCorrectionMode();
260   else AnalyseDataMode();
261 }
262
263 /******************** ANALYSIS METHODS *****************************/
264
265 //____________________________________________________________________
266 void  AliAnalysisTaskPhiCorrelations::AddSettingsTree()
267 {
268   //Write settings to output list
269   TTree *settingsTree   = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
270   settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
271   settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
272   settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
273   settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
274   settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
275   settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
276   settingsTree->Fill();
277   fListOfHistos->Add(settingsTree);
278 }  
279
280 //____________________________________________________________________
281 void  AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
282 {
283   // Run the analysis on MC to get the correction maps
284   //
285  
286   if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
287   
288   Double_t centrality = 0;
289   
290   if (fCentralityMethod.Length() > 0)
291   {
292     AliCentrality *centralityObj = 0;
293     if (fAOD)
294       centralityObj = fAOD->GetHeader()->GetCentralityP();
295     else if (fESD)
296       centralityObj = fESD->GetCentrality();
297     
298     if (centralityObj)
299     {
300       centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
301       AliInfo(Form("Centrality is %f", centrality));
302     }
303     else
304     {
305       Printf("WARNING: Centrality object is 0");
306       centrality = -1;
307      }
308   }
309   
310   // Support for ESD and AOD based analysis
311   AliVEvent* inputEvent = fAOD;
312   if (!inputEvent)
313     inputEvent = fESD;
314   
315   TObject* mc = fArrayMC;
316   if (!mc)
317     mc = fMcEvent;
318   
319   // count all events
320   fHistos->FillEvent(centrality, -1);
321   
322   // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
323   if (!fAnalyseUE->VertexSelection(fMcEvent, 0, fZVertex)) 
324     return;
325   
326   Float_t zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
327     
328   // Get MC primaries
329   TObjArray* tracksMC = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
330   
331   // (MC-true all particles)
332   // STEP 0
333   fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC);
334   
335   // Trigger selection ************************************************
336   if (fAnalyseUE->TriggerSelection(fInputHandler))
337   {  
338     // (MC-true all particles)
339     // STEP 1
340     if (!fReduceMemoryFootprint)
341       fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC);
342     else
343       fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
344       
345     if (!inputEvent) {
346       AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
347       return;
348     }
349     
350     // Vertex selection *************************************************
351     if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
352     {
353       // fill here for tracking efficiency
354       // loop over particle species
355       
356       for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
357       {
358         TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
359         TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
360         TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
361       
362         fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies, centrality);
363         
364 //      Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
365
366         delete primMCParticles;
367         delete primRecoTracksMatched;
368         delete allRecoTracksMatched;
369       }
370     
371       // (MC-true all particles)
372       // STEP 2
373       if (!fReduceMemoryFootprint)
374         fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC);
375       else
376         fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
377       
378       // Get MC primaries that match reconstructed track
379       TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
380       
381       // (RECO-matched (quantities from MC particle) primary particles)
382       // STEP 4
383       fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim);
384       
385       // Get MC primaries + secondaries that match reconstructed track
386       TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
387       
388       // (RECO-matched (quantities from MC particle) all particles)
389       // STEP 5
390       fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll);
391       
392       // Get RECO tracks
393       TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
394       
395       // (RECO all tracks)
396       // STEP 6
397       fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks);
398       
399       if (0 && !fReduceMemoryFootprint)
400       {
401         // make list of secondaries (matched with MC)
402         TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
403         for (Int_t i=0; i<tracksRecoMatchedAll->GetEntries(); i++)
404           if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
405             tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
406       
407         // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
408         fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll);
409         
410         // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
411         fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries);
412       
413         // plot delta phi vs process id of secondaries
414         // trigger particles: primaries in 4 < pT < 10
415         // associated particles: secondaries in 1 < pT < 10
416         
417         for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntries(); i++)
418         {
419           AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
420           
421           if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
422             continue;
423           
424           for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntries(); j++)
425           {
426             AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
427             
428             if (particle->Pt() < 1 || particle->Pt() > 10)
429               continue;
430             
431             if (particle->Pt() > triggerParticle->Pt())
432               continue;
433               
434             Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
435             if (deltaPhi > 1.5 * TMath::Pi()) 
436               deltaPhi -= TMath::TwoPi();
437             if (deltaPhi < -0.5 * TMath::Pi())
438               deltaPhi += TMath::TwoPi();
439               
440             Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
441               
442             ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
443           }
444         }
445         
446         delete tracksRecoMatchedSecondaries;
447       }
448   
449       delete tracksRecoMatchedPrim;
450       delete tracksRecoMatchedAll;
451       delete tracks;
452     }
453   }
454   
455   delete tracksMC;
456 }
457
458 //____________________________________________________________________
459 void  AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
460 {
461
462   // Run the analysis on DATA or MC to get raw distributions
463  
464   if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
465
466   // skip not selected events here (the AOD is not updated for those)
467   if (!fInputHandler)
468     return;
469     
470   if (!(fInputHandler->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kUserDefined)))
471     return;
472
473   Double_t centrality = 0;
474   
475   AliCentrality *centralityObj = 0;
476   if (fCentralityMethod.Length() > 0)
477   {
478     if (fAOD)
479       centralityObj = fAOD->GetHeader()->GetCentralityP();
480     else if (fESD)
481       centralityObj = fESD->GetCentrality();
482     
483     if (centralityObj)
484       centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
485       //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
486     else
487       centrality = -1;
488     AliInfo(Form("Centrality is %f", centrality));
489   }
490   
491   // Support for ESD and AOD based analysis
492   AliVEvent* inputEvent = fAOD;
493   if (!inputEvent)
494     inputEvent = fESD;
495
496   fHistos->SetRunNumber(inputEvent->GetRunNumber());
497   
498   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
499   fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
500   
501   // Trigger selection ************************************************
502   if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
503   
504   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
505   fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
506   
507   // Vertex selection *************************************************
508   if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
509   
510   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
511   fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
512  
513   // optimization
514   if (centrality < 0 && !fCompareCentralities)
515     return;
516
517   TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
518   //Printf("Accepted %d tracks", tracks->GetEntries());
519   
520   const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
521   Double_t zVtx = vertex->GetZ();
522   
523   // fill two different centralities (syst study)
524   // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
525   if (fCompareCentralities)
526     zVtx = 0;
527   
528   // Fill containers at STEP 6 (reconstructed)
529   if (centrality >= 0)
530     fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks);
531
532   // fill second time with SPD centrality
533   if (fCompareCentralities && centralityObj)
534   {
535     centrality = centralityObj->GetCentralityPercentile("CL1");
536     if (centrality >= 0)
537       fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracks);
538   }
539     
540   if (fFillMixed)
541   {
542     // event mixing
543     
544     // 1. First get an event pool corresponding in mult (cent) and
545     //    zvertex to the current event. Once initialized, the pool
546     //    should contain nMix (reduced) events. This routine does not
547     //    pre-scan the chain. The first several events of every chain
548     //    will be skipped until the needed pools are filled to the
549     //    specified depth. If the pool categories are not too rare, this
550     //    should not be a problem. If they are rare, you could lose
551     //    statistics.
552
553     // 2. Collect the whole pool's content of tracks into one TObjArray
554     //    (bgTracks), which is effectively a single background super-event.
555
556     // 3. The reduced and bgTracks arrays must both be passed into
557     //    FillCorrelations(). Also nMix should be passed in, so a weight
558     //    of 1./nMix can be applied.
559
560     AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
561     
562     if (!pool)
563       AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
564     
565     //pool->SetDebug(1);
566       
567     if (pool->IsReady()) 
568     {
569       
570       Int_t nMix = pool->GetCurrentNEvents();
571       //cout << "nMix = " << nMix << endl;
572     
573       // Fill mixed-event histos here  
574       for (Int_t jMix=0; jMix<nMix; jMix++) {
575         TObjArray* bgTracks = pool->GetEvent(jMix);
576         fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, bgTracks, 1.0 / nMix, (jMix == 0));
577       }
578     }
579     
580     // clone and give ownership to event pool
581     TObjArray* tracksClone = (TObjArray*) tracks->Clone();
582     tracksClone->SetOwner(kTRUE);
583     
584     pool->UpdatePool(tracksClone);
585     //pool->PrintInfo();
586   }
587
588   delete tracks;
589 }
590
591 //____________________________________________________________________
592 void  AliAnalysisTaskPhiCorrelations::Initialize()
593 {
594   // input handler
595   fInputHandler = (AliInputEventHandler*)
596          ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
597   // MC handler
598   fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
599 }