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 **************************************************************************/
27 #include "AliAnalysisTaskPhiCorrelations.h"
28 #include "AliAnalyseLeadingTrackUE.h"
29 #include "AliUEHistograms.h"
30 #include "AliUEHist.h"
32 #include "AliAnalysisManager.h"
33 #include "AliAODHandler.h"
34 #include "AliAODInputHandler.h"
35 #include "AliAODMCParticle.h"
36 #include "AliInputEventHandler.h"
38 #include "AliMCEventHandler.h"
39 #include "AliVParticle.h"
40 #include "AliCFContainer.h"
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44 #include "AliMultiplicity.h"
45 #include "AliCentrality.h"
47 #include "AliAODMCHeader.h"
49 #include "AliEventPoolManager.h"
52 ////////////////////////////////////////////////////////////////////////
54 // Analysis class for azimuthal correlation studies
55 // Based on the UE task from Sara Vallero and Jan Fiete
57 // This class needs input AODs.
58 // The output is a list of analysis-specific containers.
60 // The AOD can be either connected to the InputEventHandler
61 // for a chain of AOD files
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
68 // Jan Fiete Grosse-Oetringhaus
70 ////////////////////////////////////////////////////////////////////////
73 ClassImp( AliAnalysisTaskPhiCorrelations )
74 ClassImp( AliDPhiBasicParticle )
76 //____________________________________________________________________
77 AliAnalysisTaskPhiCorrelations:: AliAnalysisTaskPhiCorrelations(const char* name):
78 AliAnalysisTask(name,""),
79 // general configuration
82 fReduceMemoryFootprint(kFALSE),
85 fCompareCentralities(kFALSE),
86 fTwoTrackEfficiencyStudy(kFALSE),
88 // pointers to UE classes
92 fkTrackingEfficiency(0x0),
93 // handlers and events
101 // histogram settings
104 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
106 fCentralityMethod("V0M"),
112 fUseChargeHadrons(kFALSE),
114 fTriggerRestrictEta(-1),
115 fEtaOrdering(kFALSE),
116 fCutConversions(kFALSE),
117 fCutResonances(kFALSE),
120 // Default constructor
121 // Define input and output slots here
122 // Input slot #0 works with a TChain
123 DefineInput(0, TChain::Class());
124 // Output slot #0 writes into a TList container
125 DefineOutput(0, TList::Class());
128 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
132 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
133 delete fListOfHistos;
136 //____________________________________________________________________
137 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
140 // Connect the input data
141 if (fDebug > 1) AliInfo("ConnectInputData() ");
143 // Since AODs can either be connected to the InputEventHandler
144 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
145 // we need to get the pointer to the AODEvent correctly.
147 // Delta AODs are also accepted.
149 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
151 if( handler && handler->InheritsFrom("AliAODInputHandler") )
153 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
154 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
158 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
159 if (handler && handler->InheritsFrom("AliAODHandler") )
161 fAOD = ((AliAODHandler*)handler)->GetAOD();
162 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
166 AliWarning("I can't get any AOD Event Handler");
170 if (handler && handler->InheritsFrom("AliESDInputHandler") )
172 // pointer received per event in ::Exec
173 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
176 // Initialize common pointers
180 //____________________________________________________________________
181 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
183 // Create the output container
185 if (fDebug > 1) AliInfo("CreateOutputObjects()");
187 // Initialize class with main algorithms, event and track selection.
188 fAnalyseUE = new AliAnalyseLeadingTrackUE();
189 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
190 fAnalyseUE->SetDebug(fDebug);
191 fAnalyseUE->DefineESDCuts(fFilterBit);
193 // Initialize output list of containers
194 if (fListOfHistos != NULL){
195 delete fListOfHistos;
196 fListOfHistos = NULL;
199 fListOfHistos = new TList();
200 fListOfHistos->SetOwner(kTRUE);
203 // Initialize class to handle histograms
204 const char* histType = "4";
207 fHistos = new AliUEHistograms("AliUEHistogramsSame", histType);
208 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", histType);
210 fHistos->SetSelectCharge(fSelectCharge);
211 fHistosMixed->SetSelectCharge(fSelectCharge);
213 fHistos->SetTriggerRestrictEta(fTriggerRestrictEta);
214 fHistosMixed->SetTriggerRestrictEta(fTriggerRestrictEta);
216 fHistos->SetEtaOrdering(fEtaOrdering);
217 fHistosMixed->SetEtaOrdering(fEtaOrdering);
219 fHistos->SetPairCuts(fCutConversions, fCutResonances);
220 fHistosMixed->SetPairCuts(fCutConversions, fCutResonances);
222 // add histograms to list
223 fListOfHistos->Add(fHistos);
224 fListOfHistos->Add(fHistosMixed);
226 fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
227 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));
228 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
229 fListOfHistos->Add(new TH2F("mixedDist", ";centrality;tracks;events", 101, 0, 101, 200, 0, fMixingTracks * 1.5));
230 fListOfHistos->Add(new TH1F("pids", ";pdg;tracks", 2001, -1000.5, 1000.5));
232 PostData(0,fListOfHistos);
234 // Add task configuration to output list
238 Int_t trackDepth = fMixingTracks;
239 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemented of AliEventPool
241 Int_t nCentralityBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
242 Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
244 Int_t nZvtxBins = 7+1+7;
245 // bins for second buffer are shifted by 100 cm
246 Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7, 93, 95, 97, 99, 101, 103, 105, 107 };
247 Double_t* zvtxbin = vertexBins;
249 if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
251 nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
252 zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
255 fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
258 //____________________________________________________________________
259 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
261 // receive ESD pointer if we are not running AOD analysis
264 AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
265 if (handler && handler->InheritsFrom("AliESDInputHandler"))
266 fESD = ((AliESDInputHandler*)handler)->GetEvent();
274 fMcEvent = fMcHandler->MCEvent();
278 // array of MC particles
279 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
281 AliFatal("No array of MC particles found !!!");
284 AnalyseCorrectionMode();
286 else AnalyseDataMode();
289 /******************** ANALYSIS METHODS *****************************/
291 //____________________________________________________________________
292 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
294 //Write settings to output list
295 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
296 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
297 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
298 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
299 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
300 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
301 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
302 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
303 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
304 settingsTree->Branch("fSelectCharge", &fSelectCharge,"SelectCharge/I");
305 settingsTree->Branch("fTriggerRestrictEta", &fTriggerRestrictEta,"TriggerRestrictEta/D");
306 settingsTree->Branch("fEtaOrdering", &fEtaOrdering,"EtaOrdering/O");
307 settingsTree->Branch("fCutConversions", &fCutConversions,"CutConversions/O");
308 settingsTree->Branch("fCutResonances", &fCutResonances,"CutResonances/O");
309 settingsTree->Branch("fFillpT", &fFillpT,"FillpT/O");
310 settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
311 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
312 settingsTree->Fill();
313 fListOfHistos->Add(settingsTree);
316 //____________________________________________________________________
317 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
319 // Run the analysis on MC to get the correction maps
322 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
324 Double_t centrality = 0;
326 if (fCentralityMethod.Length() > 0)
328 AliCentrality *centralityObj = 0;
330 centralityObj = fAOD->GetHeader()->GetCentralityP();
332 centralityObj = fESD->GetCentrality();
336 centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
337 AliInfo(Form("Centrality is %f", centrality));
341 Printf("WARNING: Centrality object is 0");
346 // Support for ESD and AOD based analysis
347 AliVEvent* inputEvent = fAOD;
351 TObject* mc = fArrayMC;
356 fHistos->FillEvent(centrality, -1);
358 // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
359 TObject* vertexSupplier = fMcEvent;
361 vertexSupplier = fAOD->FindListObject(AliAODMCHeader::StdBranchName());
363 if (!fAnalyseUE->VertexSelection(vertexSupplier, 0, fZVertex))
368 zVtx = ((AliAODMCHeader*) vertexSupplier)->GetVtxZ();
370 zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
376 TObjArray* tracksMC = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
380 for (Int_t i=0; i<fArrayMC->GetEntriesFast(); i++)
381 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(((AliAODMCParticle*) fArrayMC->At(i))->PdgCode());
385 for (Int_t i=0; i<fMcEvent->GetNumberOfTracks(); i++)
386 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(fMcEvent->GetTrack(i)->PdgCode());
389 // (MC-true all particles)
391 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, 0, weight);
396 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
398 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
399 for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)
400 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, pool->GetEvent(jMix), 1.0 / pool->GetCurrentNEvents(), (jMix == 0));
401 pool->UpdatePool(CloneAndReduceTrackList(tracksMC));
404 // Trigger selection ************************************************
405 if (fAnalyseUE->TriggerSelection(fInputHandler))
407 // (MC-true all particles)
409 if (!fReduceMemoryFootprint)
410 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC, 0, weight);
412 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
415 AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
419 // Vertex selection *************************************************
420 if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
422 // fill here for tracking efficiency
423 // loop over particle species
425 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
427 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
428 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
429 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
431 fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies, centrality);
433 // Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
435 delete primMCParticles;
436 delete primRecoTracksMatched;
437 delete allRecoTracksMatched;
440 // (MC-true all particles)
442 if (!fReduceMemoryFootprint)
443 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC, 0, weight);
445 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
447 // Get MC primaries that match reconstructed track
448 TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
450 // (RECO-matched (quantities from MC particle) primary particles)
452 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim, 0, weight);
454 // Get MC primaries + secondaries that match reconstructed track
455 TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
457 // (RECO-matched (quantities from MC particle) all particles)
459 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll, 0, weight);
462 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
466 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
471 AliEventPool* pool2 = fPoolMgr->GetEventPool(centrality, zVtx + 100);
472 //pool2->PrintInfo();
473 if (pool2->IsReady() || pool2->NTracksInPool() > fMixingTracks / 10 || pool2->GetCurrentNEvents() >= 5)
474 for (Int_t jMix=0; jMix<pool2->GetCurrentNEvents(); jMix++)
475 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0));
476 pool2->UpdatePool(CloneAndReduceTrackList(tracks));
479 if (0 && !fReduceMemoryFootprint)
481 // make list of secondaries (matched with MC)
482 TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
483 for (Int_t i=0; i<tracksRecoMatchedAll->GetEntriesFast(); i++)
484 if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
485 tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
487 // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
488 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll, weight);
490 // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
491 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries, weight);
493 // plot delta phi vs process id of secondaries
494 // trigger particles: primaries in 4 < pT < 10
495 // associated particles: secondaries in 1 < pT < 10
497 for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntriesFast(); i++)
499 AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
501 if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
504 for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntriesFast(); j++)
506 AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
508 if (particle->Pt() < 1 || particle->Pt() > 10)
511 if (particle->Pt() > triggerParticle->Pt())
514 Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
515 if (deltaPhi > 1.5 * TMath::Pi())
516 deltaPhi -= TMath::TwoPi();
517 if (deltaPhi < -0.5 * TMath::Pi())
518 deltaPhi += TMath::TwoPi();
520 Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
522 ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
526 delete tracksRecoMatchedSecondaries;
529 delete tracksRecoMatchedPrim;
530 delete tracksRecoMatchedAll;
538 //____________________________________________________________________
539 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
542 // Run the analysis on DATA or MC to get raw distributions
544 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
546 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
548 // skip not selected events here (the AOD is not updated for those)
552 if (!(fInputHandler->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kUserDefined)))
555 Double_t centrality = 0;
557 AliCentrality *centralityObj = 0;
558 if (fCentralityMethod.Length() > 0)
561 centralityObj = fAOD->GetHeader()->GetCentralityP();
563 centralityObj = fESD->GetCentrality();
566 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
567 //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
576 if (fAOD->GetVZEROData())
579 for (Int_t i=0; i<64; i++)
580 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
584 AliInfo("Rejecting event due to too small V0 multiplicity");
590 AliInfo(Form("Centrality is %f", centrality));
593 // Support for ESD and AOD based analysis
594 AliVEvent* inputEvent = fAOD;
598 fHistos->SetRunNumber(inputEvent->GetRunNumber());
600 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
601 fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
603 // Trigger selection ************************************************
604 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
606 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
607 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
609 // Vertex selection *************************************************
610 if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
612 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
613 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
616 if (centrality < 0 && !fCompareCentralities)
619 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
620 //Printf("Accepted %d tracks", tracks->GetEntries());
622 const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
623 Double_t zVtx = vertex->GetZ();
629 // fill two different centralities (syst study)
630 // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
631 if (fCompareCentralities)
634 // Fill containers at STEP 6 (reconstructed)
637 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
638 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
641 // Two-track effect study
642 if (fTwoTrackEfficiencyStudy)
644 Float_t bSign = (fESD->GetMagneticField() > 0) ? 1 : -1;
645 fHistos->TwoTrackEfficiency(tracks, 0, bSign);
648 // fill second time with SPD centrality
649 if (fCompareCentralities && centralityObj)
651 centrality = centralityObj->GetCentralityPercentile("CL1");
653 fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
660 // 1. First get an event pool corresponding in mult (cent) and
661 // zvertex to the current event. Once initialized, the pool
662 // should contain nMix (reduced) events. This routine does not
663 // pre-scan the chain. The first several events of every chain
664 // will be skipped until the needed pools are filled to the
665 // specified depth. If the pool categories are not too rare, this
666 // should not be a problem. If they are rare, you could lose
669 // 2. Collect the whole pool's content of tracks into one TObjArray
670 // (bgTracks), which is effectively a single background super-event.
672 // 3. The reduced and bgTracks arrays must both be passed into
673 // FillCorrelations(). Also nMix should be passed in, so a weight
674 // of 1./nMix can be applied.
676 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
679 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
683 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
686 Int_t nMix = pool->GetCurrentNEvents();
687 // cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
689 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
690 ((TH2F*) fListOfHistos->FindObject("mixedDist"))->Fill(centrality, pool->NTracksInPool());
692 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(3);
694 // Fill mixed-event histos here
695 for (Int_t jMix=0; jMix<nMix; jMix++)
697 TObjArray* bgTracks = pool->GetEvent(jMix);
699 if (fTwoTrackEfficiencyStudy)
701 Float_t bSign = (fESD->GetMagneticField() > 0) ? 1 : -1;
702 fHistos->TwoTrackEfficiency(tracks, bgTracks, bSign);
705 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, bgTracks, 1.0 / nMix, (jMix == 0));
709 // create a list of reduced objects (to reduce memory consumption) and give ownership to event pool
710 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
711 pool->UpdatePool(tracksClone);
718 TObjArray* AliAnalysisTaskPhiCorrelations::CloneAndReduceTrackList(TObjArray* tracks)
720 // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
722 TObjArray* tracksClone = new TObjArray;
723 tracksClone->SetOwner(kTRUE);
725 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
727 AliVParticle* particle = (AliVParticle*) tracks->At(i);
728 tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
734 //____________________________________________________________________
735 void AliAnalysisTaskPhiCorrelations::Initialize()
738 fInputHandler = (AliInputEventHandler*)
739 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
741 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());