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"
53 #include "AliESDZDC.h"
54 #include "AliESDtrackCuts.h"
57 ////////////////////////////////////////////////////////////////////////
59 // Analysis class for azimuthal correlation studies
60 // Based on the UE task from Sara Vallero and Jan Fiete
62 // This class needs input AODs.
63 // The output is a list of analysis-specific containers.
65 // The AOD can be either connected to the InputEventHandler
66 // for a chain of AOD files
68 // to the OutputEventHandler
69 // for a chain of ESD files,
70 // in this case the class should be in the train after the jet-finder
73 // Jan Fiete Grosse-Oetringhaus
75 ////////////////////////////////////////////////////////////////////////
78 ClassImp( AliAnalysisTaskPhiCorrelations )
79 ClassImp( AliDPhiBasicParticle )
81 //____________________________________________________________________
82 AliAnalysisTaskPhiCorrelations:: AliAnalysisTaskPhiCorrelations(const char* name):
83 AliAnalysisTask(name,""),
84 // general configuration
87 fReduceMemoryFootprint(kFALSE),
90 fCompareCentralities(kFALSE),
91 fTwoTrackEfficiencyStudy(kFALSE),
92 fTwoTrackEfficiencyCut(0),
94 fCourseCentralityBinning(kFALSE),
96 fInjectedSignals(kFALSE),
97 // pointers to UE classes
101 fEfficiencyCorrection(0),
102 // handlers and events
110 // histogram settings
113 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
115 fCentralityMethod("V0M"),
121 fSelectBit(AliVEvent::kMB|AliVEvent::kUserDefined),
122 fUseChargeHadrons(kFALSE),
124 fTriggerSelectCharge(0),
125 fTriggerRestrictEta(-1),
126 fEtaOrdering(kFALSE),
127 fCutConversions(kFALSE),
128 fCutResonances(kFALSE),
129 fFillOnlyStep0(kFALSE),
131 fRejectCentralityOutliers(kFALSE),
134 // Default constructor
135 // Define input and output slots here
136 // Input slot #0 works with a TChain
137 DefineInput(0, TChain::Class());
138 // Output slot #0 writes into a TList container
139 DefineOutput(0, TList::Class());
142 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
146 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
147 delete fListOfHistos;
150 //____________________________________________________________________
151 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
154 // Connect the input data
155 if (fDebug > 1) AliInfo("ConnectInputData() ");
157 // Since AODs can either be connected to the InputEventHandler
158 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
159 // we need to get the pointer to the AODEvent correctly.
161 // Delta AODs are also accepted.
163 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
165 if( handler && handler->InheritsFrom("AliAODInputHandler") )
167 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
168 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
172 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
173 if (handler && handler->InheritsFrom("AliAODHandler") )
175 fAOD = ((AliAODHandler*)handler)->GetAOD();
176 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
180 AliWarning("I can't get any AOD Event Handler");
184 if (handler && handler->InheritsFrom("AliESDInputHandler") )
186 // pointer received per event in ::Exec
187 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliESDInputHandler");
190 // Initialize common pointers
194 //____________________________________________________________________
195 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
197 // Create the output container
199 if (fDebug > 1) AliInfo("CreateOutputObjects()");
201 // Initialize class with main algorithms, event and track selection.
202 fAnalyseUE = new AliAnalyseLeadingTrackUE();
203 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, fPtMin);
204 fAnalyseUE->SetDebug(fDebug);
205 fAnalyseUE->DefineESDCuts(fFilterBit);
206 fAnalyseUE->SetEventSelection(fSelectBit);
208 // Initialize output list of containers
209 if (fListOfHistos != NULL){
210 delete fListOfHistos;
211 fListOfHistos = NULL;
214 fListOfHistos = new TList();
215 fListOfHistos->SetOwner(kTRUE);
218 // Initialize class to handle histograms
219 TString histType = "4R";
220 if (fUseVtxAxis == 1)
222 else if (fUseVtxAxis == 2)
224 if (fCourseCentralityBinning)
226 fHistos = new AliUEHistograms("AliUEHistogramsSame", histType);
227 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", histType);
229 fHistos->SetSelectCharge(fSelectCharge);
230 fHistosMixed->SetSelectCharge(fSelectCharge);
232 fHistos->SetSelectTriggerCharge(fTriggerSelectCharge);
233 fHistosMixed->SetSelectTriggerCharge(fTriggerSelectCharge);
235 fHistos->SetTriggerRestrictEta(fTriggerRestrictEta);
236 fHistosMixed->SetTriggerRestrictEta(fTriggerRestrictEta);
238 fHistos->SetOnlyOneEtaSide(fOnlyOneEtaSide);
239 fHistosMixed->SetOnlyOneEtaSide(fOnlyOneEtaSide);
241 fHistos->SetEtaOrdering(fEtaOrdering);
242 fHistosMixed->SetEtaOrdering(fEtaOrdering);
244 fHistos->SetPairCuts(fCutConversions, fCutResonances);
245 fHistosMixed->SetPairCuts(fCutConversions, fCutResonances);
247 fHistos->SetEfficiencyCorrection(fEfficiencyCorrection);
248 fHistosMixed->SetEfficiencyCorrection(fEfficiencyCorrection);
250 // add histograms to list
251 fListOfHistos->Add(fHistos);
252 fListOfHistos->Add(fHistosMixed);
254 fListOfHistos->Add(new TH2F("trackletsVsV0Cent", ";L1 clusters;v0 centrality", 100, -0.5, 9999.5, 101, 0, 101));
255 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));
256 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
257 fListOfHistos->Add(new TH2F("mixedDist", ";centrality;tracks;events", 101, 0, 101, 200, 0, fMixingTracks * 1.5));
258 fListOfHistos->Add(new TH1F("pids", ";pdg;tracks", 2001, -1000.5, 1000.5));
259 fListOfHistos->Add(new TH2F("referenceMultiplicity", ";centrality;tracks;events", 101, 0, 101, 200, 0, 200));
261 PostData(0,fListOfHistos);
263 // Add task configuration to output list
267 Int_t trackDepth = fMixingTracks;
268 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemented of AliEventPoolManager
270 Int_t nCentralityBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(1);
271 Double_t* centralityBins = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(1, 0)->GetXbins()->GetArray();
273 Int_t nZvtxBins = 7+1+7;
274 // bins for second buffer are shifted by 100 cm
275 Double_t vertexBins[] = { -7, -5, -3, -1, 1, 3, 5, 7, 93, 95, 97, 99, 101, 103, 105, 107 };
276 Double_t* zvtxbin = vertexBins;
278 if (fHistos->GetUEHist(2)->GetEventHist()->GetNVar() > 2)
280 nZvtxBins = fHistos->GetUEHist(2)->GetEventHist()->GetNBins(2);
281 zvtxbin = (Double_t*) fHistos->GetUEHist(2)->GetEventHist()->GetAxis(2, 0)->GetXbins()->GetArray();
284 fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
287 //____________________________________________________________________
288 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
290 // receive ESD pointer if we are not running AOD analysis
293 AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
294 if (handler && handler->InheritsFrom("AliESDInputHandler"))
295 fESD = ((AliESDInputHandler*)handler)->GetEvent();
303 fMcEvent = fMcHandler->MCEvent();
307 // array of MC particles
308 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
310 AliFatal("No array of MC particles found !!!");
313 AnalyseCorrectionMode();
315 else AnalyseDataMode();
318 /******************** ANALYSIS METHODS *****************************/
320 //____________________________________________________________________
321 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
323 //Write settings to output list
324 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
325 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
326 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
327 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
328 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
329 settingsTree->Branch("fOnlyOneEtaSide", &fOnlyOneEtaSide,"OnlyOneEtaSide/I");
330 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
331 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
332 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
333 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
334 settingsTree->Branch("fSelectCharge", &fSelectCharge,"SelectCharge/I");
335 settingsTree->Branch("fTriggerSelectCharge", &fTriggerSelectCharge,"TriggerSelectCharge/I");
336 settingsTree->Branch("fTriggerRestrictEta", &fTriggerRestrictEta,"TriggerRestrictEta/D");
337 settingsTree->Branch("fEtaOrdering", &fEtaOrdering,"EtaOrdering/O");
338 settingsTree->Branch("fCutConversions", &fCutConversions,"CutConversions/O");
339 settingsTree->Branch("fCutResonances", &fCutResonances,"CutResonances/O");
340 settingsTree->Branch("fFillpT", &fFillpT,"FillpT/O");
341 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
342 settingsTree->Branch("fSkipTrigger", &fSkipTrigger,"SkipTrigger/O");
343 settingsTree->Branch("fInjectedSignals", &fInjectedSignals,"SkipTrigger/O");
344 settingsTree->Branch("fRejectCentralityOutliers", &fRejectCentralityOutliers,"RejectCentralityOutliers/O");
346 settingsTree->Fill();
347 fListOfHistos->Add(settingsTree);
350 //____________________________________________________________________
351 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
353 // Run the analysis on MC to get the correction maps
356 if ( fDebug > 3 ) AliInfo( " Processing event in Corrections mode ..." );
358 Double_t centrality = 0;
360 if (fCentralityMethod.Length() > 0)
362 AliCentrality *centralityObj = 0;
364 centralityObj = fAOD->GetHeader()->GetCentralityP();
366 centralityObj = fESD->GetCentrality();
370 centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
371 AliInfo(Form("Centrality is %f", centrality));
375 Printf("WARNING: Centrality object is 0");
380 // Support for ESD and AOD based analysis
381 AliVEvent* inputEvent = fAOD;
389 fHistos->SetRunNumber(inputEvent->GetRunNumber());
390 bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
393 TObject* mc = fArrayMC;
398 fHistos->FillEvent(centrality, -1);
403 // Only consider MC events within the vtx-z region used also as cut on the reconstructed vertex
404 TObject* vertexSupplier = fMcEvent;
406 vertexSupplier = fAOD->FindListObject(AliAODMCHeader::StdBranchName());
408 if (!fAnalyseUE->VertexSelection(vertexSupplier, 0, fZVertex))
413 zVtx = ((AliAODMCHeader*) vertexSupplier)->GetVtxZ();
415 zVtx = fMcEvent->GetPrimaryVertex()->GetZ();
420 // For productions with injected signals, figure out above which label to skip particles/tracks
421 Int_t skipParticlesAbove = 0;
422 if (fInjectedSignals)
424 AliGenEventHeader* eventHeader = 0;
430 AliHeader* header = (AliHeader*) fMcEvent->Header();
432 AliFatal("fInjectedSignals set but no MC header found");
434 AliGenCocktailEventHeader* cocktailHeader = dynamic_cast<AliGenCocktailEventHeader*> (header->GenEventHeader());
438 AliFatal("fInjectedSignals set but no MC cocktail header found");
441 headers = cocktailHeader->GetHeaders()->GetEntries();
442 eventHeader = dynamic_cast<AliGenEventHeader*> (cocktailHeader->GetHeaders()->First());
447 AliAODMCHeader* header = (AliAODMCHeader*) fAOD->GetList()->FindObject(AliAODMCHeader::StdBranchName());
449 AliFatal("fInjectedSignals set but no MC header found");
451 headers = header->GetNCocktailHeaders();
452 eventHeader = header->GetCocktailHeader(0);
457 // We avoid AliFatal here, because the AOD productions sometimes have events where the MC header is missing
458 // (due to unreadable Kinematics) and we don't want to loose the whole job because of a few events
459 AliError("First event header not found. Skipping this event.");
460 fHistos->FillEvent(centrality, AliUEHist::kCFStepAnaTopology);
464 skipParticlesAbove = eventHeader->NProduced();
465 AliInfo(Form("Injected signals in this event (%d headers). Keeping events of %s. Will skip particles/tracks above %d.", headers, eventHeader->ClassName(), skipParticlesAbove));
469 TObjArray* tmpList = fAnalyseUE->GetAcceptedParticles(mc, 0, kTRUE, -1, kTRUE);
470 if (fInjectedSignals)
471 fAnalyseUE->RemoveInjectedSignals(tmpList, mc, skipParticlesAbove);
472 TObjArray* tracksMC = CloneAndReduceTrackList(tmpList);
478 for (Int_t i=0; i<fArrayMC->GetEntriesFast(); i++)
479 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(((AliAODMCParticle*) fArrayMC->At(i))->PdgCode());
483 for (Int_t i=0; i<fMcEvent->GetNumberOfTracks(); i++)
484 ((TH1F*) fListOfHistos->FindObject("pids"))->Fill(fMcEvent->GetTrack(i)->PdgCode());
491 // (MC-true all particles)
493 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, 0, weight);
498 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
500 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
501 for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)
502 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepAll, tracksMC, pool->GetEvent(jMix), 1.0 / pool->GetCurrentNEvents(), (jMix == 0));
503 pool->UpdatePool(CloneAndReduceTrackList(tracksMC));
506 // Printf("trigger: %d", ((AliInputEventHandler*)fInputHandler)->IsEventSelected());
508 // Trigger selection ************************************************
509 if (!fFillOnlyStep0 && (fSkipTrigger || fAnalyseUE->TriggerSelection(fInputHandler)))
511 // (MC-true all particles)
513 if (!fReduceMemoryFootprint)
514 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTriggered, tracksMC, 0, weight);
516 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
519 AliFatal("UNEXPECTED: inputEvent is 0. Trigger selection should have failed");
523 // Vertex selection *************************************************
524 if (fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex))
526 // fill here for tracking efficiency
527 // loop over particle species
529 for (Int_t particleSpecies = 0; particleSpecies < 4; particleSpecies++)
531 TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(mc, 0x0, kTRUE, particleSpecies, kTRUE);
532 TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, particleSpecies, kTRUE);
533 TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, particleSpecies, kTRUE);
535 if (fInjectedSignals)
537 fAnalyseUE->RemoveInjectedSignals(primMCParticles, mc, skipParticlesAbove);
538 fAnalyseUE->RemoveInjectedSignals(primRecoTracksMatched, mc, skipParticlesAbove);
539 fAnalyseUE->RemoveInjectedSignals(allRecoTracksMatched, mc, skipParticlesAbove);
542 fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, particleSpecies, centrality);
544 // Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
546 delete primMCParticles;
547 delete primRecoTracksMatched;
548 delete allRecoTracksMatched;
550 TObjArray* fakeParticles = fAnalyseUE->GetFakeParticles(inputEvent, mc, kFALSE, -1, kTRUE);
551 if (fInjectedSignals)
553 fAnalyseUE->RemoveInjectedSignals((TObjArray*) fakeParticles->At(0), mc, skipParticlesAbove);
554 fAnalyseUE->RemoveInjectedSignals((TObjArray*) fakeParticles->At(1), mc, skipParticlesAbove);
556 fHistos->FillTrackingEfficiency(0, 0, 0, (TObjArray*) fakeParticles->At(2), -1, centrality);
557 fHistos->FillFakePt(fakeParticles, centrality);
558 // Printf(">>>>> %d %d %d fakes", ((TObjArray*) fakeParticles->At(0))->GetEntriesFast(), ((TObjArray*) fakeParticles->At(1))->GetEntriesFast(), ((TObjArray*) fakeParticles->At(2))->GetEntriesFast());
559 delete fakeParticles;
561 // (MC-true all particles)
563 if (!fReduceMemoryFootprint)
564 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepVertex, tracksMC, 0, weight);
566 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
568 // Get MC primaries that match reconstructed track
569 TObjArray* tracksRecoMatchedPrim = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kTRUE, -1, kTRUE);
570 if (fInjectedSignals)
571 fAnalyseUE->RemoveInjectedSignals(tracksRecoMatchedPrim, mc, skipParticlesAbove);
573 // (RECO-matched (quantities from MC particle) primary particles)
575 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTrackedOnlyPrim, tracksRecoMatchedPrim, 0, weight);
577 // Get MC primaries + secondaries that match reconstructed track
578 TObjArray* tracksRecoMatchedAll = fAnalyseUE->GetAcceptedParticles(inputEvent, mc, kFALSE, -1, kTRUE);
579 if (fInjectedSignals)
580 fAnalyseUE->RemoveInjectedSignals(tracksRecoMatchedAll, mc, skipParticlesAbove);
582 // (RECO-matched (quantities from MC particle) all particles)
584 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepTracked, tracksRecoMatchedAll, 0, weight);
587 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
588 if (fInjectedSignals)
589 fAnalyseUE->RemoveInjectedSignals(tracks, mc, skipParticlesAbove);
594 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, 0, weight);
596 // two track cut, STEP 8
597 if (fTwoTrackEfficiencyCut > 0)
598 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
601 if (fFillMixed && !fSkipStep6)
603 AliEventPool* pool2 = fPoolMgr->GetEventPool(centrality, zVtx + 100);
604 //pool2->PrintInfo();
605 if (pool2->IsReady() || pool2->NTracksInPool() > fMixingTracks / 10 || pool2->GetCurrentNEvents() >= 5)
606 for (Int_t jMix=0; jMix<pool2->GetCurrentNEvents(); jMix++)
609 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0));
611 // two track cut, STEP 8
612 if (fTwoTrackEfficiencyCut > 0)
613 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracks, pool2->GetEvent(jMix), 1.0 / pool2->GetCurrentNEvents(), (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
615 pool2->UpdatePool(CloneAndReduceTrackList(tracks));
618 if (0 && !fReduceMemoryFootprint)
620 // make list of secondaries (matched with MC)
621 TObjArray* tracksRecoMatchedSecondaries = new TObjArray;
622 for (Int_t i=0; i<tracksRecoMatchedAll->GetEntriesFast(); i++)
623 if (((AliAODMCParticle*)tracksRecoMatchedAll->At(i))->IsPhysicalPrimary() == kFALSE)
624 tracksRecoMatchedSecondaries->Add(tracksRecoMatchedAll->At(i));
626 // Study: Use only secondaries as trigger particles and plot the correlation vs. all particles; store in step 9
627 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy2, tracksRecoMatchedSecondaries, tracksRecoMatchedAll, weight);
629 // Study: Use only primaries as trigger particles and plot the correlation vs. secondaries; store in step 8
630 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksRecoMatchedPrim, tracksRecoMatchedSecondaries, weight);
632 // plot delta phi vs process id of secondaries
633 // trigger particles: primaries in 4 < pT < 10
634 // associated particles: secondaries in 1 < pT < 10
636 for (Int_t i=0; i<tracksRecoMatchedPrim->GetEntriesFast(); i++)
638 AliVParticle* triggerParticle = (AliVParticle*) tracksRecoMatchedPrim->At(i);
640 if (triggerParticle->Pt() < 4 || triggerParticle->Pt() > 10)
643 for (Int_t j=0; j<tracksRecoMatchedSecondaries->GetEntriesFast(); j++)
645 AliAODMCParticle* particle = (AliAODMCParticle*) tracksRecoMatchedSecondaries->At(j);
647 if (particle->Pt() < 1 || particle->Pt() > 10)
650 if (particle->Pt() > triggerParticle->Pt())
653 Double_t deltaPhi = triggerParticle->Phi() - particle->Phi();
654 if (deltaPhi > 1.5 * TMath::Pi())
655 deltaPhi -= TMath::TwoPi();
656 if (deltaPhi < -0.5 * TMath::Pi())
657 deltaPhi += TMath::TwoPi();
659 Int_t processID = fMcEvent->Stack()->Particle(particle->GetLabel())->GetUniqueID();
661 ((TH2F*) fListOfHistos->FindObject("processIDs"))->Fill(deltaPhi, processID);
665 delete tracksRecoMatchedSecondaries;
668 delete tracksRecoMatchedPrim;
669 delete tracksRecoMatchedAll;
677 //____________________________________________________________________
678 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
681 // Run the analysis on DATA or MC to get raw distributions
683 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
685 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
690 // skip not selected events here (the AOD is not updated for those)
691 if (!fSkipTrigger && !(fInputHandler->IsEventSelected() & fSelectBit))
694 // Support for ESD and AOD based analysis
695 AliVEvent* inputEvent = fAOD;
699 Double_t centrality = 0;
701 AliCentrality *centralityObj = 0;
702 if (fCentralityMethod.Length() > 0)
704 if (fCentralityMethod == "ZNA_MANUAL")
707 for(Int_t j = 0; j < 4; ++j) {
708 if (fESD->GetZDCData()->GetZDCTDCData(12,j) != 0) {
713 // Printf("%d %f", zna, fZNAtower[0]);
716 // code from Chiara O (23.10.12)
717 const Double_t *fZNAtower = fESD->GetZDCData()->GetZN2TowerEnergy();
718 Float_t znacut[4] = {681., 563., 413., 191.};
720 if(fZNAtower[0]>znacut[0]) centrality = 1;
721 else if(fZNAtower[0]>znacut[1]) centrality = 21;
722 else if(fZNAtower[0]>znacut[2]) centrality = 41;
723 else if(fZNAtower[0]>znacut[3]) centrality = 61;
724 else centrality = 81;
729 else if (fCentralityMethod == "TRACKS_MANUAL")
732 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
733 centrality = tracks->GetEntriesFast();
736 // Printf("%d %f", tracks->GetEntriesFast(), centrality);
742 centralityObj = fAOD->GetHeader()->GetCentralityP();
744 centralityObj = fESD->GetCentrality();
747 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
748 //centrality = centralityObj->GetCentralityPercentileUnchecked(fCentralityMethod);
757 if (fAOD->GetVZEROData())
760 for (Int_t i=0; i<64; i++)
761 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
765 AliInfo("Rejecting event due to too small V0 multiplicity");
772 AliInfo(Form("Centrality is %f", centrality));
775 Float_t bSign = (inputEvent->GetMagneticField() > 0) ? 1 : -1;
777 fHistos->SetRunNumber(inputEvent->GetRunNumber());
779 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
780 fHistos->FillEvent(centrality, AliUEHist::kCFStepAll);
782 // Trigger selection ************************************************
783 if (!fSkipTrigger && !fAnalyseUE->TriggerSelection(fInputHandler)) return;
785 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
786 fHistos->FillEvent(centrality, AliUEHist::kCFStepTriggered);
788 // Vertex selection *************************************************
789 if(!fAnalyseUE->VertexSelection(inputEvent, fnTracksVertex, fZVertex)) return;
791 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
792 fHistos->FillEvent(centrality, AliUEHist::kCFStepVertex);
795 if (centrality < 0 && !fCompareCentralities)
798 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(inputEvent, 0, kTRUE, -1, kTRUE);
800 // check for outlier in centrality vs number of tracks (rough constants extracted from correlation histgram)
801 Bool_t reject = kFALSE;
802 if (fRejectCentralityOutliers)
804 if (centrality > 40 && centrality <= 50 && tracks->GetEntriesFast() > 1160)
806 if (centrality > 50 && centrality <= 60 && tracks->GetEntriesFast() > 650)
808 if (centrality > 60 && centrality <= 70 && tracks->GetEntriesFast() > 370)
810 if (centrality > 70 && centrality <= 80 && tracks->GetEntriesFast() > 220)
812 if (centrality > 80 && centrality <= 90 && tracks->GetEntriesFast() > 130)
814 if (centrality > 90 && tracks->GetEntriesFast() > 75)
820 AliInfo(Form("Rejecting event due to centrality vs tracks correlation: %f %d", centrality, tracks->GetEntriesFast()));
821 fHistos->FillEvent(centrality, AliUEHist::kCFStepAnaTopology);
826 // reference multiplicity
827 Int_t referenceMultiplicity = -1;
829 referenceMultiplicity = AliESDtrackCuts::GetReferenceMultiplicity(fESD);
831 ((TH2F*) fListOfHistos->FindObject("referenceMultiplicity"))->Fill(centrality, referenceMultiplicity);
833 // create a list of reduced objects. This speeds up processing and reduces memory consumption for the event pool
834 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
837 //Printf("Accepted %d tracks", tracks->GetEntries());
839 const AliVVertex* vertex = inputEvent->GetPrimaryVertex();
840 Double_t zVtx = vertex->GetZ();
846 // fill two different centralities (syst study)
847 // the zvtx axis is used to distinguish the results of both centralities: configured centrality in zvtx = 0, SPD in zvtx = 2
848 if (fCompareCentralities)
851 // Fill containers at STEP 6 (reconstructed)
855 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
857 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
859 if (fTwoTrackEfficiencyCut > 0)
860 fHistos->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, 0, weight, kTRUE, kTRUE, bSign, fTwoTrackEfficiencyCut);
863 // fill second time with SPD centrality
864 if (fCompareCentralities && centralityObj)
866 centrality = centralityObj->GetCentralityPercentile("CL1");
867 if (centrality >= 0 && !fSkipStep6)
868 fHistos->FillCorrelations(centrality, 2, AliUEHist::kCFStepReconstructed, tracksClone, 0, weight);
875 // 1. First get an event pool corresponding in mult (cent) and
876 // zvertex to the current event. Once initialized, the pool
877 // should contain nMix (reduced) events. This routine does not
878 // pre-scan the chain. The first several events of every chain
879 // will be skipped until the needed pools are filled to the
880 // specified depth. If the pool categories are not too rare, this
881 // should not be a problem. If they are rare, you could lose
884 // 2. Collect the whole pool's content of tracks into one TObjArray
885 // (bgTracks), which is effectively a single background super-event.
887 // 3. The reduced and bgTracks arrays must both be passed into
888 // FillCorrelations(). Also nMix should be passed in, so a weight
889 // of 1./nMix can be applied.
891 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
894 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
898 if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5)
901 Int_t nMix = pool->GetCurrentNEvents();
902 // cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
904 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
905 ((TH2F*) fListOfHistos->FindObject("mixedDist"))->Fill(centrality, pool->NTracksInPool());
907 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(3);
909 // Fill mixed-event histos here
910 for (Int_t jMix=0; jMix<nMix; jMix++)
912 TObjArray* bgTracks = pool->GetEvent(jMix);
915 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepReconstructed, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0));
917 if (fTwoTrackEfficiencyCut > 0)
918 fHistosMixed->FillCorrelations(centrality, zVtx, AliUEHist::kCFStepBiasStudy, tracksClone, bgTracks, 1.0 / nMix, (jMix == 0), kTRUE, bSign, fTwoTrackEfficiencyCut);
922 // ownership is with the pool now
923 pool->UpdatePool(tracksClone);
930 TObjArray* AliAnalysisTaskPhiCorrelations::CloneAndReduceTrackList(TObjArray* tracks)
932 // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
934 TObjArray* tracksClone = new TObjArray;
935 tracksClone->SetOwner(kTRUE);
937 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
939 AliVParticle* particle = (AliVParticle*) tracks->At(i);
940 tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
946 //____________________________________________________________________
947 void AliAnalysisTaskPhiCorrelations::Initialize()
950 fInputHandler = (AliInputEventHandler*)
951 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
953 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());