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