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