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