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