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