]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskLeadingTrackUE.cxx
Fixed dependencies
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskLeadingTrackUE.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 //#include <TVector3.h>
27
28 #include "AliAnalysisTaskLeadingTrackUE.h"
29 #include "AliAnalyseLeadingTrackUE.h"
30 #include "AliUEHistograms.h"
31 #include "AliUEHist.h"
32
33 #include "AliAnalysisHelperJetTasks.h"
34 #include "AliAnalysisManager.h"
35 #include "AliAODHandler.h"
36 #include "AliAODInputHandler.h"
37 #include "AliAODMCParticle.h"
38 #include "AliGenPythiaEventHeader.h"
39 #include "AliInputEventHandler.h"
40 #include "AliLog.h"
41 #include "AliMCEventHandler.h"
42 #include "AliVParticle.h"
43
44
45 ////////////////////////////////////////////////////////////////////////
46 //
47 // Analysis class for Underlying Event studies w.r.t. leading track
48 //
49 // Look for correlations on the tranverse regions w.r.t
50 // the leading track in the event
51 //
52 // This class needs input AODs.
53 // The output is a list of analysis-specific containers.
54 //
55 // The AOD can be either connected to the InputEventHandler  
56 // for a chain of AOD files 
57 // or 
58 // to the OutputEventHandler
59 // for a chain of ESD files,
60 // in this case the class should be in the train after the jet-finder
61 //
62 //    Authors:
63 //    Arian Abrahantes Quintana 
64 //    Jan Fiete Grosse-Oetringhaus
65 //    Ernesto Lopez Torres
66 //    Sara Vallero
67 // 
68 ////////////////////////////////////////////////////////////////////////
69
70
71 ClassImp( AliAnalysisTaskLeadingTrackUE )
72
73 // Define global pointer
74 AliAnalysisTaskLeadingTrackUE* AliAnalysisTaskLeadingTrackUE::fgTaskLeadingTrackUE=NULL;
75
76 //____________________________________________________________________
77 AliAnalysisTaskLeadingTrackUE:: AliAnalysisTaskLeadingTrackUE(const char* name):
78 AliAnalysisTask(name,""),
79 // general configuration
80 fDebug(0),
81 fMode(0),
82 fReduceMemoryFootprint(kFALSE),
83 // pointers to UE classes
84 fAnalyseUE(0x0),
85 fHistosUE(0x0),
86 fkTrackingEfficiency(0x0),
87 // handlers and events
88 fAOD(0x0),           
89 fArrayMC(0x0),
90 fInputHandler(0x0),
91 fMcEvent(0x0),
92 fMcHandler(0x0),
93 // histogram settings
94 fListOfHistos(0x0), 
95 fBinsPtInHist(30),     
96 fMinJetPtInHist(0.),
97 fMaxJetPtInHist(300.), 
98 // event QA
99 fnTracksVertex(1),  // QA tracks pointing to principal vertex (= 3 default) 
100 fZVertex(10.),
101 // track cuts
102 fTrackEtaCut(0.8),
103 fLeadingTrackEtaCut(0.8),
104 fFilterBit(0xFF),
105 fSelectBit(0),
106 fUseChargeHadrons(kFALSE),
107 //For MC weighting
108 fAvgTrials(1)
109 {
110   // Default constructor
111   // Define input and output slots here
112   // Input slot #0 works with a TChain
113   DefineInput(0, TChain::Class());
114   // Output slot #0 writes into a TList container
115   DefineOutput(0, TList::Class());
116
117 }
118
119 AliAnalysisTaskLeadingTrackUE::~AliAnalysisTaskLeadingTrackUE() 
120
121   // destructor
122   
123   if (fListOfHistos  && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) 
124     delete fListOfHistos;
125 }
126
127 /************** INTERFACE METHODS *****************************/
128
129 //______________________________________________________________
130 Bool_t AliAnalysisTaskLeadingTrackUE::Notify()
131 {
132   
133   // Implemented Notify() to read the cross sections
134   // and number of trials from pyxsec.root.
135   // This will be used when merging different MC samples.
136   // (Copied from AliAnalysisTaskJFSystematics)
137   
138   fAvgTrials = 1;
139   TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
140   Float_t xsection = 0;
141   Float_t trials  = 1;
142   if(tree){
143         TFile *curfile = tree->GetCurrentFile();
144         if (!curfile) {
145                 Error("Notify","No current file");
146                 return kFALSE;
147                 }
148         
149         AliAnalysisHelperJetTasks::PythiaInfoFromFile(curfile->GetName(),xsection,trials); 
150         
151         //TO-DO
152         //fHistosUE->GetXsec()->Fill("<#sigma>",xsection);
153
154         // construct average trials
155         Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
156         if(trials>=nEntries && nEntries>0.)fAvgTrials = trials/nEntries;
157         }
158   
159   return kTRUE;
160 }
161
162 //____________________________________________________________________
163 void AliAnalysisTaskLeadingTrackUE::ConnectInputData(Option_t* /*option*/)
164 {
165   
166   // Connect the input data  
167   if (fDebug > 1) AliInfo("ConnectInputData() ");
168   
169   // Since AODs can either be connected to the InputEventHandler
170   // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
171   // we need to get the pointer to the AODEvent correctly.
172   
173   // Delta AODs are also accepted.
174   
175   TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
176   
177   if( handler && handler->InheritsFrom("AliAODInputHandler") ) { // input AOD
178         fAOD = ((AliAODInputHandler*)handler)->GetEvent();
179         if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
180   } else {  //output AOD
181         handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
182         if( handler && handler->InheritsFrom("AliAODHandler") ) {
183                 fAOD = ((AliAODHandler*)handler)->GetAOD();
184                 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
185         } else {  // no AOD
186                 AliFatal("I can't get any AOD Event Handler");
187                 return;
188                 }
189         }       
190   
191   // Initialize common pointers
192   Initialize();
193    
194 }
195
196 //____________________________________________________________________
197 void  AliAnalysisTaskLeadingTrackUE::CreateOutputObjects()
198 {
199   // Create the output container
200   
201   if (fDebug > 1) AliInfo("CreateOutputObjects()");
202    
203   // Initialize class with main algorithms, event and track selection. 
204   fAnalyseUE = new AliAnalyseLeadingTrackUE();
205   fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fLeadingTrackEtaCut);
206   fAnalyseUE->SetDebug(fDebug); 
207
208   // Initialize output list of containers
209   if (fListOfHistos != NULL){
210         delete fListOfHistos;
211         fListOfHistos = NULL;
212         }
213   if (!fListOfHistos){
214         fListOfHistos = new TList();
215         fListOfHistos->SetOwner(kTRUE); 
216         }
217
218   // Initialize class to handle histograms 
219   fHistosUE = new AliUEHistograms("AliUEHistograms", "123");
220   
221   // add histograms to list
222   fListOfHistos->Add(fHistosUE);
223   
224   //fListOfHistos->Add(new TH2F("multVsLeadStep5", ";multiplicity;leading pT", 100, -0.5, 99.5, 20, 0, 10));
225   
226   // Add task configuration to output list 
227   AddSettingsTree();
228   
229
230   PostData(0,fListOfHistos);
231 }
232
233 //____________________________________________________________________
234 void  AliAnalysisTaskLeadingTrackUE::Exec(Option_t */*option*/)
235 {
236   // array of MC particles
237   if (fMcHandler){
238         fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
239         if (!fArrayMC)AliFatal("No array of MC particles found !!!");
240         }
241
242   // Get number of trials from MC header
243   
244 //   Float_t nTrials = 1;
245   if (fMcHandler) {
246         fMcEvent = fMcHandler->MCEvent();
247 //      if (fMcEvent) {
248                 // TO-DO: extend to PHOJET
249 //              AliGenPythiaEventHeader*  pythiaGenHeader = AliAnalysisHelperJetTasks::GetPythiaEventHeader(fMcEvent);
250 //              if(pythiaGenHeader){
251 //                      nTrials = pythiaGenHeader->Trials();
252 //                      }
253 //              }
254         }
255
256   // TO-DO      
257   //fHistosUE->GetTrials()->Fill("#sum{ntrials}",fAvgTrials);
258   
259   // Analyse the event
260   if (fMode) AnalyseCorrectionMode();
261   else AnalyseDataMode();
262
263   PostData(0,fListOfHistos);
264 }
265
266 //____________________________________________________________________
267 void  AliAnalysisTaskLeadingTrackUE::Terminate(Option_t */*option*/)
268 {
269   
270   // Terminate analysis
271   if( fDebug > 1 ) AliInfo("End analysis");
272
273   if (!gROOT->IsBatch()){
274         fListOfHistos = dynamic_cast<TList*> (GetOutputData(0));
275         if (!fListOfHistos){
276                 AliError("Histogram List is not available");
277                 return;
278         }else{
279         // Draw something
280         }
281   } else {
282         AliInfo(" Batch mode, not histograms will be shown...");
283         }
284
285 }
286
287
288 /******************** ANALYSIS METHODS *****************************/
289
290 //____________________________________________________________________
291 void  AliAnalysisTaskLeadingTrackUE::AddSettingsTree()
292 {
293   //Write settings to output list
294   TTree *settingsTree   = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
295   settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
296   settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
297   settingsTree->Branch("fLeadingTrackEtaCut", &fLeadingTrackEtaCut, "LeadingTrackEtaCut/D");
298   settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
299   settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
300   settingsTree->Fill();
301   fListOfHistos->Add(settingsTree);
302 }  
303
304 //____________________________________________________________________
305 void  AliAnalysisTaskLeadingTrackUE::AnalyseCorrectionMode()
306 {
307   // Run the analysis on MC to get the correction maps
308   //
309   // if fReduceMemoryFootprint step 3,4,5,7,9 are not filled
310
311   PostData(0,fListOfHistos);
312   
313   if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
314   
315   //PROCESS TYPE (ND,SD,DD)
316   AliAnalysisHelperJetTasks::MCProcessType eventId = AliAnalysisHelperJetTasks::kInvalidProcess;
317   AliGenEventHeader* genHeader = fMcEvent->GenEventHeader();
318   eventId = AliAnalysisHelperJetTasks::GetPythiaEventProcessType(genHeader,kFALSE);
319   if (eventId<0)
320     eventId = AliAnalysisHelperJetTasks::GetDPMjetEventProcessType(genHeader,kFALSE);
321   if (eventId<0 && fDebug>1)
322     AliInfo("No Pythia or Phojet Header retrived!");
323     
324   Int_t fillId=-1;
325   if (eventId == AliAnalysisHelperJetTasks::kND)fillId = 0; 
326   if (eventId == AliAnalysisHelperJetTasks::kSD)fillId = 1; 
327   if (eventId == AliAnalysisHelperJetTasks::kDD)fillId = 2; 
328   
329   // count all events
330   fHistosUE->FillEvent(fillId, -1);
331   
332   // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
333   if (!fAnalyseUE->VertexSelection(fMcEvent, 0, fZVertex)) 
334     return;
335   
336   // Get MC-true leading particle (but do not cut out events!)
337   TObjArray *ltMC = (TObjArray*)fAnalyseUE->FindLeadingObjects(fArrayMC);
338   AliVParticle* leadingMC = 0;
339   if (ltMC)
340     leadingMC = (AliVParticle*) ltMC->At(0);
341     
342   // it can happen that there is no MC leading particle in the acceptance required (|eta|<0.8)
343   // and we do not want to base the event slection on MC information
344   
345   // Sort MC-true charged particles
346   // as output you get an array of 3 lists of  particles belonging to different regions:
347   // - at 0: towards
348   // - at 1: away
349   // - at 2: transverse MIN
350   // - at 3: transverse MAX
351   TObjArray *regionSortedParticlesMC = (TObjArray*)fAnalyseUE->SortRegions(leadingMC, fArrayMC, 0x0); 
352   TObjArray *regionsMinMaxMC = (TObjArray*)fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesMC->At(2),(TList*)regionSortedParticlesMC->At(3));
353   // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
354   // (MC-true leading particle and MC-true all particles)
355   // STEP 0
356   fHistosUE->Fill(fillId,0,AliUEHist::kCFStepAll,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
357   
358   // Trigger selection ************************************************
359   if (fAnalyseUE->TriggerSelection(fInputHandler))
360   {
361     // PILEUP-CUT 
362     Bool_t select = kFALSE;
363     if (fSelectBit) select = AliAnalysisHelperJetTasks::TestSelectInfo(fSelectBit);
364     if (select)
365       fHistosUE->FillEvent(fillId, -2);
366     else
367     {
368     
369             // Count events that pass AliPhysicsSelection
370     
371             // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
372             // (MC-true leading particle and MC-true all particles)
373             // STEP 1
374             fHistosUE->Fill(fillId,0,AliUEHist::kCFStepTriggered,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
375   
376             // count number of MC tracks above 150 MeV/c
377             Int_t nMCTracks = 0; 
378             if (leadingMC && leadingMC->Pt() > 0.15)
379               nMCTracks++;
380             for (Int_t i=0; i<4; i++)
381               for (Int_t j=0; j<((TList*)regionSortedParticlesMC->At(i))->GetEntries(); j++)
382                 if (((AliVParticle*) ((TList*)regionSortedParticlesMC->At(i))->At(j))->Pt() > 0.15)
383                   nMCTracks++;
384       
385             //((TH2F*)fListOfHistos->FindObject("multVsLeadStep5"))->Fill(nMCTracks, leadingMC->Pt());
386     
387             // Vertex selection *************************************************
388             if (fAnalyseUE->VertexSelection(fAOD, fnTracksVertex, fZVertex))
389             {
390               // Count events that pass Vertex selection
391             
392               // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
393               // (MC-true leading particle and MC-true all particles)
394               // STEP 2
395               fHistosUE->Fill(fillId,0,AliUEHist::kCFStepVertex,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
396     
397               // fill here for tracking efficiency
398               // loop over particle species
399               for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
400               {
401                 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(fArrayMC, 0x0, kTRUE, particleSpecies);
402                 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(fAOD, fArrayMC, kTRUE, particleSpecies);
403                 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(fAOD, fArrayMC, kFALSE, particleSpecies);
404
405                 fHistosUE->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, 0, 0, particleSpecies);
406               
407                 delete primMCParticles;
408                 delete primRecoTracksMatched;
409                 delete allRecoTracksMatched;
410               }
411       
412               // Get Reconstructed leading particle *******************************
413               TObjArray *ltRECO = fAnalyseUE->FindLeadingObjects(fAOD);
414               if (ltRECO)
415               {
416                 // Count events where a reconstructed track was found in |eta|<0.8
417                 // the pT cut will be set when projecting output containers
418                 // for leading particle correlation plots
419                 if (leadingMC) {
420                       fHistosUE->Fill(leadingMC, (AliVParticle*)ltRECO->At(0));
421                       }
422               
423                 // If there is no MC leading track the container is not filled, so the number of entries in the container might be different  
424                 // from the number of events after the selection, since the selection is based on RECO tracks
425                 // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
426                 // (MC-true leading particle and MC-true all particles)
427                 // STEP 3
428                 if (!fReduceMemoryFootprint)
429                   fHistosUE->Fill(fillId,0,AliUEHist::kCFStepAnaTopology,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
430       
431                 //Sort RECO particles w.r.t. MC-leading and return matched (primary) MC particle 
432                 // (you cannot sort tracks w.r.t. RECO-leading and plot it vs. MC-leading ...)
433                 TObjArray *regionSortedParticlesRECOLTMC = (TObjArray*)fAnalyseUE->SortRegions(leadingMC, fAOD, fArrayMC, kTRUE); 
434                 TObjArray *regionsMinMaxRECOLTMC = (TObjArray*)fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECOLTMC->At(2),(TList*)regionSortedParticlesRECOLTMC->At(3));
435                 // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
436                 // (MC leading particle and RECO-matched (quantities from MC particle)  all particles)
437                 // STEP 4
438                 //if (!fReduceMemoryFootprint)
439                   fHistosUE->Fill(fillId,0,AliUEHist::kCFStepTrackedOnlyPrim,leadingMC,(TList*)regionSortedParticlesRECOLTMC->At(0),(TList*)regionSortedParticlesRECOLTMC->At(1),(TList*)regionsMinMaxRECOLTMC->At(0),(TList*)regionsMinMaxRECOLTMC->At(1));
440                 // comparing this step with step 3 (for all-tracks observables) you get the tracking efficiency
441         
442                 //Sort RECO particles w.r.t. MC-leading and return matched (primary+secondary) MC particle 
443                 // (you cannot sort tracks w.r.t. RECO-leading and plot it vs. MC-leading ...)
444                 TObjArray *regionSortedParticlesRECOLTMC2 = (TObjArray*)fAnalyseUE->SortRegions(leadingMC, fAOD, fArrayMC,kFALSE); 
445                 TObjArray *regionsMinMaxRECOLTMC2 = (TObjArray*)fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECOLTMC2->At(2),(TList*)regionSortedParticlesRECOLTMC2->At(3));
446         
447                 // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
448                 // (MC leading particle and RECO-matched (quantities from MC particle)  all particles)
449                 // STEP 5
450                 //if (!fReduceMemoryFootprint)
451                   fHistosUE->Fill(fillId,0,AliUEHist::kCFStepTracked,leadingMC,(TList*)regionSortedParticlesRECOLTMC2->At(0),(TList*)regionSortedParticlesRECOLTMC2->At(1),(TList*)regionsMinMaxRECOLTMC2->At(0),(TList*)regionsMinMaxRECOLTMC2->At(1));
452                 // comparing this step with step 3 (for all-tracks observables) you get the tracking efficiency
453           
454                 // SWITCH TO RECONSTRUCTED TRACKS  ************************************
455                 // The next steps correspond to track selections
456                 // Sort RECO particles w.r.t. RECO-leading and return RECO particle  
457                 TObjArray *regionSortedParticlesRECO = (TObjArray*)fAnalyseUE->SortRegions((AliVParticle*)ltRECO->At(0), fAOD,0); 
458                 TObjArray *regionsMinMaxRECO = (TObjArray*)fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECO->At(2),(TList*)regionSortedParticlesRECO->At(3));
459                 // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
460                 // (RECO leading particle and RECO  all particles)
461                 // STEP 6
462                 fHistosUE->Fill(fillId,0,AliUEHist::kCFStepReconstructed,(AliVParticle*)ltRECO->At(0),(TList*)regionSortedParticlesRECO->At(0),(TList*)regionSortedParticlesRECO->At(1),(TList*)regionsMinMaxRECO->At(0),(TList*)regionsMinMaxRECO->At(1));
463         
464                 // STEP 8 for reduced efficiency study
465                 FillReducedEfficiency(fillId, AliUEHist::kCFStepBiasStudy, ltRECO, kFALSE);
466                 if (!fReduceMemoryFootprint)
467                   FillReducedEfficiency(fillId, AliUEHist::kCFStepBiasStudy2, ltRECO, kTRUE);
468         
469                 // count number of reco tracks above 150 MeV/c
470                 Int_t nRecoTracks = 0; 
471                 if (((AliVParticle*) ltRECO->At(0))->Pt() > 0.15)
472                   nRecoTracks++;
473                 for (Int_t i=0; i<4; i++)
474                   for (Int_t j=0; j<((TList*)regionSortedParticlesRECO->At(i))->GetEntries(); j++)
475                     if (((AliVParticle*) ((TList*)regionSortedParticlesRECO->At(i))->At(j))->Pt() > 0.15)
476                       nRecoTracks++;
477                 
478                 if (leadingMC && leadingMC->Pt() > 0.5)
479                   fHistosUE->GetCorrelationMultiplicity()->Fill(nMCTracks, nRecoTracks);
480         
481                 if (leadingMC)
482                 {
483                   // Match reco leading track with true *********************************
484                   Int_t recoLabel = ((AliAODTrack*)ltRECO->At(0))->GetLabel();
485                   Int_t mcLabel   = ((AliAODMCParticle*)leadingMC)->GetLabel();
486                   if (recoLabel != mcLabel) 
487                     return;
488           
489                   // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
490                   // (RECO-MATCHED leading particle and RECO all particles)
491                   // STEP 7
492                   fHistosUE->Fill(fillId,0,AliUEHist::kCFStepRealLeading,(AliVParticle*)ltRECO->At(0),(TList*)regionSortedParticlesRECO->At(0),(TList*)regionSortedParticlesRECO->At(1),(TList*)regionsMinMaxRECO->At(0),(TList*)regionsMinMaxRECO->At(1));
493                   // comparing this step with step 6 (for leading-track observables) you get the efficiency to reconstruct the leading track
494                   // comparing this step with step 6 (for all-tracks observables) you see how leading-track misidentification affects the final distributions
495                 }
496         
497                 delete regionSortedParticlesRECOLTMC;
498                 delete regionsMinMaxRECOLTMC;
499                 delete regionSortedParticlesRECOLTMC2;
500                 delete regionsMinMaxRECOLTMC2;
501                 delete regionSortedParticlesRECO;
502                 delete regionsMinMaxRECO;
503                 delete ltRECO;
504               } // lt reco
505         } // vertex
506     } // pileup
507   } //phyiscs selection
508   
509   if (ltMC)
510     delete ltMC;
511   delete regionSortedParticlesMC;
512   delete regionsMinMaxMC;
513 }
514
515 //____________________________________________________________________
516 void  AliAnalysisTaskLeadingTrackUE::AnalyseDataMode()
517 {
518
519   // Run the analysis on DATA or MC to get raw distributions
520  
521   PostData(0,fListOfHistos);
522   
523   if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
524   Int_t eventId = 0;
525
526   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
527   fHistosUE->FillEvent(eventId, AliUEHist::kCFStepAll);
528
529   // Trigger selection ************************************************
530   if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
531   // PILEUP-CUT 
532   Bool_t select = kFALSE;
533   if (fSelectBit) select = AliAnalysisHelperJetTasks::TestSelectInfo(fSelectBit);
534   if (select) 
535   {
536     fHistosUE->FillEvent(eventId, -2);
537     return;
538   }
539   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
540   fHistosUE->FillEvent(eventId, AliUEHist::kCFStepTriggered);
541   
542   // Vertex selection *************************************************
543   if(!fAnalyseUE->VertexSelection(fAOD, fnTracksVertex, fZVertex)) return;
544   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
545   fHistosUE->FillEvent(eventId, AliUEHist::kCFStepVertex);
546   // comparing this step with previous one you get the vertex selection efficiency from data (?)
547  
548
549   // Get Reconstructed leading particle *******************************
550   TObjArray *ltRECO = fAnalyseUE->FindLeadingObjects(fAOD);
551   if (!ltRECO) return;
552   
553   // fill control distributions
554   fHistosUE->Fill(0, (AliVParticle*)ltRECO->At(0));
555   
556   // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
557   fHistosUE->FillEvent(eventId, AliUEHist::kCFStepAnaTopology);
558
559   // Switch to reconstructed tracks ************************************
560   // Sort RECO particles w.r.t. RECO-leading and return RECO particle  
561   TObjArray *regionSortedParticlesRECO = (TObjArray*)fAnalyseUE->SortRegions((AliVParticle*)ltRECO->At(0), fAOD,0); 
562   TObjArray *regionsMinMaxRECO = (TObjArray*)fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECO->At(2),(TList*)regionSortedParticlesRECO->At(3));
563   // Fill UE containers (step, leading track, towards particles, away particles, transverse MIN and MAX particles)
564   // (RECO leading particle and RECO  all particles)
565   // STEP 6
566   fHistosUE->Fill(eventId,0,AliUEHist::kCFStepReconstructed,(AliVParticle*)ltRECO->At(0),(TList*)regionSortedParticlesRECO->At(0),(TList*)regionSortedParticlesRECO->At(1),(TList*)regionsMinMaxRECO->At(0),(TList*)regionsMinMaxRECO->At(1));
567   
568   // STEP 8 and 9 for reduced efficiency study
569   FillReducedEfficiency(eventId, AliUEHist::kCFStepBiasStudy, ltRECO, kFALSE);
570   FillReducedEfficiency(eventId, AliUEHist::kCFStepBiasStudy2, ltRECO, kTRUE);
571   
572   delete ltRECO;
573   delete regionSortedParticlesRECO;
574   delete regionsMinMaxRECO;
575 }
576
577 //____________________________________________________________________
578 void AliAnalysisTaskLeadingTrackUE::FillReducedEfficiency(Int_t eventId, AliUEHist::CFStep step, const TObjArray* ltRECO, Bool_t twoStep)
579 {
580   // remove leading particle using fkTrackingEfficiency and use subleading particle to fill the histograms
581   //
582   // if twoStep is kTRUE, do a two step procedure where in each step only 50% of the loss due to the tracking efficiency is applied 
583   
584   if (!fkTrackingEfficiency)
585     return;
586     
587   TObjArray* particleList =  new TObjArray(*ltRECO);
588   AliVParticle* leading = (AliVParticle*) particleList->At(0);
589   if (!leading)
590   {
591     delete particleList;
592     return;
593   }
594   
595   // remove particles depending on tracking efficiency
596   Int_t count = (twoStep) ? 2 : 1;
597   
598   for (Int_t i=0; i<count; i++)
599   {
600     const TAxis * xax = fkTrackingEfficiency->GetXaxis();
601     Float_t trackingEff = fkTrackingEfficiency->GetBinContent(xax->FindFixBin(leading->Pt()));
602     if (twoStep)
603       trackingEff = 0.5 * (trackingEff + 1);
604       
605     if (gRandom->Uniform() > trackingEff)
606     {
607       //Printf("LOWEFF: Removing leading particle");
608       particleList->RemoveAt(0);
609       particleList->Compress();
610     }
611       
612     if (particleList->GetEntries() == 0)
613     {
614       delete particleList;
615       return;
616     }
617     
618     leading = (AliVParticle*) particleList->At(0);
619   }
620   
621   TObjArray *regionSortedParticlesRECOLowEff = fAnalyseUE->SortRegions(leading, particleList, 0); 
622   TObjArray *regionsMinMaxRECOLowEff = fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECOLowEff->At(2), (TList*)regionSortedParticlesRECOLowEff->At(3));
623     
624   fHistosUE->Fill(eventId,0,step,leading,(TList*)regionSortedParticlesRECOLowEff->At(0), (TList*)regionSortedParticlesRECOLowEff->At(1), (TList*)regionsMinMaxRECOLowEff->At(0), (TList*)regionsMinMaxRECOLowEff->At(1));
625
626   delete regionSortedParticlesRECOLowEff;
627   delete regionsMinMaxRECOLowEff;
628   delete particleList;
629 }
630
631 //____________________________________________________________________
632 void  AliAnalysisTaskLeadingTrackUE::Initialize()
633 {
634   // input handler
635   fInputHandler = (AliInputEventHandler*)
636          ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
637   // MC handler
638   fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
639 }
640
641
642
643
644
645