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),
87 fTwoTrackEfficiencyCut(0),
90 // pointers to UE classes
94 fkTrackingEfficiency(0x0),
95 // handlers and events
103 // histogram settings
106 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
108 fCentralityMethod("V0M"),
115 fUseChargeHadrons(kFALSE),
117 fTriggerRestrictEta(-1),
118 fEtaOrdering(kFALSE),
119 fCutConversions(kFALSE),
120 fCutResonances(kFALSE),
121 fFillOnlyStep0(kFALSE),
125 // Default constructor
126 // Define input and output slots here
127 // Input slot #0 works with a TChain
128 DefineInput(0, TChain::Class());
129 // Output slot #0 writes into a TList container
130 DefineOutput(0, TList::Class());
133 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
137 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
138 delete fListOfHistos;
141 //____________________________________________________________________
142 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
145 // Connect the input data
146 if (fDebug > 1) AliInfo("ConnectInputData() ");
148 // Since AODs can either be connected to the InputEventHandler
149 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
150 // we need to get the pointer to the AODEvent correctly.
152 // Delta AODs are also accepted.
154 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
156 if( handler && handler->InheritsFrom("AliAODInputHandler") )
158 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
159 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
163 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
164 if (handler && handler->InheritsFrom("AliAODHandler") )
166 fAOD = ((AliAODHandler*)handler)->GetAOD();
167 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
171 AliWarning("I can't get any AOD Event Handler");
175 if (handler && handler->InheritsFrom("AliESDInputHandler") )
177 // pointer received per event in ::Exec
178 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
181 // Initialize common pointers
185 //____________________________________________________________________
186 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
188 // Create the output container
190 if (fDebug > 1) AliInfo("CreateOutputObjects()");
192 // Initialize class with main algorithms, event and track selection.
193 fAnalyseUE = new AliAnalyseLeadingTrackUE();
194 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
195 fAnalyseUE->SetDebug(fDebug);
196 fAnalyseUE->DefineESDCuts(fFilterBit);
198 // Initialize output list of containers
199 if (fListOfHistos != NULL){
200 delete fListOfHistos;
201 fListOfHistos = NULL;
204 fListOfHistos = new TList();
205 fListOfHistos->SetOwner(kTRUE);
208 // Initialize class to handle histograms
209 const char* histType = "4";
212 fHistos = new AliUEHistograms("AliUEHistogramsSame", histType);
213 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", histType);
215 fHistos->SetSelectCharge(fSelectCharge);
216 fHistosMixed->SetSelectCharge(fSelectCharge);
218 fHistos->SetTriggerRestrictEta(fTriggerRestrictEta);
219 fHistosMixed->SetTriggerRestrictEta(fTriggerRestrictEta);
221 fHistos->SetOnlyOneEtaSide(fOnlyOneEtaSide);
222 fHistosMixed->SetOnlyOneEtaSide(fOnlyOneEtaSide);
224 fHistos->SetEtaOrdering(fEtaOrdering);
225 fHistosMixed->SetEtaOrdering(fEtaOrdering);
227 fHistos->SetPairCuts(fCutConversions, fCutResonances);
228 fHistosMixed->SetPairCuts(fCutConversions, fCutResonances);
230 // add histograms to list
231 fListOfHistos->Add(fHistos);
232 fListOfHistos->Add(fHistosMixed);
234 fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
235 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));
236 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
237 fListOfHistos->Add(new TH2F("mixedDist", ";centrality;tracks;events", 101, 0, 101, 200, 0, fMixingTracks * 1.5));
238 fListOfHistos->Add(new TH1F("pids", ";pdg;tracks", 2001, -1000.5, 1000.5));
240 PostData(0,fListOfHistos);
242 // Add task configuration to output list
246 Int_t trackDepth = fMixingTracks;
247 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemented of AliEventPoolManager
249 Int_t nCentralityBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
250 Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
252 Int_t nZvtxBins = 7+1+7;
253 // bins for second buffer are shifted by 100 cm
254 Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7, 93, 95, 97, 99, 101, 103, 105, 107 };
255 Double_t* zvtxbin = vertexBins;
257 if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
259 nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
260 zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
263 fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
266 //____________________________________________________________________
267 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
269 // receive ESD pointer if we are not running AOD analysis
272 AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
273 if (handler && handler->InheritsFrom("AliESDInputHandler"))
274 fESD = ((AliESDInputHandler*)handler)->GetEvent();
282 fMcEvent = fMcHandler->MCEvent();
286 // array of MC particles
287 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
289 AliFatal("No array of MC particles found !!!");
292 AnalyseCorrectionMode();
294 else AnalyseDataMode();
297 /******************** ANALYSIS METHODS *****************************/
299 //____________________________________________________________________
300 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
302 //Write settings to output list
303 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
304 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
305 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
306 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
307 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
308 settingsTree->Branch("fOnlyOneEtaSide", &fOnlyOneEtaSide,"OnlyOneEtaSide/I");
309 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
310 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
311 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
312 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
313 settingsTree->Branch("fSelectCharge", &fSelectCharge,"SelectCharge/I");
314 settingsTree->Branch("fTriggerRestrictEta", &fTriggerRestrictEta,"TriggerRestrictEta/D");
315 settingsTree->Branch("fEtaOrdering", &fEtaOrdering,"EtaOrdering/O");
316 settingsTree->Branch("fCutConversions", &fCutConversions,"CutConversions/O");
317 settingsTree->Branch("fCutResonances", &fCutResonances,"CutResonances/O");
318 settingsTree->Branch("fFillpT", &fFillpT,"FillpT/O");
319 settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
320 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
321 settingsTree->Branch("fSkipTrigger", &fSkipTrigger,"SkipTrigger/O");
323 settingsTree->Fill();
324 fListOfHistos->Add(settingsTree);
327 //____________________________________________________________________
328 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
330 // Run the analysis on MC to get the correction maps
333 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
335 Double_t centrality = 0;
337 if (fCentralityMethod.Length() > 0)
339 AliCentrality *centralityObj = 0;
341 centralityObj = fAOD->GetHeader()->GetCentralityP();
343 centralityObj = fESD->GetCentrality();
347 centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
348 AliInfo(Form("Centrality is %f", centrality));
352 Printf("WARNING: Centrality object is 0");
357 // Support for ESD and AOD based analysis
358 AliVEvent* inputEvent = fAOD;
363 fHistos->SetRunNumber(inputEvent->GetRunNumber());
365 TObject* mc = fArrayMC;
370 fHistos->FillEvent(centrality, -1);
372 // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
373 TObject* vertexSupplier = fMcEvent;
375 vertexSupplier = fAOD->FindListObject(AliAODMCHeader::StdBranchName());
377 if (!fAnalyseUE->VertexSelection(vertexSupplier, 0, fZVertex))
382 zVtx = ((AliAODMCHeader*) vertexSupplier)->GetVtxZ();
384 zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
390 TObjArray* tracksMC = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
394 for (Int_t i=0; i<fArrayMC->GetEntriesFast(); i++)
395 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(((AliAODMCParticle*) fArrayMC->At(i))->PdgCode());
399 for (Int_t i=0; i<fMcEvent->GetNumberOfTracks(); i++)
400 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(fMcEvent->GetTrack(i)->PdgCode());
403 // (MC-true all particles)
405 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, 0, weight);
410 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
412 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
413 for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)
414 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, pool->GetEvent(jMix), 1.0 / pool->GetCurrentNEvents(), (jMix == 0));
415 pool->UpdatePool(CloneAndReduceTrackList(tracksMC));
418 // Printf("trigger: %d", ((AliInputEventHandler*)fInputHandler)->IsEventSelected());
420 // Trigger selection ************************************************
421 if (!fFillOnlyStep0 && (fSkipTrigger || fAnalyseUE->TriggerSelection(fInputHandler)))
423 // (MC-true all particles)
425 if (!fReduceMemoryFootprint)
426 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC, 0, weight);
428 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
431 AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
435 // Vertex selection *************************************************
436 if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
438 // fill here for tracking efficiency
439 // loop over particle species
441 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
443 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
444 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
445 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
447 fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies, centrality);
449 // Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
451 delete primMCParticles;
452 delete primRecoTracksMatched;
453 delete allRecoTracksMatched;
456 // (MC-true all particles)
458 if (!fReduceMemoryFootprint)
459 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC, 0, weight);
461 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
463 // Get MC primaries that match reconstructed track
464 TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
466 // (RECO-matched (quantities from MC particle) primary particles)
468 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim, 0, weight);
470 // Get MC primaries + secondaries that match reconstructed track
471 TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
473 // (RECO-matched (quantities from MC particle) all particles)
475 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll, 0, weight);
478 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
483 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
486 if (fFillMixed && !fSkipStep6)
488 AliEventPool* pool2 = fPoolMgr->GetEventPool(centrality, zVtx + 100);
489 //pool2->PrintInfo();
490 if (pool2->IsReady() || pool2->NTracksInPool() > fMixingTracks / 10 || pool2->GetCurrentNEvents() >= 5)
491 for (Int_t jMix=0; jMix<pool2->GetCurrentNEvents(); jMix++)
492 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0));
493 pool2->UpdatePool(CloneAndReduceTrackList(tracks));
496 if (0 && !fReduceMemoryFootprint)
498 // make list of secondaries (matched with MC)
499 TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
500 for (Int_t i=0; i<tracksRecoMatchedAll->GetEntriesFast(); i++)
501 if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
502 tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
504 // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
505 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll, weight);
507 // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
508 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries, weight);
510 // plot delta phi vs process id of secondaries
511 // trigger particles: primaries in 4 < pT < 10
512 // associated particles: secondaries in 1 < pT < 10
514 for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntriesFast(); i++)
516 AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
518 if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
521 for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntriesFast(); j++)
523 AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
525 if (particle->Pt() < 1 || particle->Pt() > 10)
528 if (particle->Pt() > triggerParticle->Pt())
531 Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
532 if (deltaPhi > 1.5 * TMath::Pi())
533 deltaPhi -= TMath::TwoPi();
534 if (deltaPhi < -0.5 * TMath::Pi())
535 deltaPhi += TMath::TwoPi();
537 Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
539 ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
543 delete tracksRecoMatchedSecondaries;
546 delete tracksRecoMatchedPrim;
547 delete tracksRecoMatchedAll;
555 //____________________________________________________________________
556 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
559 // Run the analysis on DATA or MC to get raw distributions
561 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
563 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
565 // skip not selected events here (the AOD is not updated for those)
569 if (!(fInputHandler->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kUserDefined)))
572 Double_t centrality = 0;
574 AliCentrality *centralityObj = 0;
575 if (fCentralityMethod.Length() > 0)
578 centralityObj = fAOD->GetHeader()->GetCentralityP();
580 centralityObj = fESD->GetCentrality();
583 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
584 //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
593 if (fAOD->GetVZEROData())
596 for (Int_t i=0; i<64; i++)
597 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
601 AliInfo("Rejecting event due to too small V0 multiplicity");
607 AliInfo(Form("Centrality is %f", centrality));
610 // Support for ESD and AOD based analysis
611 AliVEvent* inputEvent = fAOD;
615 Float_t bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
617 fHistos->SetRunNumber(inputEvent->GetRunNumber());
619 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
620 fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
622 // Trigger selection ************************************************
623 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
625 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
626 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
628 // Vertex selection *************************************************
629 if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
631 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
632 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
635 if (centrality < 0 && !fCompareCentralities)
638 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
639 // create a list of reduced objects. This speeds up processing and reduces memory consumption for the event pool
640 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
643 //Printf("Accepted %d tracks", tracks->GetEntries());
645 const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
646 Double_t zVtx = vertex->GetZ();
652 // fill two different centralities (syst study)
653 // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
654 if (fCompareCentralities)
657 // Fill containers at STEP 6 (reconstructed)
661 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
663 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
665 if (fTwoTrackEfficiencyCut > 0)
666 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
669 // fill second time with SPD centrality
670 if (fCompareCentralities && centralityObj)
672 centrality = centralityObj->GetCentralityPercentile("CL1");
673 if (centrality >= 0 && !fSkipStep6)
674 fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
681 // 1. First get an event pool corresponding in mult (cent) and
682 // zvertex to the current event. Once initialized, the pool
683 // should contain nMix (reduced) events. This routine does not
684 // pre-scan the chain. The first several events of every chain
685 // will be skipped until the needed pools are filled to the
686 // specified depth. If the pool categories are not too rare, this
687 // should not be a problem. If they are rare, you could lose
690 // 2. Collect the whole pool's content of tracks into one TObjArray
691 // (bgTracks), which is effectively a single background super-event.
693 // 3. The reduced and bgTracks arrays must both be passed into
694 // FillCorrelations(). Also nMix should be passed in, so a weight
695 // of 1./nMix can be applied.
697 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
700 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
704 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
707 Int_t nMix = pool->GetCurrentNEvents();
708 // cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
710 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
711 ((TH2F*) fListOfHistos->FindObject("mixedDist"))->Fill(centrality, pool->NTracksInPool());
713 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(3);
715 // Fill mixed-event histos here
716 for (Int_t jMix=0; jMix<nMix; jMix++)
718 TObjArray* bgTracks = pool->GetEvent(jMix);
721 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0));
723 if (fTwoTrackEfficiencyCut > 0)
724 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
728 // ownership is with the pool now
729 pool->UpdatePool(tracksClone);
736 TObjArray* AliAnalysisTaskPhiCorrelations::CloneAndReduceTrackList(TObjArray* tracks)
738 // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
740 TObjArray* tracksClone = new TObjArray;
741 tracksClone->SetOwner(kTRUE);
743 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
745 AliVParticle* particle = (AliVParticle*) tracks->At(i);
746 tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
752 //____________________________________________________________________
753 void AliAnalysisTaskPhiCorrelations::Initialize()
756 fInputHandler = (AliInputEventHandler*)
757 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
759 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());