1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
26 //#include <TVector3.h>
28 #include "AliAnalysisTaskLeadingTrackUE.h"
29 #include "AliAnalyseLeadingTrackUE.h"
30 #include "AliUEHistograms.h"
31 #include "AliUEHist.h"
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"
41 #include "AliMCEventHandler.h"
42 #include "AliVParticle.h"
45 ////////////////////////////////////////////////////////////////////////
47 // Analysis class for Underlying Event studies w.r.t. leading track
49 // Look for correlations on the tranverse regions w.r.t
50 // the leading track in the event
52 // This class needs input AODs.
53 // The output is a list of analysis-specific containers.
55 // The AOD can be either connected to the InputEventHandler
56 // for a chain of AOD files
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
63 // Arian Abrahantes Quintana
64 // Jan Fiete Grosse-Oetringhaus
65 // Ernesto Lopez Torres
68 ////////////////////////////////////////////////////////////////////////
71 ClassImp( AliAnalysisTaskLeadingTrackUE )
73 // Define global pointer
74 AliAnalysisTaskLeadingTrackUE* AliAnalysisTaskLeadingTrackUE::fgTaskLeadingTrackUE=NULL;
76 //____________________________________________________________________
77 AliAnalysisTaskLeadingTrackUE:: AliAnalysisTaskLeadingTrackUE(const char* name):
78 AliAnalysisTask(name,""),
79 // general configuration
82 fReduceMemoryFootprint(kFALSE),
83 // pointers to UE classes
86 fkTrackingEfficiency(0x0),
87 // handlers and events
97 fMaxJetPtInHist(300.),
99 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
103 fLeadingTrackEtaCut(0.8),
106 fUseChargeHadrons(kFALSE),
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());
119 AliAnalysisTaskLeadingTrackUE::~AliAnalysisTaskLeadingTrackUE()
123 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
124 delete fListOfHistos;
127 /************** INTERFACE METHODS *****************************/
129 //______________________________________________________________
130 Bool_t AliAnalysisTaskLeadingTrackUE::Notify()
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)
139 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
140 Float_t xsection = 0;
143 TFile *curfile = tree->GetCurrentFile();
145 Error("Notify","No current file");
149 AliAnalysisHelperJetTasks::PythiaInfoFromFile(curfile->GetName(),xsection,trials);
152 //fHistosUE->GetXsec()->Fill("<#sigma>",xsection);
154 // construct average trials
155 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
156 if(trials>=nEntries && nEntries>0.)fAvgTrials = trials/nEntries;
162 //____________________________________________________________________
163 void AliAnalysisTaskLeadingTrackUE::ConnectInputData(Option_t* /*option*/)
166 // Connect the input data
167 if (fDebug > 1) AliInfo("ConnectInputData() ");
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.
173 // Delta AODs are also accepted.
175 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
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");
186 AliFatal("I can't get any AOD Event Handler");
191 // Initialize common pointers
196 //____________________________________________________________________
197 void AliAnalysisTaskLeadingTrackUE::CreateOutputObjects()
199 // Create the output container
201 if (fDebug > 1) AliInfo("CreateOutputObjects()");
203 // Initialize class with main algorithms, event and track selection.
204 fAnalyseUE = new AliAnalyseLeadingTrackUE();
205 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fLeadingTrackEtaCut);
206 fAnalyseUE->SetDebug(fDebug);
208 // Initialize output list of containers
209 if (fListOfHistos != NULL){
210 delete fListOfHistos;
211 fListOfHistos = NULL;
214 fListOfHistos = new TList();
215 fListOfHistos->SetOwner(kTRUE);
218 // Initialize class to handle histograms
219 fHistosUE = new AliUEHistograms("AliUEHistograms", "123");
221 // add histograms to list
222 fListOfHistos->Add(fHistosUE);
224 //fListOfHistos->Add(new TH2F("multVsLeadStep5", ";multiplicity;leading pT", 100, -0.5, 99.5, 20, 0, 10));
226 // Add task configuration to output list
230 PostData(0,fListOfHistos);
233 //____________________________________________________________________
234 void AliAnalysisTaskLeadingTrackUE::Exec(Option_t */*option*/)
236 // array of MC particles
238 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
239 if (!fArrayMC)AliFatal("No array of MC particles found !!!");
242 // Get number of trials from MC header
244 // Float_t nTrials = 1;
246 fMcEvent = fMcHandler->MCEvent();
248 // TO-DO: extend to PHOJET
249 // AliGenPythiaEventHeader* pythiaGenHeader = AliAnalysisHelperJetTasks::GetPythiaEventHeader(fMcEvent);
250 // if(pythiaGenHeader){
251 // nTrials = pythiaGenHeader->Trials();
257 //fHistosUE->GetTrials()->Fill("#sum{ntrials}",fAvgTrials);
260 if (fMode) AnalyseCorrectionMode();
261 else AnalyseDataMode();
263 PostData(0,fListOfHistos);
266 //____________________________________________________________________
267 void AliAnalysisTaskLeadingTrackUE::Terminate(Option_t */*option*/)
270 // Terminate analysis
271 if( fDebug > 1 ) AliInfo("End analysis");
273 if (!gROOT->IsBatch()){
274 fListOfHistos = dynamic_cast<TList*> (GetOutputData(0));
276 AliError("Histogram List is not available");
282 AliInfo(" Batch mode, not histograms will be shown...");
288 /******************** ANALYSIS METHODS *****************************/
290 //____________________________________________________________________
291 void AliAnalysisTaskLeadingTrackUE::AddSettingsTree()
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);
304 //____________________________________________________________________
305 void AliAnalysisTaskLeadingTrackUE::AnalyseCorrectionMode()
307 // Run the analysis on MC to get the correction maps
309 // if fReduceMemoryFootprint step 3,4,5,7,9 are not filled
311 PostData(0,fListOfHistos);
313 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
315 //PROCESS TYPE (ND,SD,DD)
316 AliAnalysisHelperJetTasks::MCProcessType eventId = AliAnalysisHelperJetTasks::kInvalidProcess;
317 AliGenEventHeader* genHeader = fMcEvent->GenEventHeader();
318 eventId = AliAnalysisHelperJetTasks::GetPythiaEventProcessType(genHeader,kFALSE);
320 eventId = AliAnalysisHelperJetTasks::GetDPMjetEventProcessType(genHeader,kFALSE);
321 if (eventId<0 && fDebug>1)
322 AliInfo("No Pythia or Phojet Header retrived!");
325 if (eventId == AliAnalysisHelperJetTasks::kND)fillId = 0;
326 if (eventId == AliAnalysisHelperJetTasks::kSD)fillId = 1;
327 if (eventId == AliAnalysisHelperJetTasks::kDD)fillId = 2;
330 fHistosUE->FillEvent(fillId, -1);
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))
336 // Get MC-true leading particle (but do not cut out events!)
337 TObjArray *ltMC = (TObjArray*)fAnalyseUE->FindLeadingObjects(fArrayMC);
338 AliVParticle* leadingMC = 0;
340 leadingMC = (AliVParticle*) ltMC->At(0);
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
345 // Sort MC-true charged particles
346 // as output you get an array of 3 lists of particles belonging to different regions:
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)
356 fHistosUE->Fill(fillId,0,AliUEHist::kCFStepAll,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
358 // Trigger selection ************************************************
359 if (fAnalyseUE->TriggerSelection(fInputHandler))
362 Bool_t select = kFALSE;
363 if (fSelectBit) select = AliAnalysisHelperJetTasks::TestSelectInfo(fSelectBit);
365 fHistosUE->FillEvent(fillId, -2);
369 // Count events that pass AliPhysicsSelection
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)
374 fHistosUE->Fill(fillId,0,AliUEHist::kCFStepTriggered,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
376 // count number of MC tracks above 150 MeV/c
378 if (leadingMC && leadingMC->Pt() > 0.15)
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)
385 //((TH2F*)fListOfHistos->FindObject("multVsLeadStep5"))->Fill(nMCTracks, leadingMC->Pt());
387 // Vertex selection *************************************************
388 if (fAnalyseUE->VertexSelection(fAOD, fnTracksVertex, fZVertex))
390 // Count events that pass Vertex selection
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)
395 fHistosUE->Fill(fillId,0,AliUEHist::kCFStepVertex,leadingMC,(TList*)regionSortedParticlesMC->At(0),(TList*)regionSortedParticlesMC->At(1),(TList*)regionsMinMaxMC->At(0),(TList*)regionsMinMaxMC->At(1));
397 // fill here for tracking efficiency
398 // loop over particle species
399 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
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);
405 fHistosUE->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, 0, 0, particleSpecies);
407 delete primMCParticles;
408 delete primRecoTracksMatched;
409 delete allRecoTracksMatched;
412 // Get Reconstructed leading particle *******************************
413 TObjArray *ltRECO = fAnalyseUE->FindLeadingObjects(fAOD);
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
420 fHistosUE->Fill(leadingMC, (AliVParticle*)ltRECO->At(0));
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)
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));
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)
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
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));
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)
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
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)
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));
464 // STEP 8 for reduced efficiency study
465 FillReducedEfficiency(fillId, AliUEHist::kCFStepBiasStudy, ltRECO, kFALSE);
466 if (!fReduceMemoryFootprint)
467 FillReducedEfficiency(fillId, AliUEHist::kCFStepBiasStudy2, ltRECO, kTRUE);
469 // count number of reco tracks above 150 MeV/c
470 Int_t nRecoTracks = 0;
471 if (((AliVParticle*) ltRECO->At(0))->Pt() > 0.15)
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)
478 if (leadingMC && leadingMC->Pt() > 0.5)
479 fHistosUE->GetCorrelationMultiplicity()->Fill(nMCTracks, nRecoTracks);
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)
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)
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
497 delete regionSortedParticlesRECOLTMC;
498 delete regionsMinMaxRECOLTMC;
499 delete regionSortedParticlesRECOLTMC2;
500 delete regionsMinMaxRECOLTMC2;
501 delete regionSortedParticlesRECO;
502 delete regionsMinMaxRECO;
507 } //phyiscs selection
511 delete regionSortedParticlesMC;
512 delete regionsMinMaxMC;
515 //____________________________________________________________________
516 void AliAnalysisTaskLeadingTrackUE::AnalyseDataMode()
519 // Run the analysis on DATA or MC to get raw distributions
521 PostData(0,fListOfHistos);
523 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
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);
529 // Trigger selection ************************************************
530 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
532 Bool_t select = kFALSE;
533 if (fSelectBit) select = AliAnalysisHelperJetTasks::TestSelectInfo(fSelectBit);
536 fHistosUE->FillEvent(eventId, -2);
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);
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 (?)
549 // Get Reconstructed leading particle *******************************
550 TObjArray *ltRECO = fAnalyseUE->FindLeadingObjects(fAOD);
553 // fill control distributions
554 fHistosUE->Fill(0, (AliVParticle*)ltRECO->At(0));
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);
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)
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));
568 // STEP 8 and 9 for reduced efficiency study
569 FillReducedEfficiency(eventId, AliUEHist::kCFStepBiasStudy, ltRECO, kFALSE);
570 FillReducedEfficiency(eventId, AliUEHist::kCFStepBiasStudy2, ltRECO, kTRUE);
573 delete regionSortedParticlesRECO;
574 delete regionsMinMaxRECO;
577 //____________________________________________________________________
578 void AliAnalysisTaskLeadingTrackUE::FillReducedEfficiency(Int_t eventId, AliUEHist::CFStep step, const TObjArray* ltRECO, Bool_t twoStep)
580 // remove leading particle using fkTrackingEfficiency and use subleading particle to fill the histograms
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
584 if (!fkTrackingEfficiency)
587 TObjArray* particleList = new TObjArray(*ltRECO);
588 AliVParticle* leading = (AliVParticle*) particleList->At(0);
595 // remove particles depending on tracking efficiency
596 Int_t count = (twoStep) ? 2 : 1;
598 for (Int_t i=0; i<count; i++)
600 Float_t trackingEff = fkTrackingEfficiency->GetBinContent(fkTrackingEfficiency->GetXaxis()->FindBin(leading->Pt()));
602 trackingEff = 0.5 * (trackingEff + 1);
604 if (gRandom->Uniform() > trackingEff)
606 //Printf("LOWEFF: Removing leading particle");
607 particleList->RemoveAt(0);
608 particleList->Compress();
611 if (particleList->GetEntries() == 0)
617 leading = (AliVParticle*) particleList->At(0);
620 TObjArray *regionSortedParticlesRECOLowEff = fAnalyseUE->SortRegions(leading, particleList, 0);
621 TObjArray *regionsMinMaxRECOLowEff = fAnalyseUE->GetMinMaxRegion((TList*)regionSortedParticlesRECOLowEff->At(2), (TList*)regionSortedParticlesRECOLowEff->At(3));
623 fHistosUE->Fill(eventId,0,step,leading,(TList*)regionSortedParticlesRECOLowEff->At(0), (TList*)regionSortedParticlesRECOLowEff->At(1), (TList*)regionsMinMaxRECOLowEff->At(0), (TList*)regionsMinMaxRECOLowEff->At(1));
625 delete regionSortedParticlesRECOLowEff;
626 delete regionsMinMaxRECOLowEff;
630 //____________________________________________________________________
631 void AliAnalysisTaskLeadingTrackUE::Initialize()
634 fInputHandler = (AliInputEventHandler*)
635 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
637 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());