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