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