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 fTriggerSelectCharge(0),
118 fTriggerRestrictEta(-1),
119 fEtaOrdering(kFALSE),
120 fCutConversions(kFALSE),
121 fCutResonances(kFALSE),
122 fFillOnlyStep0(kFALSE),
124 fRejectCentralityOutliers(kFALSE),
127 // Default constructor
128 // Define input and output slots here
129 // Input slot #0 works with a TChain
130 DefineInput(0, TChain::Class());
131 // Output slot #0 writes into a TList container
132 DefineOutput(0, TList::Class());
135 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
139 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
140 delete fListOfHistos;
143 //____________________________________________________________________
144 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
147 // Connect the input data
148 if (fDebug > 1) AliInfo("ConnectInputData() ");
150 // Since AODs can either be connected to the InputEventHandler
151 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
152 // we need to get the pointer to the AODEvent correctly.
154 // Delta AODs are also accepted.
156 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
158 if( handler && handler->InheritsFrom("AliAODInputHandler") )
160 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
161 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
165 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
166 if (handler && handler->InheritsFrom("AliAODHandler") )
168 fAOD = ((AliAODHandler*)handler)->GetAOD();
169 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
173 AliWarning("I can't get any AOD Event Handler");
177 if (handler && handler->InheritsFrom("AliESDInputHandler") )
179 // pointer received per event in ::Exec
180 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
183 // Initialize common pointers
187 //____________________________________________________________________
188 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
190 // Create the output container
192 if (fDebug > 1) AliInfo("CreateOutputObjects()");
194 // Initialize class with main algorithms, event and track selection.
195 fAnalyseUE = new AliAnalyseLeadingTrackUE();
196 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
197 fAnalyseUE->SetDebug(fDebug);
198 fAnalyseUE->DefineESDCuts(fFilterBit);
200 // Initialize output list of containers
201 if (fListOfHistos != NULL){
202 delete fListOfHistos;
203 fListOfHistos = NULL;
206 fListOfHistos = new TList();
207 fListOfHistos->SetOwner(kTRUE);
210 // Initialize class to handle histograms
211 const char* histType = "4";
214 fHistos = new AliUEHistograms("AliUEHistogramsSame", histType);
215 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", histType);
217 fHistos->SetSelectCharge(fSelectCharge);
218 fHistosMixed->SetSelectCharge(fSelectCharge);
220 fHistos->SetSelectCharge(fTriggerSelectCharge);
221 fHistosMixed->SetSelectCharge(fTriggerSelectCharge);
223 fHistos->SetTriggerRestrictEta(fTriggerRestrictEta);
224 fHistosMixed->SetTriggerRestrictEta(fTriggerRestrictEta);
226 fHistos->SetOnlyOneEtaSide(fOnlyOneEtaSide);
227 fHistosMixed->SetOnlyOneEtaSide(fOnlyOneEtaSide);
229 fHistos->SetEtaOrdering(fEtaOrdering);
230 fHistosMixed->SetEtaOrdering(fEtaOrdering);
232 fHistos->SetPairCuts(fCutConversions, fCutResonances);
233 fHistosMixed->SetPairCuts(fCutConversions, fCutResonances);
235 // add histograms to list
236 fListOfHistos->Add(fHistos);
237 fListOfHistos->Add(fHistosMixed);
239 fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
240 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));
241 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
242 fListOfHistos->Add(new TH2F("mixedDist", ";centrality;tracks;events", 101, 0, 101, 200, 0, fMixingTracks * 1.5));
243 fListOfHistos->Add(new TH1F("pids", ";pdg;tracks", 2001, -1000.5, 1000.5));
245 PostData(0,fListOfHistos);
247 // Add task configuration to output list
251 Int_t trackDepth = fMixingTracks;
252 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemented of AliEventPoolManager
254 Int_t nCentralityBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
255 Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
257 Int_t nZvtxBins = 7+1+7;
258 // bins for second buffer are shifted by 100 cm
259 Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7, 93, 95, 97, 99, 101, 103, 105, 107 };
260 Double_t* zvtxbin = vertexBins;
262 if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
264 nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
265 zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
268 fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
271 //____________________________________________________________________
272 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
274 // receive ESD pointer if we are not running AOD analysis
277 AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
278 if (handler && handler->InheritsFrom("AliESDInputHandler"))
279 fESD = ((AliESDInputHandler*)handler)->GetEvent();
287 fMcEvent = fMcHandler->MCEvent();
291 // array of MC particles
292 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
294 AliFatal("No array of MC particles found !!!");
297 AnalyseCorrectionMode();
299 else AnalyseDataMode();
302 /******************** ANALYSIS METHODS *****************************/
304 //____________________________________________________________________
305 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
307 //Write settings to output list
308 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
309 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
310 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
311 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
312 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
313 settingsTree->Branch("fOnlyOneEtaSide", &fOnlyOneEtaSide,"OnlyOneEtaSide/I");
314 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
315 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
316 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
317 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
318 settingsTree->Branch("fSelectCharge", &fSelectCharge,"SelectCharge/I");
319 settingsTree->Branch("fTriggerSelectCharge", &fTriggerSelectCharge,"TriggerSelectCharge/I");
320 settingsTree->Branch("fTriggerRestrictEta", &fTriggerRestrictEta,"TriggerRestrictEta/D");
321 settingsTree->Branch("fEtaOrdering", &fEtaOrdering,"EtaOrdering/O");
322 settingsTree->Branch("fCutConversions", &fCutConversions,"CutConversions/O");
323 settingsTree->Branch("fCutResonances", &fCutResonances,"CutResonances/O");
324 settingsTree->Branch("fFillpT", &fFillpT,"FillpT/O");
325 settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
326 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
327 settingsTree->Branch("fSkipTrigger", &fSkipTrigger,"SkipTrigger/O");
328 settingsTree->Branch("fRejectCentralityOutliers", &fRejectCentralityOutliers,"RejectCentralityOutliers/O");
330 settingsTree->Fill();
331 fListOfHistos->Add(settingsTree);
334 //____________________________________________________________________
335 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
337 // Run the analysis on MC to get the correction maps
340 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
342 Double_t centrality = 0;
344 if (fCentralityMethod.Length() > 0)
346 AliCentrality *centralityObj = 0;
348 centralityObj = fAOD->GetHeader()->GetCentralityP();
350 centralityObj = fESD->GetCentrality();
354 centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
355 AliInfo(Form("Centrality is %f", centrality));
359 Printf("WARNING: Centrality object is 0");
364 // Support for ESD and AOD based analysis
365 AliVEvent* inputEvent = fAOD;
373 fHistos->SetRunNumber(inputEvent->GetRunNumber());
374 bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
377 TObject* mc = fArrayMC;
382 fHistos->FillEvent(centrality, -1);
384 // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
385 TObject* vertexSupplier = fMcEvent;
387 vertexSupplier = fAOD->FindListObject(AliAODMCHeader::StdBranchName());
389 if (!fAnalyseUE->VertexSelection(vertexSupplier, 0, fZVertex))
394 zVtx = ((AliAODMCHeader*) vertexSupplier)->GetVtxZ();
396 zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
402 TObjArray* tracksMC = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
406 for (Int_t i=0; i<fArrayMC->GetEntriesFast(); i++)
407 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(((AliAODMCParticle*) fArrayMC->At(i))->PdgCode());
411 for (Int_t i=0; i<fMcEvent->GetNumberOfTracks(); i++)
412 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(fMcEvent->GetTrack(i)->PdgCode());
418 // (MC-true all particles)
420 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, 0, weight);
425 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
427 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
428 for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)
429 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, pool->GetEvent(jMix), 1.0 / pool->GetCurrentNEvents(), (jMix == 0));
430 pool->UpdatePool(CloneAndReduceTrackList(tracksMC));
433 // Printf("trigger: %d", ((AliInputEventHandler*)fInputHandler)->IsEventSelected());
435 // Trigger selection ************************************************
436 if (!fFillOnlyStep0 && (fSkipTrigger || fAnalyseUE->TriggerSelection(fInputHandler)))
438 // (MC-true all particles)
440 if (!fReduceMemoryFootprint)
441 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC, 0, weight);
443 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
446 AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
450 // Vertex selection *************************************************
451 if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
453 // fill here for tracking efficiency
454 // loop over particle species
456 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
458 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
459 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
460 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
462 fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies, centrality);
464 // Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
466 delete primMCParticles;
467 delete primRecoTracksMatched;
468 delete allRecoTracksMatched;
471 // (MC-true all particles)
473 if (!fReduceMemoryFootprint)
474 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC, 0, weight);
476 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
478 // Get MC primaries that match reconstructed track
479 TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
481 // (RECO-matched (quantities from MC particle) primary particles)
483 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim, 0, weight);
485 // Get MC primaries + secondaries that match reconstructed track
486 TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
488 // (RECO-matched (quantities from MC particle) all particles)
490 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll, 0, weight);
493 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
498 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
500 // two track cut, STEP 8
501 if (fTwoTrackEfficiencyCut > 0)
502 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
505 if (fFillMixed && !fSkipStep6)
507 AliEventPool* pool2 = fPoolMgr->GetEventPool(centrality, zVtx + 100);
508 //pool2->PrintInfo();
509 if (pool2->IsReady() || pool2->NTracksInPool() > fMixingTracks / 10 || pool2->GetCurrentNEvents() >= 5)
510 for (Int_t jMix=0; jMix<pool2->GetCurrentNEvents(); jMix++)
513 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0));
515 // two track cut, STEP 8
516 if (fTwoTrackEfficiencyCut > 0)
517 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
519 pool2->UpdatePool(CloneAndReduceTrackList(tracks));
522 if (0 && !fReduceMemoryFootprint)
524 // make list of secondaries (matched with MC)
525 TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
526 for (Int_t i=0; i<tracksRecoMatchedAll->GetEntriesFast(); i++)
527 if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
528 tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
530 // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
531 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll, weight);
533 // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
534 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries, weight);
536 // plot delta phi vs process id of secondaries
537 // trigger particles: primaries in 4 < pT < 10
538 // associated particles: secondaries in 1 < pT < 10
540 for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntriesFast(); i++)
542 AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
544 if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
547 for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntriesFast(); j++)
549 AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
551 if (particle->Pt() < 1 || particle->Pt() > 10)
554 if (particle->Pt() > triggerParticle->Pt())
557 Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
558 if (deltaPhi > 1.5 * TMath::Pi())
559 deltaPhi -= TMath::TwoPi();
560 if (deltaPhi < -0.5 * TMath::Pi())
561 deltaPhi += TMath::TwoPi();
563 Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
565 ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
569 delete tracksRecoMatchedSecondaries;
572 delete tracksRecoMatchedPrim;
573 delete tracksRecoMatchedAll;
581 //____________________________________________________________________
582 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
585 // Run the analysis on DATA or MC to get raw distributions
587 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
589 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
594 // skip not selected events here (the AOD is not updated for those)
595 if (!(fInputHandler->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kUserDefined)))
598 Double_t centrality = 0;
600 AliCentrality *centralityObj = 0;
601 if (fCentralityMethod.Length() > 0)
604 centralityObj = fAOD->GetHeader()->GetCentralityP();
606 centralityObj = fESD->GetCentrality();
609 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
610 //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
619 if (fAOD->GetVZEROData())
622 for (Int_t i=0; i<64; i++)
623 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
627 AliInfo("Rejecting event due to too small V0 multiplicity");
633 AliInfo(Form("Centrality is %f", centrality));
636 // Support for ESD and AOD based analysis
637 AliVEvent* inputEvent = fAOD;
641 Float_t bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
643 fHistos->SetRunNumber(inputEvent->GetRunNumber());
645 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
646 fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
648 // Trigger selection ************************************************
649 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
651 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
652 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
654 // Vertex selection *************************************************
655 if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
657 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
658 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
661 if (centrality < 0 && !fCompareCentralities)
664 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
666 // check for outlier in centrality vs number of tracks (rough constants extracted from correlation histgram)
667 Bool_t reject = kFALSE;
668 if (fRejectCentralityOutliers)
670 if (centrality > 40 && centrality <= 50 && tracks->GetEntriesFast() > 1160)
672 if (centrality > 50 && centrality <= 60 && tracks->GetEntriesFast() > 650)
674 if (centrality > 60 && centrality <= 70 && tracks->GetEntriesFast() > 370)
676 if (centrality > 70 && centrality <= 80 && tracks->GetEntriesFast() > 220)
678 if (centrality > 80 && centrality <= 90 && tracks->GetEntriesFast() > 130)
680 if (centrality > 90 && tracks->GetEntriesFast() > 75)
686 AliInfo(Form("Rejecting event due to centrality vs tracks correlation: %f %d", centrality, tracks->GetEntriesFast()));
687 fHistos->FillEvent(centrality, AliUEHist::kCFStepAnaTopology);
693 // create a list of reduced objects. This speeds up processing and reduces memory consumption for the event pool
694 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
697 //Printf("Accepted %d tracks", tracks->GetEntries());
699 const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
700 Double_t zVtx = vertex->GetZ();
706 // fill two different centralities (syst study)
707 // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
708 if (fCompareCentralities)
711 // Fill containers at STEP 6 (reconstructed)
715 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
717 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
719 if (fTwoTrackEfficiencyCut > 0)
720 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
723 // fill second time with SPD centrality
724 if (fCompareCentralities && centralityObj)
726 centrality = centralityObj->GetCentralityPercentile("CL1");
727 if (centrality >= 0 && !fSkipStep6)
728 fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
735 // 1. First get an event pool corresponding in mult (cent) and
736 // zvertex to the current event. Once initialized, the pool
737 // should contain nMix (reduced) events. This routine does not
738 // pre-scan the chain. The first several events of every chain
739 // will be skipped until the needed pools are filled to the
740 // specified depth. If the pool categories are not too rare, this
741 // should not be a problem. If they are rare, you could lose
744 // 2. Collect the whole pool's content of tracks into one TObjArray
745 // (bgTracks), which is effectively a single background super-event.
747 // 3. The reduced and bgTracks arrays must both be passed into
748 // FillCorrelations(). Also nMix should be passed in, so a weight
749 // of 1./nMix can be applied.
751 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
754 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
758 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
761 Int_t nMix = pool->GetCurrentNEvents();
762 // cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
764 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
765 ((TH2F*) fListOfHistos->FindObject("mixedDist"))->Fill(centrality, pool->NTracksInPool());
767 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(3);
769 // Fill mixed-event histos here
770 for (Int_t jMix=0; jMix<nMix; jMix++)
772 TObjArray* bgTracks = pool->GetEvent(jMix);
775 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0));
777 if (fTwoTrackEfficiencyCut > 0)
778 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
782 // ownership is with the pool now
783 pool->UpdatePool(tracksClone);
790 TObjArray* AliAnalysisTaskPhiCorrelations::CloneAndReduceTrackList(TObjArray* tracks)
792 // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
794 TObjArray* tracksClone = new TObjArray;
795 tracksClone->SetOwner(kTRUE);
797 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
799 AliVParticle* particle = (AliVParticle*) tracks->At(i);
800 tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
806 //____________________________________________________________________
807 void AliAnalysisTaskPhiCorrelations::Initialize()
810 fInputHandler = (AliInputEventHandler*)
811 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
813 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());