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"
48 #include "AliGenCocktailEventHeader.h"
49 #include "AliGenEventHeader.h"
51 #include "AliEventPoolManager.h"
54 ////////////////////////////////////////////////////////////////////////
56 // Analysis class for azimuthal correlation studies
57 // Based on the UE task from Sara Vallero and Jan Fiete
59 // This class needs input AODs.
60 // The output is a list of analysis-specific containers.
62 // The AOD can be either connected to the InputEventHandler
63 // for a chain of AOD files
65 // to the OutputEventHandler
66 // for a chain of ESD files,
67 // in this case the class should be in the train after the jet-finder
70 // Jan Fiete Grosse-Oetringhaus
72 ////////////////////////////////////////////////////////////////////////
75 ClassImp( AliAnalysisTaskPhiCorrelations )
76 ClassImp( AliDPhiBasicParticle )
78 //____________________________________________________________________
79 AliAnalysisTaskPhiCorrelations:: AliAnalysisTaskPhiCorrelations(const char* name):
80 AliAnalysisTask(name,""),
81 // general configuration
84 fReduceMemoryFootprint(kFALSE),
87 fCompareCentralities(kFALSE),
88 fTwoTrackEfficiencyStudy(kFALSE),
89 fTwoTrackEfficiencyCut(0),
92 fInjectedSignals(kFALSE),
93 // pointers to UE classes
97 fkTrackingEfficiency(0x0),
98 // handlers and events
106 // histogram settings
109 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
111 fCentralityMethod("V0M"),
117 fSelectBit(AliVEvent::kMB|AliVEvent::kUserDefined),
118 fUseChargeHadrons(kFALSE),
120 fTriggerSelectCharge(0),
121 fTriggerRestrictEta(-1),
122 fEtaOrdering(kFALSE),
123 fCutConversions(kFALSE),
124 fCutResonances(kFALSE),
125 fFillOnlyStep0(kFALSE),
127 fRejectCentralityOutliers(kFALSE),
130 // Default constructor
131 // Define input and output slots here
132 // Input slot #0 works with a TChain
133 DefineInput(0, TChain::Class());
134 // Output slot #0 writes into a TList container
135 DefineOutput(0, TList::Class());
138 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
142 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
143 delete fListOfHistos;
146 //____________________________________________________________________
147 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
150 // Connect the input data
151 if (fDebug > 1) AliInfo("ConnectInputData() ");
153 // Since AODs can either be connected to the InputEventHandler
154 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
155 // we need to get the pointer to the AODEvent correctly.
157 // Delta AODs are also accepted.
159 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
161 if( handler && handler->InheritsFrom("AliAODInputHandler") )
163 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
164 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
168 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
169 if (handler && handler->InheritsFrom("AliAODHandler") )
171 fAOD = ((AliAODHandler*)handler)->GetAOD();
172 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
176 AliWarning("I can't get any AOD Event Handler");
180 if (handler && handler->InheritsFrom("AliESDInputHandler") )
182 // pointer received per event in ::Exec
183 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
186 // Initialize common pointers
190 //____________________________________________________________________
191 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
193 // Create the output container
195 if (fDebug > 1) AliInfo("CreateOutputObjects()");
197 // Initialize class with main algorithms, event and track selection.
198 fAnalyseUE = new AliAnalyseLeadingTrackUE();
199 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
200 fAnalyseUE->SetDebug(fDebug);
201 fAnalyseUE->DefineESDCuts(fFilterBit);
202 fAnalyseUE->SetEventSelection(fSelectBit);
204 // Initialize output list of containers
205 if (fListOfHistos != NULL){
206 delete fListOfHistos;
207 fListOfHistos = NULL;
210 fListOfHistos = new TList();
211 fListOfHistos->SetOwner(kTRUE);
214 // Initialize class to handle histograms
215 const char* histType = "4";
218 fHistos = new AliUEHistograms("AliUEHistogramsSame", histType);
219 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", histType);
221 fHistos->SetSelectCharge(fSelectCharge);
222 fHistosMixed->SetSelectCharge(fSelectCharge);
224 fHistos->SetSelectCharge(fTriggerSelectCharge);
225 fHistosMixed->SetSelectCharge(fTriggerSelectCharge);
227 fHistos->SetTriggerRestrictEta(fTriggerRestrictEta);
228 fHistosMixed->SetTriggerRestrictEta(fTriggerRestrictEta);
230 fHistos->SetOnlyOneEtaSide(fOnlyOneEtaSide);
231 fHistosMixed->SetOnlyOneEtaSide(fOnlyOneEtaSide);
233 fHistos->SetEtaOrdering(fEtaOrdering);
234 fHistosMixed->SetEtaOrdering(fEtaOrdering);
236 fHistos->SetPairCuts(fCutConversions, fCutResonances);
237 fHistosMixed->SetPairCuts(fCutConversions, fCutResonances);
239 // add histograms to list
240 fListOfHistos->Add(fHistos);
241 fListOfHistos->Add(fHistosMixed);
243 fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
244 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));
245 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
246 fListOfHistos->Add(new TH2F("mixedDist", ";centrality;tracks;events", 101, 0, 101, 200, 0, fMixingTracks * 1.5));
247 fListOfHistos->Add(new TH1F("pids", ";pdg;tracks", 2001, -1000.5, 1000.5));
249 PostData(0,fListOfHistos);
251 // Add task configuration to output list
255 Int_t trackDepth = fMixingTracks;
256 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemented of AliEventPoolManager
258 Int_t nCentralityBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
259 Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
261 Int_t nZvtxBins = 7+1+7;
262 // bins for second buffer are shifted by 100 cm
263 Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7, 93, 95, 97, 99, 101, 103, 105, 107 };
264 Double_t* zvtxbin = vertexBins;
266 if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
268 nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
269 zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
272 fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
275 //____________________________________________________________________
276 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
278 // receive ESD pointer if we are not running AOD analysis
281 AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
282 if (handler && handler->InheritsFrom("AliESDInputHandler"))
283 fESD = ((AliESDInputHandler*)handler)->GetEvent();
291 fMcEvent = fMcHandler->MCEvent();
295 // array of MC particles
296 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
298 AliFatal("No array of MC particles found !!!");
301 AnalyseCorrectionMode();
303 else AnalyseDataMode();
306 /******************** ANALYSIS METHODS *****************************/
308 //____________________________________________________________________
309 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
311 //Write settings to output list
312 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
313 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
314 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
315 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
316 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
317 settingsTree->Branch("fOnlyOneEtaSide", &fOnlyOneEtaSide,"OnlyOneEtaSide/I");
318 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
319 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
320 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
321 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
322 settingsTree->Branch("fSelectCharge", &fSelectCharge,"SelectCharge/I");
323 settingsTree->Branch("fTriggerSelectCharge", &fTriggerSelectCharge,"TriggerSelectCharge/I");
324 settingsTree->Branch("fTriggerRestrictEta", &fTriggerRestrictEta,"TriggerRestrictEta/D");
325 settingsTree->Branch("fEtaOrdering", &fEtaOrdering,"EtaOrdering/O");
326 settingsTree->Branch("fCutConversions", &fCutConversions,"CutConversions/O");
327 settingsTree->Branch("fCutResonances", &fCutResonances,"CutResonances/O");
328 settingsTree->Branch("fFillpT", &fFillpT,"FillpT/O");
329 settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
330 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
331 settingsTree->Branch("fSkipTrigger", &fSkipTrigger,"SkipTrigger/O");
332 settingsTree->Branch("fInjectedSignals", &fInjectedSignals,"SkipTrigger/O");
333 settingsTree->Branch("fRejectCentralityOutliers", &fRejectCentralityOutliers,"RejectCentralityOutliers/O");
335 settingsTree->Fill();
336 fListOfHistos->Add(settingsTree);
339 //____________________________________________________________________
340 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
342 // Run the analysis on MC to get the correction maps
345 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
347 Double_t centrality = 0;
349 if (fCentralityMethod.Length() > 0)
351 AliCentrality *centralityObj = 0;
353 centralityObj = fAOD->GetHeader()->GetCentralityP();
355 centralityObj = fESD->GetCentrality();
359 centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
360 AliInfo(Form("Centrality is %f", centrality));
364 Printf("WARNING: Centrality object is 0");
369 // Support for ESD and AOD based analysis
370 AliVEvent* inputEvent = fAOD;
378 fHistos->SetRunNumber(inputEvent->GetRunNumber());
379 bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
382 TObject* mc = fArrayMC;
387 fHistos->FillEvent(centrality, -1);
389 // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
390 TObject* vertexSupplier = fMcEvent;
392 vertexSupplier = fAOD->FindListObject(AliAODMCHeader::StdBranchName());
394 if (!fAnalyseUE->VertexSelection(vertexSupplier, 0, fZVertex))
399 zVtx = ((AliAODMCHeader*) vertexSupplier)->GetVtxZ();
401 zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
406 // For productions with injected signals, figure out above which label to skip particles/tracks
407 Int_t skipParticlesAbove = 0;
408 if (fInjectedSignals)
410 AliGenEventHeader* eventHeader = 0;
416 AliHeader* header = (AliHeader*) fMcEvent->Header();
418 AliFatal("fInjectedSignals set but no MC header found");
420 AliGenCocktailEventHeader* cocktailHeader = dynamic_cast<AliGenCocktailEventHeader*> (header->GenEventHeader());
424 AliFatal("fInjectedSignals set but no MC cocktail header found");
427 headers = cocktailHeader->GetHeaders()->GetEntries();
428 eventHeader = dynamic_cast<AliGenEventHeader*> (cocktailHeader->GetHeaders()->First());
433 AliAODMCHeader* header = (AliAODMCHeader*) fAOD->GetList()->FindObject(AliAODMCHeader::StdBranchName());
435 AliFatal("fInjectedSignals set but no MC header found");
437 headers = header->GetNCocktailHeaders();
438 eventHeader = header->GetCocktailHeader(0);
442 AliFatal("First event header not found");
444 skipParticlesAbove = eventHeader->NProduced();
445 AliInfo(Form("Injected signals in this event (%d headers). Keeping events of %s. Will skip particles/tracks above %d.", headers, eventHeader->ClassName(), skipParticlesAbove));
449 TObjArray* tmpList = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
450 if (fInjectedSignals)
451 fAnalyseUE->RemoveInjectedSignals(tmpList, mc, skipParticlesAbove);
452 TObjArray* tracksMC = CloneAndReduceTrackList(tmpList);
458 for (Int_t i=0; i<fArrayMC->GetEntriesFast(); i++)
459 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(((AliAODMCParticle*) fArrayMC->At(i))->PdgCode());
463 for (Int_t i=0; i<fMcEvent->GetNumberOfTracks(); i++)
464 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(fMcEvent->GetTrack(i)->PdgCode());
471 // (MC-true all particles)
473 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, 0, weight);
478 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
480 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
481 for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)
482 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, pool->GetEvent(jMix), 1.0 / pool->GetCurrentNEvents(), (jMix == 0));
483 pool->UpdatePool(CloneAndReduceTrackList(tracksMC));
486 // Printf("trigger: %d", ((AliInputEventHandler*)fInputHandler)->IsEventSelected());
488 // Trigger selection ************************************************
489 if (!fFillOnlyStep0 && (fSkipTrigger || fAnalyseUE->TriggerSelection(fInputHandler)))
491 // (MC-true all particles)
493 if (!fReduceMemoryFootprint)
494 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC, 0, weight);
496 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
499 AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
503 // Vertex selection *************************************************
504 if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
506 // fill here for tracking efficiency
507 // loop over particle species
509 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
511 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
512 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
513 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
515 if (fInjectedSignals)
517 fAnalyseUE->RemoveInjectedSignals(primMCParticles, mc, skipParticlesAbove);
518 fAnalyseUE->RemoveInjectedSignals(primRecoTracksMatched, mc, skipParticlesAbove);
519 fAnalyseUE->RemoveInjectedSignals(allRecoTracksMatched, mc, skipParticlesAbove);
522 fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, particleSpecies, centrality);
524 // Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
526 delete primMCParticles;
527 delete primRecoTracksMatched;
528 delete allRecoTracksMatched;
530 TObjArray* fakeParticles = fAnalyseUE->GetFakeParticles(inputEvent, mc, kFALSE, -1, kTRUE);
531 fHistos->FillTrackingEfficiency(0, 0, 0, (TObjArray*) fakeParticles->At(2), -1, centrality);
532 fHistos->FillFakePt(fakeParticles, centrality);
533 delete fakeParticles;
535 // (MC-true all particles)
537 if (!fReduceMemoryFootprint)
538 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC, 0, weight);
540 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
542 // Get MC primaries that match reconstructed track
543 TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
544 if (fInjectedSignals)
545 fAnalyseUE->RemoveInjectedSignals(tracksRecoMatchedPrim, mc, skipParticlesAbove);
547 // (RECO-matched (quantities from MC particle) primary particles)
549 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim, 0, weight);
551 // Get MC primaries + secondaries that match reconstructed track
552 TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
553 if (fInjectedSignals)
554 fAnalyseUE->RemoveInjectedSignals(tracksRecoMatchedAll, mc, skipParticlesAbove);
556 // (RECO-matched (quantities from MC particle) all particles)
558 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll, 0, weight);
561 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
562 if (fInjectedSignals)
563 fAnalyseUE->RemoveInjectedSignals(tracks, mc, skipParticlesAbove);
568 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
570 // two track cut, STEP 8
571 if (fTwoTrackEfficiencyCut > 0)
572 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
575 if (fFillMixed && !fSkipStep6)
577 AliEventPool* pool2 = fPoolMgr->GetEventPool(centrality, zVtx + 100);
578 //pool2->PrintInfo();
579 if (pool2->IsReady() || pool2->NTracksInPool() > fMixingTracks / 10 || pool2->GetCurrentNEvents() >= 5)
580 for (Int_t jMix=0; jMix<pool2->GetCurrentNEvents(); jMix++)
583 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0));
585 // two track cut, STEP 8
586 if (fTwoTrackEfficiencyCut > 0)
587 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
589 pool2->UpdatePool(CloneAndReduceTrackList(tracks));
592 if (0 && !fReduceMemoryFootprint)
594 // make list of secondaries (matched with MC)
595 TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
596 for (Int_t i=0; i<tracksRecoMatchedAll->GetEntriesFast(); i++)
597 if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
598 tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
600 // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
601 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll, weight);
603 // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
604 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries, weight);
606 // plot delta phi vs process id of secondaries
607 // trigger particles: primaries in 4 < pT < 10
608 // associated particles: secondaries in 1 < pT < 10
610 for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntriesFast(); i++)
612 AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
614 if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
617 for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntriesFast(); j++)
619 AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
621 if (particle->Pt() < 1 || particle->Pt() > 10)
624 if (particle->Pt() > triggerParticle->Pt())
627 Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
628 if (deltaPhi > 1.5 * TMath::Pi())
629 deltaPhi -= TMath::TwoPi();
630 if (deltaPhi < -0.5 * TMath::Pi())
631 deltaPhi += TMath::TwoPi();
633 Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
635 ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
639 delete tracksRecoMatchedSecondaries;
642 delete tracksRecoMatchedPrim;
643 delete tracksRecoMatchedAll;
651 //____________________________________________________________________
652 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
655 // Run the analysis on DATA or MC to get raw distributions
657 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
659 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
664 // skip not selected events here (the AOD is not updated for those)
665 if (!(fInputHandler->IsEventSelected() & fSelectBit))
668 Double_t centrality = 0;
670 AliCentrality *centralityObj = 0;
671 if (fCentralityMethod.Length() > 0)
674 centralityObj = fAOD->GetHeader()->GetCentralityP();
676 centralityObj = fESD->GetCentrality();
679 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
680 //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
689 if (fAOD->GetVZEROData())
692 for (Int_t i=0; i<64; i++)
693 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
697 AliInfo("Rejecting event due to too small V0 multiplicity");
703 AliInfo(Form("Centrality is %f", centrality));
706 // Support for ESD and AOD based analysis
707 AliVEvent* inputEvent = fAOD;
711 Float_t bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
713 fHistos->SetRunNumber(inputEvent->GetRunNumber());
715 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
716 fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
718 // Trigger selection ************************************************
719 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
721 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
722 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
724 // Vertex selection *************************************************
725 if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
727 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
728 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
731 if (centrality < 0 && !fCompareCentralities)
734 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
736 // check for outlier in centrality vs number of tracks (rough constants extracted from correlation histgram)
737 Bool_t reject = kFALSE;
738 if (fRejectCentralityOutliers)
740 if (centrality > 40 && centrality <= 50 && tracks->GetEntriesFast() > 1160)
742 if (centrality > 50 && centrality <= 60 && tracks->GetEntriesFast() > 650)
744 if (centrality > 60 && centrality <= 70 && tracks->GetEntriesFast() > 370)
746 if (centrality > 70 && centrality <= 80 && tracks->GetEntriesFast() > 220)
748 if (centrality > 80 && centrality <= 90 && tracks->GetEntriesFast() > 130)
750 if (centrality > 90 && tracks->GetEntriesFast() > 75)
756 AliInfo(Form("Rejecting event due to centrality vs tracks correlation: %f %d", centrality, tracks->GetEntriesFast()));
757 fHistos->FillEvent(centrality, AliUEHist::kCFStepAnaTopology);
763 // create a list of reduced objects. This speeds up processing and reduces memory consumption for the event pool
764 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
767 //Printf("Accepted %d tracks", tracks->GetEntries());
769 const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
770 Double_t zVtx = vertex->GetZ();
776 // fill two different centralities (syst study)
777 // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
778 if (fCompareCentralities)
781 // Fill containers at STEP 6 (reconstructed)
785 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
787 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
789 if (fTwoTrackEfficiencyCut > 0)
790 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
793 // fill second time with SPD centrality
794 if (fCompareCentralities && centralityObj)
796 centrality = centralityObj->GetCentralityPercentile("CL1");
797 if (centrality >= 0 && !fSkipStep6)
798 fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
805 // 1. First get an event pool corresponding in mult (cent) and
806 // zvertex to the current event. Once initialized, the pool
807 // should contain nMix (reduced) events. This routine does not
808 // pre-scan the chain. The first several events of every chain
809 // will be skipped until the needed pools are filled to the
810 // specified depth. If the pool categories are not too rare, this
811 // should not be a problem. If they are rare, you could lose
814 // 2. Collect the whole pool's content of tracks into one TObjArray
815 // (bgTracks), which is effectively a single background super-event.
817 // 3. The reduced and bgTracks arrays must both be passed into
818 // FillCorrelations(). Also nMix should be passed in, so a weight
819 // of 1./nMix can be applied.
821 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
824 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
828 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
831 Int_t nMix = pool->GetCurrentNEvents();
832 // cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
834 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
835 ((TH2F*) fListOfHistos->FindObject("mixedDist"))->Fill(centrality, pool->NTracksInPool());
837 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(3);
839 // Fill mixed-event histos here
840 for (Int_t jMix=0; jMix<nMix; jMix++)
842 TObjArray* bgTracks = pool->GetEvent(jMix);
845 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0));
847 if (fTwoTrackEfficiencyCut > 0)
848 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
852 // ownership is with the pool now
853 pool->UpdatePool(tracksClone);
860 TObjArray* AliAnalysisTaskPhiCorrelations::CloneAndReduceTrackList(TObjArray* tracks)
862 // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
864 TObjArray* tracksClone = new TObjArray;
865 tracksClone->SetOwner(kTRUE);
867 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
869 AliVParticle* particle = (AliVParticle*) tracks->At(i);
870 tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
876 //____________________________________________________________________
877 void AliAnalysisTaskPhiCorrelations::Initialize()
880 fInputHandler = (AliInputEventHandler*)
881 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
883 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());