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 "AliAnalysisHelperJetTasks.h"
33 #include "AliAnalysisManager.h"
34 #include "AliAODHandler.h"
35 #include "AliAODInputHandler.h"
36 #include "AliAODMCParticle.h"
37 #include "AliGenPythiaEventHeader.h"
38 #include "AliInputEventHandler.h"
40 #include "AliMCEventHandler.h"
41 #include "AliVParticle.h"
43 #include "AliESDEvent.h"
44 #include "AliESDInputHandler.h"
45 #include "AliMultiplicity.h"
47 #include "EventMixing/AliMixEventInputHandler.h"
50 ////////////////////////////////////////////////////////////////////////
52 // Analysis class for azimuthal correlation studies
53 // Based on the UE task from Sara Vallero and Jan Fiete
55 // This class needs input AODs.
56 // The output is a list of analysis-specific containers.
58 // The AOD can be either connected to the InputEventHandler
59 // for a chain of AOD files
61 // to the OutputEventHandler
62 // for a chain of ESD files,
63 // in this case the class should be in the train after the jet-finder
66 // Jan Fiete Grosse-Oetringhaus
68 ////////////////////////////////////////////////////////////////////////
71 ClassImp( AliAnalysisTaskPhiCorrelations )
73 //____________________________________________________________________
74 AliAnalysisTaskPhiCorrelations:: AliAnalysisTaskPhiCorrelations(const char* name):
75 AliAnalysisTask(name,""),
76 // general configuration
79 fReduceMemoryFootprint(kFALSE),
80 // pointers to UE classes
84 fkTrackingEfficiency(0x0),
85 // handlers and events
94 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
101 fUseChargeHadrons(kFALSE)
103 // Default constructor
104 // Define input and output slots here
105 // Input slot #0 works with a TChain
106 DefineInput(0, TChain::Class());
107 // Output slot #0 writes into a TList container
108 DefineOutput(0, TList::Class());
112 AliAnalysisTaskPhiCorrelations::~AliAnalysisTaskPhiCorrelations()
116 if (fListOfHistos && !AliAnalysisManager::GetAnalysisManager()->IsProofMode())
117 delete fListOfHistos;
120 //____________________________________________________________________
121 void AliAnalysisTaskPhiCorrelations::ConnectInputData(Option_t* /*option*/)
124 // Connect the input data
125 if (fDebug > 1) AliInfo("ConnectInputData() ");
127 // Since AODs can either be connected to the InputEventHandler
128 // or to the OutputEventHandler ( the AOD is created by a previus task in the train )
129 // we need to get the pointer to the AODEvent correctly.
131 // Delta AODs are also accepted.
133 TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
135 if( handler && handler->InheritsFrom("AliAODInputHandler") ) { // input AOD
136 fAOD = ((AliAODInputHandler*)handler)->GetEvent();
137 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
138 } else { //output AOD
139 handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
140 if( handler && handler->InheritsFrom("AliAODHandler") ) {
141 fAOD = ((AliAODHandler*)handler)->GetAOD();
142 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
144 AliFatal("I can't get any AOD Event Handler");
149 // Initialize common pointers
154 //____________________________________________________________________
155 void AliAnalysisTaskPhiCorrelations::CreateOutputObjects()
157 // Create the output container
159 if (fDebug > 1) AliInfo("CreateOutputObjects()");
161 // Initialize class with main algorithms, event and track selection.
162 fAnalyseUE = new AliAnalyseLeadingTrackUE();
163 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, fUseChargeHadrons, fTrackEtaCut, 1.0);
164 fAnalyseUE->SetDebug(fDebug);
165 fAnalyseUE->DefineESDCuts(0);
167 // Initialize output list of containers
168 if (fListOfHistos != NULL){
169 delete fListOfHistos;
170 fListOfHistos = NULL;
173 fListOfHistos = new TList();
174 fListOfHistos->SetOwner(kTRUE);
177 // Initialize class to handle histograms
178 fHistos = new AliUEHistograms("AliUEHistogramsSame", "4");
179 fHistosMixed = new AliUEHistograms("AliUEHistogramsMixed", "4");
181 // add histograms to list
182 fListOfHistos->Add(fHistos);
183 fListOfHistos->Add(fHistosMixed);
185 //fListOfHistos->Add(new TH2F("multVsLeadStep5", ";multiplicity;leading pT", 100, -0.5, 99.5, 20, 0, 10));
187 // Add task configuration to output list
190 PostData(0,fListOfHistos);
193 //____________________________________________________________________
194 void AliAnalysisTaskPhiCorrelations::Exec(Option_t */*option*/)
196 // array of MC particles
198 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
199 if (!fArrayMC)AliFatal("No array of MC particles found !!!");
203 if (fMode) AnalyseCorrectionMode();
204 else AnalyseDataMode();
206 PostData(0,fListOfHistos);
209 /******************** ANALYSIS METHODS *****************************/
211 //____________________________________________________________________
212 void AliAnalysisTaskPhiCorrelations::AddSettingsTree()
214 //Write settings to output list
215 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
216 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
217 settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
218 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
219 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
220 settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
221 settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
222 settingsTree->Fill();
223 fListOfHistos->Add(settingsTree);
226 //____________________________________________________________________
227 void AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
229 // Run the analysis on MC to get the correction maps
234 //____________________________________________________________________
235 void AliAnalysisTaskPhiCorrelations::AnalyseDataMode()
238 // Run the analysis on DATA or MC to get raw distributions
240 PostData(0,fListOfHistos);
242 if ( fDebug > 3 ) AliInfo( " Processing event in Data mode ..." );
245 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
246 fHistos->FillEvent(eventId, AliUEHist::kCFStepAll);
248 // Trigger selection ************************************************
249 if (!fAnalyseUE->TriggerSelection(fInputHandler)) return;
251 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
252 fHistos->FillEvent(eventId, AliUEHist::kCFStepTriggered);
254 // Vertex selection *************************************************
255 if(!fAnalyseUE->VertexSelection(fAOD, fnTracksVertex, fZVertex)) return;
257 // Fill the "event-counting-container", it is needed to get the number of events remaining after each event-selection cut
258 fHistos->FillEvent(eventId, AliUEHist::kCFStepVertex);
260 // TODO HACK centrality should be retrieved from AOD, this needs the ESD
261 AliESDInputHandler* handler = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
262 AliESDEvent* esd = handler->GetEvent();
263 Int_t centrality = esd->GetMultiplicity()->GetNumberOfITSClusters(1);
265 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(fAOD, 0, kTRUE, -1, kTRUE);
267 // reduce to pt and eta range
269 TObjArray* reduced = new TObjArray;
270 for (Int_t i=0; i<tracks->GetEntries(); i++)
272 AliVParticle* particle = (AliVParticle*) tracks->At(i);
273 if (TMath::Abs(particle->Eta()) < fTrackEtaCut && particle->Pt() > fPtMin)
274 reduced->Add(particle);
278 // Fill containers at STEP 6 (reconstructed)
279 fHistos->FillCorrelations(eventId, centrality, AliUEHist::kCFStepReconstructed, tracks);
281 AliInfo("**Main Event");
282 AliInfo(Form(" Tracklets %d", fAOD->GetTracklets()->GetNumberOfTracklets()));
283 AliInfo(Form(" Vz %f", fAOD->GetPrimaryVertex()->GetZ()));
285 AliMixEventInputHandler* mixEH = (AliMixEventInputHandler*) fInputHandler->MixingHandler();
289 AliInfo(Form("num mixed %d", mixEH->MixedEventNumber()));
291 if (mixEH->MixedEventNumber() > 0)
293 // TODO in principle need to apply same quality cuts as above
294 // or somewhere before on the level of the pool (would be more efficient)
295 for (Int_t i = 0; i < mixEH->BufferSize(); i++)
297 AliESDEvent *eventMix = (AliESDEvent*) mixEH->InputEventHandler(i)->GetEvent();
299 AliError(Form("Could not retrieve event %d", i));
303 AliInfo(Form("**Mixed Event %d", i));
304 AliInfo(Form(" Tracklets %d", eventMix->GetMultiplicity()->GetNumberOfTracklets()));
305 AliInfo(Form(" Vz %f", eventMix->GetPrimaryVertex()->GetZ()));
307 TObjArray* tracksMixed = fAnalyseUE->GetAcceptedParticles(eventMix, 0, kTRUE, -1, kTRUE);
309 fHistosMixed->FillCorrelations(eventId, centrality, AliUEHist::kCFStepReconstructed, tracks, tracksMixed);
311 tracksMixed->SetOwner(); // contains tpc only tracks, that is why we have to delete the content as well
323 //____________________________________________________________________
324 void AliAnalysisTaskPhiCorrelations::Initialize()
327 fInputHandler = (AliInputEventHandler*)
328 ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
330 fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());