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