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