]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/DiHadronPID/AliAnalysisTaskDiHadronPID.cxx
cd9a109b221fc8f5f970ebce9a187b6f1ca0e12c
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / DiHadronPID / AliAnalysisTaskDiHadronPID.cxx
1 /************************************************************************* 
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. * 
3 *                                                                        * 
4 * Author: The ALICE Off-line Project.                                    * 
5 * Contributors are mentioned in the code where appropriate.              * 
6 *                                                                        * 
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 **************************************************************************/
15
16 // -----------------------------------------------------------------------
17 //  This analysis task fills histograms with PID information of tracks 
18 //  associated to a high p_T trigger.
19 // -----------------------------------------------------------------------
20 //  Author: Misha Veldhoen (misha.veldhoen@cern.ch)
21
22 #include <iostream>
23
24 // Basic Includes
25 #include "TH1F.h"
26 #include "TH2F.h"
27 #include "TH3F.h"
28 #include "THn.h"
29 #include "TFile.h"
30 #include "TChain.h"
31 #include "TObject.h"
32 #include "TObjArray.h"
33
34 // Manager/ Handler
35 #include "AliAnalysisManager.h"
36 #include "AliInputEventHandler.h"
37
38 // Event pool includes.
39 #include "AliEventPoolManager.h"
40
41 // PID includes.
42 #include "AliPIDResponse.h"
43
44 // AOD includes.
45 #include "AliAODEvent.h"
46 #include "AliAODTrack.h"
47 #include "AliAODHandler.h"
48 #include "AliAODVertex.h"
49 #include "AliAODInputHandler.h"
50 #include "AliAODMCParticle.h"
51 #include "AliAODMCHeader.h"
52
53 // Additional includes.
54 #include "AliTrackDiHadronPID.h"
55 #include "AliAODTrackCutsDiHadronPID.h"
56 #include "AliAODEventCutsDiHadronPID.h"
57 #include "AliHistToolsDiHadronPID.h"
58 #include "AliFunctionsDiHadronPID.h"
59
60 // AnalysisTask headers.
61 #include "AliAnalysisTaskSE.h"
62 #include "AliAnalysisTaskDiHadronPID.h"
63
64 using namespace std;
65
66 ClassImp(AliAnalysisTaskDiHadronPID);
67
68 // -----------------------------------------------------------------------
69 AliAnalysisTaskDiHadronPID::AliAnalysisTaskDiHadronPID():
70         AliAnalysisTaskSE(),
71         fPIDResponse(0x0),
72         fEventCuts(0x0),
73         fTrackCutsTrigger(0x0),
74         fTrackCutsAssociated(0x0),
75         fPoolMgr(0x0),
76         fTriggerTracks(0x0),
77         fAssociatedTracks(0x0),
78         fCurrentAODEvent(0x0),
79         fOutputList(0x0),
80         fPtSpectrumTOFbins(0x0),
81         fCorrelationsTOFbins(0x0),
82         fMixedEventsTOFbins(0x0),
83         fPtSpectrumTOFTPCbins(0x0),
84         fCorrelationsTOFTPCbins(0x0),
85         fMixedEventsTOFTPCbins(0x0),
86         fMixedEventsTOFTPCbinsPID(0x0), 
87         fTOFhistos(0x0),
88         fTOFmismatch(0x0),
89         fTOFPtAxis(0x0),
90         fTOFTPChistos(0x0),
91         fTOFTPCmismatch(0x0),
92         fTOFTPCPtAxis(0x0),
93         fNDEtaBins(32),
94         fNDPhiBins(32), 
95         fMinNEventsForMixing(5),
96         fPoolTrackDepth(2000),
97         fPoolSize(1000),
98         fMixEvents(kTRUE),
99         fMixTriggers(kFALSE),
100         fCalculateMismatch(kTRUE),
101         fT0Fill(0x0),
102         fLvsEta(0x0),
103         fLvsEtaProjections(0x0),        
104         fMakeTOFcorrelations(kTRUE),
105         fMakeTOFTPCcorrelationsPi(kFALSE),
106         fMakeTOFTPCcorrelationsKa(kFALSE),
107         fMakeTOFTPCcorrelationsPr(kFALSE),      
108         fTOFIntervalFactorTOFTPC(1.),
109         fExtendPtAxis(kFALSE)
110
111 {
112
113         //
114         // Default Constructor.
115         //
116
117         if (fDebug > 0) {AliInfo("AliAnalysisTaskDiHadronPID Default Constructor.");}           
118
119 }
120
121 // -----------------------------------------------------------------------
122 AliAnalysisTaskDiHadronPID::AliAnalysisTaskDiHadronPID(const char* name):
123         AliAnalysisTaskSE(name),
124         fPIDResponse(0x0),
125         fEventCuts(0x0),
126         fTrackCutsTrigger(0x0),
127         fTrackCutsAssociated(0x0),
128         fPoolMgr(0x0),
129         fTriggerTracks(0x0),
130         fAssociatedTracks(0x0), 
131         fCurrentAODEvent(0x0),
132         fOutputList(0x0),
133         fPtSpectrumTOFbins(0x0),
134         fCorrelationsTOFbins(0x0),
135         fMixedEventsTOFbins(0x0),
136         fPtSpectrumTOFTPCbins(0x0),
137         fCorrelationsTOFTPCbins(0x0),
138         fMixedEventsTOFTPCbins(0x0),
139         fMixedEventsTOFTPCbinsPID(0x0),                 
140         fTOFhistos(0x0),
141         fTOFmismatch(0x0),
142         fTOFPtAxis(0x0),
143         fTOFTPChistos(0x0),
144         fTOFTPCmismatch(0x0),
145         fTOFTPCPtAxis(0x0),     
146         fNDEtaBins(32),
147         fNDPhiBins(32),
148         fMinNEventsForMixing(5),
149         fPoolTrackDepth(2000),
150         fPoolSize(1000),
151         fMixEvents(kTRUE),      
152         fMixTriggers(kFALSE),
153         fCalculateMismatch(kTRUE),
154         fT0Fill(0x0),
155         fLvsEta(0x0),
156         fLvsEtaProjections(0x0),
157         fMakeTOFcorrelations(kTRUE),
158         fMakeTOFTPCcorrelationsPi(kFALSE),
159         fMakeTOFTPCcorrelationsKa(kFALSE),
160         fMakeTOFTPCcorrelationsPr(kFALSE),      
161         fTOFIntervalFactorTOFTPC(1.),
162         fExtendPtAxis(kFALSE)
163
164 {
165
166         //
167         // Named Constructor.
168         //
169
170         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
171
172         DefineInput(0,TChain::Class());
173         DefineOutput(1,TList::Class());
174
175 }
176
177 // -----------------------------------------------------------------------
178 AliAnalysisTaskDiHadronPID::~AliAnalysisTaskDiHadronPID() {
179
180         //
181         // Destructor.
182         //
183
184         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
185
186         if (fPoolMgr) {delete fPoolMgr; fPoolMgr = 0x0;}
187         if (fOutputList) {delete fOutputList; fOutputList = 0x0;}
188
189 }
190
191 // -----------------------------------------------------------------------
192 void AliAnalysisTaskDiHadronPID::UserCreateOutputObjects() {
193
194         //
195         // Create Output objects.
196         //
197
198         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
199
200         AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
201         if (!manager) {AliFatal("Could not obtain analysis manager.");} 
202         AliInputEventHandler* inputHandler = dynamic_cast<AliInputEventHandler*> (manager->GetInputEventHandler());
203         if (!inputHandler) {AliFatal("Could not obtain input handler.");}       
204
205         // Getting the pointer to the PID response object.
206         fPIDResponse = inputHandler->GetPIDResponse();  
207         if (!fPIDResponse) {AliFatal("Could not obtain PID response.");}
208
209         // For now we don't bin in multiplicity for pp.
210         TArrayD* centralityBins = 0x0;
211         if (fEventCuts->GetIsPbPb()) {
212                 Double_t tmp[] = {0., 1., 2., 3., 4., 5., 10., 20., 30., 40., 50., 60., 70., 80., 90., 100.1 };
213                 centralityBins = new TArrayD(15, tmp);
214         } else {
215                 Double_t tmp[] = {0.,1.};
216                 centralityBins = new TArrayD(2, tmp);
217         }
218
219         Int_t nZvtxBins  = 7;
220         Double_t vertexBins[] = {-7., -5., -3., -1., 1., 3., 5., 7.};
221
222         fPoolMgr = new AliEventPoolManager(fPoolSize, fPoolTrackDepth, centralityBins->GetSize(), centralityBins->GetArray(), nZvtxBins, (Double_t*) vertexBins);
223     
224         delete centralityBins;
225
226         // Create the output list.
227         fOutputList = new TList();
228         fOutputList->SetOwner(kTRUE);
229
230         // Creating all requested histograms locally.
231         fEventCuts->CreateHistos();
232         fOutputList->Add(fEventCuts);
233
234         fTrackCutsTrigger->CreateHistos();
235         fOutputList->Add(fTrackCutsTrigger);
236
237         fTrackCutsAssociated->CreateHistos();
238         fOutputList->Add(fTrackCutsAssociated);
239
240         TString speciesname[] = {"Pion","Kaon","Proton"};
241
242         // Create TOF correlations histograms (DPhi,DEta,TOF).
243         if (fMakeTOFcorrelations) {
244
245                 // Get the pT axis for the TOF PID correlations.
246                 Double_t* ptaxis = fTrackCutsAssociated->GetPtAxisPID();
247                 Int_t nptbins = fTrackCutsAssociated->GetNPtBinsPID();
248
249                 // Create Pt spectrum histogram.
250                 fPtSpectrumTOFbins = new TH1F("fPtSpectrumTOFbins","p_{T} Spectrum;p_{T} (GeV/c);Count",nptbins,ptaxis);
251                 fOutputList->Add(fPtSpectrumTOFbins);
252
253                 // Create unidentified correlations histogram.
254                 fCorrelationsTOFbins = AliHistToolsDiHadronPID::MakeHist3D("fCorrelationsTOFbins","Correlations;#Delta#phi;#Delta#eta;p_{T} (GeV/c)",
255                         fNDPhiBins,-TMath::Pi()/2.,3.*TMath::Pi()/2.,
256                         fNDEtaBins,-1.6,1.6,
257                         nptbins, ptaxis);
258                 fOutputList->Add(fCorrelationsTOFbins);
259
260                 // Create unidentified mixed events histogram.
261                 fMixedEventsTOFbins = AliHistToolsDiHadronPID::MakeHist3D("fMixedEventsTOFbins","Mixed Events;#Delta#phi;#Delta#eta;p_{T} (GeV/c)",
262                         fNDPhiBins,-TMath::Pi()/2.,3.*TMath::Pi()/2.,
263                         fNDEtaBins,-1.6,1.6,
264                         nptbins, ptaxis);
265                 fOutputList->Add(fMixedEventsTOFbins);
266
267                 // Create TOFPtaxis.
268                 fTOFPtAxis = new TAxis(nptbins, ptaxis);
269                 fTOFPtAxis->SetName("fTOFPtAxis");
270                 fTOFPtAxis->SetTitle("p_{T} GeV/c");
271
272                 // Create PID histograms.
273                 fTOFhistos = new TObjArray(3);
274                 fTOFhistos->SetOwner(kTRUE);    
275                 fTOFhistos->SetName("CorrelationsTOF");
276
277                 if (fCalculateMismatch) {
278                         fTOFmismatch = new TObjArray(3);
279                         fTOFmismatch->SetOwner(kTRUE);
280                         fTOFmismatch->SetName("MismatchTOF");
281                 }
282
283                 for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
284
285                         TObjArray* TOFhistosTmp = new TObjArray(fTOFPtAxis->GetNbins());
286                         TOFhistosTmp->SetOwner(kTRUE);
287                         TOFhistosTmp->SetName(speciesname[iSpecies].Data());
288
289                         TObjArray* TOFmismatchTmp = 0x0;
290                         if (fCalculateMismatch) {
291                                 TOFmismatchTmp = new TObjArray(fTOFPtAxis->GetNbins());
292                                 TOFmismatchTmp->SetOwner(kTRUE);
293                                 TOFmismatchTmp->SetName(speciesname[iSpecies].Data());
294                         }
295
296                         for (Int_t iBinPt = 1; iBinPt < (fTOFPtAxis->GetNbins() + 1); iBinPt++) {
297
298                                 Int_t iPtClass = fTrackCutsAssociated->GetPtClass(iBinPt);
299                                 if (iPtClass == -1) {AliFatal("Not valid pT class."); continue;}
300
301                                 Int_t NBinsTOF = fTrackCutsAssociated->GetNTOFbins(iPtClass,iSpecies);
302                                 Double_t TOFmin = fTrackCutsAssociated->GetTOFmin(iPtClass,iSpecies);
303                                 Double_t TOFmax = fTrackCutsAssociated->GetTOFmax(iPtClass,iSpecies);
304
305                                 //cout << "ptbin: "<< iBinPt << " class: " << iPtClass << " TOFBins: " << NBinsTOF << " min: " << TOFmin << " max: " << TOFmax << endl; 
306
307                                 // Correlation histogram.
308                                 TH3F* htmp = new TH3F(Form("fCorrelationsTOF_%i",iBinPt),
309                                         Form("%5.3f < p_{T} < %5.3f; #Delta#phi; #Delta#eta; t_{TOF} (ps)", fTOFPtAxis->GetBinLowEdge(iBinPt), fTOFPtAxis->GetBinUpEdge(iBinPt)), 
310                                         fNDPhiBins, -TMath::Pi()/2., 3.*TMath::Pi()/2.,
311                                         fNDEtaBins, -1.6, 1.6, NBinsTOF, TOFmin, TOFmax);
312                                 htmp->SetDirectory(0);
313
314                                 TOFhistosTmp->Add(htmp);
315
316                                 if (fCalculateMismatch) {
317                                         // Mismatch histogram.
318                                         TH1F* htmp2 = new TH1F(Form("fMismatchTOF_%i",iBinPt),
319                                                 Form("%5.3f < p_{T} < %5.3f; t_{TOF} (ps)", fTOFPtAxis->GetBinLowEdge(iBinPt), fTOFPtAxis->GetBinUpEdge(iBinPt)),
320                                                 NBinsTOF, TOFmin, TOFmax);
321                                         htmp2->SetDirectory(0);
322
323                                         TOFmismatchTmp->Add(htmp2);     
324                                 }       
325                         }
326
327                         fTOFhistos->Add(TOFhistosTmp);
328                         if (fCalculateMismatch) {fTOFmismatch->Add(TOFmismatchTmp);}
329
330                 }
331
332                 fOutputList->Add(fTOFhistos);
333                 if (fCalculateMismatch) {fOutputList->Add(fTOFmismatch);}
334
335         }
336
337         // Create TOF/TPC correlation histograms. (DPhi,DEta,TOF,TPC).
338         if (fMakeTOFTPCcorrelationsPi || fMakeTOFTPCcorrelationsKa || fMakeTOFTPCcorrelationsPr) {
339
340                 Double_t ptarrayTOFTPC[16] = {2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 
341                                                                           2.8, 3.0, 3.2, 3.4, 3.6, 3.8, 
342                                                                           4.2, 4.6, 5.0};
343                 const Int_t nptbins = 15;
344                 Double_t ptarrayTOFTPCext[26] = {1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 
345                                                                           2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 
346                                                                           2.8, 3.0, 3.2, 3.4, 3.6, 3.8, 
347                                                                           4.2, 4.6, 5.0};
348                 const Int_t nptbinsext = 25;
349
350                 fTOFTPCPtAxis = new TAxis(fExtendPtAxis ? nptbinsext : nptbins, fExtendPtAxis ? ptarrayTOFTPCext : ptarrayTOFTPC);
351                 fTOFTPCPtAxis->SetName("fTOFTPCPtAxis");
352                 fTOFTPCPtAxis->SetTitle("p_{T} GeV/c");
353
354                 // Create Pt spectrum histogram.
355                 fPtSpectrumTOFTPCbins = new TH1F("fPtSpectrumTOFTPCbins","p_{T} Spectrum;p_{T} (GeV/c);Count",nptbins,ptarrayTOFTPC);
356                 fOutputList->Add(fPtSpectrumTOFTPCbins);
357
358                 // Create unidentified correlations histogram.
359                 fCorrelationsTOFTPCbins = AliHistToolsDiHadronPID::MakeHist3D("fCorrelationsTOFTPCbins","Correlations;#Delta#phi;#Delta#eta;p_{T} (GeV/c)",
360                         fNDPhiBins,-TMath::Pi()/2.,3.*TMath::Pi()/2.,
361                         fNDEtaBins,-1.6,1.6,fExtendPtAxis ? nptbinsext : nptbins, fExtendPtAxis ? ptarrayTOFTPCext : ptarrayTOFTPC);
362                 fOutputList->Add(fCorrelationsTOFTPCbins);
363
364                 // Create unidentified mixed events histogram.
365                 fMixedEventsTOFTPCbins = AliHistToolsDiHadronPID::MakeHist3D("fMixedEventsTOFTPCbins","Mixed Events;#Delta#phi;#Delta#eta;p_{T} (GeV/c)",
366                         fNDPhiBins,-TMath::Pi()/2.,3.*TMath::Pi()/2.,
367                         fNDEtaBins,-1.6,1.6,fExtendPtAxis ? nptbinsext : nptbins, fExtendPtAxis ? ptarrayTOFTPCext : ptarrayTOFTPC);
368                 fOutputList->Add(fMixedEventsTOFTPCbins);
369
370                 fTOFTPChistos = new TObjArray(3);
371                 fTOFTPChistos->SetOwner(kTRUE);
372                 fTOFTPChistos->SetName("CorrelationsTOFTPC");
373
374                 if (fCalculateMismatch) {
375                         fTOFTPCmismatch = new TObjArray(3);
376                         fTOFTPCmismatch->SetOwner(kTRUE);
377                         fTOFTPCmismatch->SetName("MismatchTOFTPC");
378                 }
379
380                 fMixedEventsTOFTPCbinsPID = new TObjArray(3);
381                 fMixedEventsTOFTPCbinsPID->SetOwner(kTRUE);
382                 fMixedEventsTOFTPCbinsPID->SetName("MixedEventsTOFTPC");
383
384                 for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
385
386                         // Create Mixed events with PID.
387                         TH3F* mixedeventsPID = AliHistToolsDiHadronPID::MakeHist3D(Form("fMixedEventsTOFTPC%s", speciesname[iSpecies].Data()),
388                         Form("Mixed Events %s;#Delta#phi;#Delta#eta;p_{T} (GeV/c)", speciesname[iSpecies].Data()),
389                         fNDPhiBins,-TMath::Pi()/2.,3.*TMath::Pi()/2.,
390                         fNDEtaBins,-1.6,1.6,fExtendPtAxis ? nptbinsext : nptbins, fExtendPtAxis ? ptarrayTOFTPCext : ptarrayTOFTPC);
391                         fMixedEventsTOFTPCbinsPID->Add(mixedeventsPID);
392                 
393                         // Create the directory structure Pion, Kaon, Proton, regardless
394                         // of wether the histograms are created (to keep the order.)
395                         TObjArray* TOFTPChistosTmp = new TObjArray(fTOFTPCPtAxis->GetNbins());
396                         TOFTPChistosTmp->SetOwner(kTRUE);
397                         TOFTPChistosTmp->SetName(speciesname[iSpecies].Data());
398
399                         TObjArray* TOFTPCmismatchTmp = 0x0;
400                         if (fCalculateMismatch) { 
401                                 TOFTPCmismatchTmp = new TObjArray(fTOFTPCPtAxis->GetNbins());
402                                 TOFTPCmismatchTmp->SetOwner(kTRUE);
403                                 TOFTPCmismatchTmp->SetName(speciesname[iSpecies].Data());       
404                         }
405
406                         // Only Create the TOF/TPC histograms when requested.
407                         Bool_t MakeTOFTPCcorrelations[3] = {fMakeTOFTPCcorrelationsPi, fMakeTOFTPCcorrelationsKa, fMakeTOFTPCcorrelationsPr};
408                         if (MakeTOFTPCcorrelations[iSpecies]) {
409                                 for (Int_t iBinPt = 1; iBinPt < (fTOFTPCPtAxis->GetNbins() + 1); iBinPt++) {
410                         
411                                         // Approximate resolutions of TOF and TPC detector.
412                                         const Double_t sTOFest = 110.;
413                                         const Double_t sTPCest = 4.5;
414                                         
415                                         // Set range +/- 5 sigma of main peak. (+ 10 sigma for TOF max, for mismatches.)
416                                         Double_t TOFmin = -5. * sTOFest;
417                                         Double_t TOFmax = 10. * sTOFest;
418                                         Double_t TPCmin = -4. * sTPCest;
419                                         Double_t TPCmax = 4. * sTPCest;
420
421                                         Double_t TOFexp = AliFunctionsDiHadronPID::TOFExpTime(fTOFTPCPtAxis->GetBinLowEdge(iBinPt), 0.4,  AliFunctionsDiHadronPID::M(iSpecies));
422                                         Double_t TPCexp = AliFunctionsDiHadronPID::TPCExpdEdX(fTOFTPCPtAxis->GetBinLowEdge(iBinPt), 0.4,  AliFunctionsDiHadronPID::M(iSpecies));
423
424                                         for (Int_t jSpecies = 0; jSpecies < 3; jSpecies++) {
425
426                                                 if (iSpecies == jSpecies) {continue;}
427
428                                                 Double_t TOFexpOther = AliFunctionsDiHadronPID::TOFExpTime(fTOFTPCPtAxis->GetBinLowEdge(iBinPt), 0.4,  AliFunctionsDiHadronPID::M(jSpecies));
429                                                 Double_t TPCexpOther = AliFunctionsDiHadronPID::TPCExpdEdX(fTOFTPCPtAxis->GetBinLowEdge(iBinPt), 0.4,  AliFunctionsDiHadronPID::M(jSpecies));
430                                         
431                                                 // If any peak is within +/- 7 sigma, then also add this peak.
432                                                 if ( (TMath::Abs(TOFexp - TOFexpOther) < 7. * sTOFest) ||
433                                                          (TMath::Abs(TPCexp - TPCexpOther) < 7. * sTPCest) ) {
434
435                                                         TOFmin = TMath::Min(TOFmin, (TOFexpOther - TOFexp - 2. * sTOFest) );
436                                                         TOFmax = TMath::Max(TOFmax, (TOFexpOther - TOFexp + 10. * sTOFest) );
437                                                         TPCmin = TMath::Min(TPCmin, (TPCexpOther - TPCexp - 2. * sTPCest) );
438                                                         TPCmax = TMath::Max(TPCmax, (TPCexpOther - TPCexp + 2. * sTPCest) );                                            
439
440                                                 }
441
442                                         }
443
444                                         // With the standard TOF range, fitting the deuterons and the TOF mismatches is
445                                         // hard. This flag doubles the range of the TOF axis in the TOF/TPC histograms,
446                                         // while leaving the resolution the same. Turning on this flag will greatly increase
447                                         // the memory consumption of the task, to the point that it's probably too much 
448                                         // to save a Buffer with all three species included.
449                                         Double_t TOFreach = TOFmax - TOFmin;
450                                         TOFmax += (TOFreach * (fTOFIntervalFactorTOFTPC - 1.));
451                                         Int_t TOFbins = (Int_t)(60. * fTOFIntervalFactorTOFTPC);
452
453                                         Int_t NBinsTOFTPC[4] = {32, 32, TOFbins, 40};
454                                         Double_t minTOFTPC[4] = {-TMath::Pi()/2., -1.6, TOFmin, TPCmin};
455                                         Double_t maxTOFTPC[4] = {3.*TMath::Pi()/2., 1.6, TOFmax, TPCmax};
456
457                                         THnF* htmp = new THnF(Form("fCorrelationsTOFTPC_%i",iBinPt),
458                                                 Form("%5.3f < p_{T} < %5.3f", fTOFTPCPtAxis->GetBinLowEdge(iBinPt), fTOFTPCPtAxis->GetBinUpEdge(iBinPt)), 
459                                                 4, NBinsTOFTPC, minTOFTPC, maxTOFTPC);
460
461                                         (htmp->GetAxis(0))->SetTitle("#Delta#phi");
462                                         (htmp->GetAxis(1))->SetTitle("#Delta#eta");
463                                         (htmp->GetAxis(2))->SetTitle("t_{TOF} (ps)");
464                                         (htmp->GetAxis(3))->SetTitle("dE/dx (a.u.)");
465
466                                         TOFTPChistosTmp->Add(htmp);
467
468                                         if (fCalculateMismatch) { 
469                                                 // Mismatch histogram.
470                                                 TH2F* htmp2 = new TH2F(Form("fMismatchTOFTPC_%i",iBinPt),
471                                                         Form("%5.3f < p_{T} < %5.3f; t_{TOF} (ps); dE/dx (a.u.)", fTOFTPCPtAxis->GetBinLowEdge(iBinPt), fTOFTPCPtAxis->GetBinUpEdge(iBinPt)), 
472                                                         NBinsTOFTPC[2], TOFmin, TOFmax, NBinsTOFTPC[3], TPCmin, TPCmax);
473                                                 htmp2->SetDirectory(0);
474
475                                                 TOFTPCmismatchTmp->Add(htmp2);
476                 
477                                         }
478                                 } // End loop over pT bins.
479                         } // End species if.
480
481                         fTOFTPChistos->Add(TOFTPChistosTmp);
482                         if (fCalculateMismatch) {fTOFTPCmismatch->Add(TOFTPCmismatchTmp);}
483
484                 }
485
486                 fOutputList->Add(fTOFTPChistos);
487                 if (fCalculateMismatch) {fOutputList->Add(fTOFTPCmismatch);}
488                 fOutputList->Add(fMixedEventsTOFTPCbinsPID);
489
490         }
491
492         // Load external TOF histograms if flag is set.
493         if (fCalculateMismatch) {LoadExtMismatchHistos();}
494
495         PostData(1,fOutputList);
496
497 }
498
499 // -----------------------------------------------------------------------
500 void AliAnalysisTaskDiHadronPID::LocalInit() {
501
502         //
503         // Initialize on the this computer. 
504         //
505
506         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
507
508 }
509
510 // -----------------------------------------------------------------------
511 void AliAnalysisTaskDiHadronPID::UserExec(Option_t*) {
512
513         //
514         // Main Loop.
515         //
516
517         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
518
519         // Input Current Event.
520         fCurrentAODEvent = dynamic_cast<AliAODEvent*>(InputEvent());
521         if (!fCurrentAODEvent) AliFatal("No Event Found!");
522
523         if (!fEventCuts->IsSelected(fCurrentAODEvent)) {return;}
524
525         // Fill the global tracks array. - NOT NEEDED I THINK, since we're not using
526         // bit 1<<7 for the associated tracks!
527
528         // Let the track cut objects know that a new event will start.
529         fTrackCutsTrigger->StartNewEvent();
530         fTrackCutsAssociated->StartNewEvent();
531
532         // Create arrays for trigger/associated particles.
533         fTriggerTracks = new TObjArray(350);
534         fTriggerTracks->SetOwner(kTRUE);
535
536         fAssociatedTracks = new TObjArray(3500);
537         fAssociatedTracks->SetOwner(kTRUE);
538
539         for (Int_t iTrack = 0; iTrack < fCurrentAODEvent->GetNumberOfTracks(); iTrack++) {
540
541                 AliAODTrack* track = (AliAODTrack*)fCurrentAODEvent->GetTrack(iTrack);
542                 AliTrackDiHadronPID* pidtrack = new AliTrackDiHadronPID(track,0x0,0x0,fPIDResponse);
543                 pidtrack->ForgetAboutPointers();
544                 pidtrack->SetDebugLevel(fDebug);
545
546                 Double_t rndhittime = -1.e21;
547                 if (fCalculateMismatch) rndhittime = GenerateRandomHit(pidtrack->Eta());
548
549                 // Fill the trigger/associated tracks array.
550                 if (fTrackCutsTrigger->IsSelectedData(pidtrack,rndhittime)) {fTriggerTracks->AddLast(pidtrack);}
551                 else if (fTrackCutsAssociated->IsSelectedData(pidtrack,rndhittime)) {
552                         
553                         fAssociatedTracks->AddLast(pidtrack);
554
555                         // Fill p_T spectrum.
556                         if (fPtSpectrumTOFbins) fPtSpectrumTOFbins->Fill(pidtrack->Pt());
557                         if (fPtSpectrumTOFTPCbins) fPtSpectrumTOFTPCbins->Fill(pidtrack->Pt());
558
559                         // Fill mismatch histograms with associateds.
560                         if (fCalculateMismatch && (rndhittime > -1.e20)) {
561
562                                 Double_t apt = pidtrack->Pt();
563
564                                 if (fMakeTOFcorrelations) {
565
566                                         for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
567
568                                                 TObjArray* atmp = (TObjArray*)fTOFmismatch->At(iSpecies);
569                                                 Int_t ptbin = fTOFPtAxis->FindBin(apt);
570
571                                                 // Only fill if histogram exists in fTOFmismatch.
572                                                 if ( !(ptbin < 1) && !(ptbin > fTOFPtAxis->GetNbins()) ) {
573
574                                                         TH1F* htmp = (TH1F*)atmp->At(ptbin - 1);
575                                                         htmp->Fill(rndhittime - pidtrack->GetTOFsignalExpected(iSpecies));
576
577                                                 }
578                                         }
579                                 }
580
581                                 Bool_t MakeTOFTPCcorrelations[3] = {fMakeTOFTPCcorrelationsPi, fMakeTOFTPCcorrelationsKa, fMakeTOFTPCcorrelationsPr};
582                                 if (fMakeTOFTPCcorrelationsPi || fMakeTOFTPCcorrelationsKa || fMakeTOFTPCcorrelationsPr) { 
583
584                                         for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
585
586                                                 if (!MakeTOFTPCcorrelations[iSpecies]) {continue;}
587
588                                                 TObjArray* atmp = (TObjArray*)fTOFTPCmismatch->At(iSpecies);
589                                                 Int_t ptbin = fTOFTPCPtAxis->FindBin(apt);
590
591                                                 // Only fill if histogram exists in fTOFTPCmismatch.
592                                                 if ( !(ptbin < 1) && !(ptbin > fTOFTPCPtAxis->GetNbins()) ) {
593
594                                                         TH2F* htmp = (TH2F*)atmp->At(ptbin - 1);
595                                                         htmp->Fill(rndhittime - pidtrack->GetTOFsignalExpected(iSpecies), pidtrack->GetTPCsignalMinusExpected(iSpecies));
596
597                                                 }
598                                         }
599                                 }
600
601                         }
602
603                 } 
604                 else {delete pidtrack;}
605
606         }
607
608         // Fill Correlation histograms.
609         for (Int_t iTrigger = 0; iTrigger < fTriggerTracks->GetEntriesFast(); iTrigger++) {
610                 AliTrackDiHadronPID* triggertrack = (AliTrackDiHadronPID*)fTriggerTracks->At(iTrigger);
611
612                 for (Int_t iAssociated = 0; iAssociated < fAssociatedTracks->GetEntriesFast(); iAssociated++) {
613                         AliTrackDiHadronPID* associatedtrack = (AliTrackDiHadronPID*)fAssociatedTracks->At(iAssociated);
614
615                         Double_t DPhi = triggertrack->Phi() - associatedtrack->Phi();
616                         if (DPhi < -TMath::Pi()/2.) {DPhi += 2.*TMath::Pi();}
617                         if (DPhi > 3.*TMath::Pi()/2.) {DPhi -= 2.*TMath::Pi();}
618
619                         Double_t DEta = triggertrack->Eta() - associatedtrack->Eta();
620                         if (fCorrelationsTOFbins) fCorrelationsTOFbins->Fill(DPhi,DEta,associatedtrack->Pt());
621                         if (fCorrelationsTOFTPCbins) fCorrelationsTOFTPCbins->Fill(DPhi,DEta,associatedtrack->Pt());
622
623                         Double_t apt = associatedtrack->Pt();
624
625                         // Fill TOF correlations.
626                         if (fMakeTOFcorrelations) {
627                                 
628                                 for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
629
630                                         TObjArray* atmp = (TObjArray*)fTOFhistos->At(iSpecies);
631                                         Int_t ptbin = fTOFPtAxis->FindBin(apt);
632
633                                         // Only fill if histogram exists in fTOFhistos.
634                                         if ( !(ptbin < 1) && !(ptbin > fTOFPtAxis->GetNbins()) ) {
635
636                                                 TH3F* htmp = (TH3F*)atmp->At(ptbin - 1);
637                                                 htmp->Fill(DPhi, DEta, associatedtrack->GetTOFsignalMinusExpected(iSpecies));
638
639                                         }
640                                 }
641                         }
642
643                         // Fill TOF/ TPC Correlations.
644                         Bool_t MakeTOFTPCcorrelations[3] = {fMakeTOFTPCcorrelationsPi, fMakeTOFTPCcorrelationsKa, fMakeTOFTPCcorrelationsPr};
645                         if (fMakeTOFTPCcorrelationsPi || fMakeTOFTPCcorrelationsKa || fMakeTOFTPCcorrelationsPr) { 
646
647                                 for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
648
649                                         if (!MakeTOFTPCcorrelations[iSpecies]) {continue;}
650
651                                         TObjArray* atmp = (TObjArray*)fTOFTPChistos->At(iSpecies);
652                                         Int_t ptbin = fTOFTPCPtAxis->FindBin(apt);
653
654                                         // Only fill if histogram exists in fTOFhistos.
655                                         if ( !(ptbin < 1) && !(ptbin > fTOFTPCPtAxis->GetNbins()) ) {
656
657                                                 THnF* htmp = (THnF*)atmp->At(ptbin - 1);
658                                                 Double_t TOFTPCfill[4] = {DPhi, DEta, 
659                                                         associatedtrack->GetTOFsignalMinusExpected(iSpecies), associatedtrack->GetTPCsignalMinusExpected(iSpecies)}; 
660
661                                                 htmp->Fill(TOFTPCfill);
662
663                                         }                               
664                                 }       
665                         }               
666                 }
667         }
668
669         //cout<<"Triggers: "<<fTriggerTracks->GetEntriesFast()<<" Associateds: "<<fAssociatedTracks->GetEntriesFast()<<endl;    
670
671         // Determine vtxz of current event.
672         AliAODVertex* currentprimaryvertex = fCurrentAODEvent->GetPrimaryVertex();
673         Double_t vtxz = currentprimaryvertex->GetZ();
674
675         // Determine centrality of current event (for PbPb).
676         AliEventPool* poolin = 0x0;
677         Float_t percentile = -1.;
678         if (fEventCuts->GetIsPbPb()) {
679                 TString centralityestimator = fEventCuts->GetCentralityEstimator();
680                 AliCentrality* currentcentrality = fCurrentAODEvent->GetCentrality();
681                 percentile = currentcentrality->GetCentralityPercentile(centralityestimator.Data());
682
683                 poolin = fPoolMgr->GetEventPool(percentile, vtxz); 
684                 if (!poolin) {AliFatal(Form("No pool found for centrality = %f, vtxz = %f", percentile, vtxz));}
685         } else {
686                 poolin = fPoolMgr->GetEventPool(0.5, vtxz);     // There are no multiplicity bins for pp yet.  
687                 if (!poolin) {AliFatal(Form("No pool found for vtxz = %f", vtxz));}
688         }
689
690         // TObjArray* fGlobalTracksArray; 
691
692         // Give a print out of the pool manager's contents.
693         if (fDebug > 0) PrintPoolManagerContents();
694
695         // Mix events if there are enough events in the pool.
696         if (poolin->GetCurrentNEvents() >= fMinNEventsForMixing) {
697                 //{cout << "Mixing Events." << endl;}
698
699                 // Loop over all events in the event pool.
700                 for (Int_t iMixEvent = 0; iMixEvent < poolin->GetCurrentNEvents(); iMixEvent++) {
701                 TObjArray* mixtracks = poolin->GetEvent(iMixEvent);
702
703                 // Mix either the triggers or the associateds.
704                 if (fMixTriggers) {
705
706                                 // Loop over all associateds in this event.
707                                 for (Int_t iAssociated = 0; iAssociated < fAssociatedTracks->GetEntriesFast(); iAssociated++) {
708                                         AliTrackDiHadronPID* associatedtrack = (AliTrackDiHadronPID*)fAssociatedTracks->At(iAssociated);
709
710                                 // Loop over all mixed tracks.
711                                 for (Int_t iMixTrack = 0; iMixTrack < mixtracks->GetEntriesFast(); iMixTrack++) {
712                                         AliTrackDiHadronPID* mixtrack = (AliTrackDiHadronPID*)mixtracks->At(iMixTrack);
713                                                         
714                                                 Double_t DPhi = mixtrack->Phi() - associatedtrack->Phi();
715                                                 if (DPhi < -TMath::Pi()/2.) {DPhi += 2.*TMath::Pi();}
716                                                 if (DPhi > 3.*TMath::Pi()/2.) {DPhi -= 2.*TMath::Pi();}
717
718                                                 Double_t DEta = mixtrack->Eta() - associatedtrack->Eta();
719                                                 if (fMixedEventsTOFbins) fMixedEventsTOFbins->Fill(DPhi,DEta,associatedtrack->Pt());
720                                                 if (fMixedEventsTOFTPCbins) fMixedEventsTOFTPCbins->Fill(DPhi,DEta,associatedtrack->Pt());
721
722                                                 // Fill the mixed event histograms with a 1 sigma PID cut.
723                                                 if (fMixedEventsTOFTPCbinsPID) {
724
725                                                         for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
726
727                                                                 TH3F* mixedeventhist = (TH3F*)fMixedEventsTOFTPCbinsPID->At(iSpecies);
728
729                                                                 // Check the nSigma of the associated tracks.
730                                                                 Double_t nSigmaTOFTPC = TMath::Sqrt( 
731                                                                 associatedtrack->GetNumberOfSigmasTOF(iSpecies) * associatedtrack->GetNumberOfSigmasTOF(iSpecies) +
732                                                                 associatedtrack->GetNumberOfSigmasTPC(iSpecies) * associatedtrack->GetNumberOfSigmasTPC(iSpecies));
733
734                                                                 if (nSigmaTOFTPC < 1.) {mixedeventhist->Fill(DPhi,DEta,associatedtrack->Pt());}
735
736                                                         }
737                                                 }
738
739                                 }
740                                 }
741
742                         } else {
743
744                                 // Loop over all triggers in this event.
745                                 for (Int_t iTrigger = 0; iTrigger < fTriggerTracks->GetEntriesFast(); iTrigger++) {
746                                         AliTrackDiHadronPID* triggertrack = (AliTrackDiHadronPID*)fTriggerTracks->At(iTrigger);
747
748                                 // Loop over all mixed tracks.
749                                 for (Int_t iMixTrack = 0; iMixTrack < mixtracks->GetEntriesFast(); iMixTrack++) {
750                                         AliTrackDiHadronPID* mixtrack = (AliTrackDiHadronPID*)mixtracks->At(iMixTrack);
751                                                         
752                                                 Double_t DPhi = triggertrack->Phi() - mixtrack->Phi();
753                                                 if (DPhi < -TMath::Pi()/2.) {DPhi += 2.*TMath::Pi();}
754                                                 if (DPhi > 3.*TMath::Pi()/2.) {DPhi -= 2.*TMath::Pi();}
755
756                                                 Double_t DEta = triggertrack->Eta() - mixtrack->Eta();
757                                                 if (fMixedEventsTOFbins) fMixedEventsTOFbins->Fill(DPhi,DEta,mixtrack->Pt());
758                                                 if (fMixedEventsTOFTPCbins) fMixedEventsTOFTPCbins->Fill(DPhi,DEta,mixtrack->Pt());
759                                 
760                                                 // Fill the mixed event histograms with a 1 sigma PID cut.
761                                                 if (fMixedEventsTOFTPCbinsPID) {
762
763                                                         for (Int_t iSpecies = 0; iSpecies < 3; iSpecies++) {
764
765                                                                 TH3F* mixedeventhist = (TH3F*)fMixedEventsTOFTPCbinsPID->At(iSpecies);
766
767                                                                 // Check the nSigma of the associated tracks.
768                                                                 Double_t nSigmaTOFTPC = TMath::Sqrt( 
769                                                                 mixtrack->GetNumberOfSigmasTOF(iSpecies) * mixtrack->GetNumberOfSigmasTOF(iSpecies) +
770                                                                 mixtrack->GetNumberOfSigmasTPC(iSpecies) * mixtrack->GetNumberOfSigmasTPC(iSpecies));
771
772                                                                 if (nSigmaTOFTPC < 1.) {mixedeventhist->Fill(DPhi,DEta,mixtrack->Pt());}
773
774                                                         }
775                                                 }
776
777                                 }
778                                 }
779
780                         } // End if     
781                 }
782         }       
783
784         // Update the event pool.
785         AliEventPool* poolout = 0x0;
786         if (fEventCuts->GetIsPbPb()) {
787                 poolout = fPoolMgr->GetEventPool(percentile, vtxz); // Get the buffer associated with the current centrality and z-vtx
788                 if (!poolout) AliFatal(Form("No pool found for centrality = %f, vtx_z = %f", percentile, vtxz));
789         } else {
790                 poolout = fPoolMgr->GetEventPool(0.5, vtxz); // Get the buffer associated with the current centrality and z-vtx
791                 if (!poolout) AliFatal(Form("No pool found for vtx_z = %f", vtxz));
792         }
793
794
795         // Q: is it a problem that the fAssociatedTracks array can be bigger than the number of tracks inside?
796         if (fMixTriggers) {
797                 poolout->UpdatePool(fTriggerTracks);
798                 fAssociatedTracks->Delete();
799                 delete fAssociatedTracks;
800         }
801         else {
802                 poolout->UpdatePool(fAssociatedTracks);
803                 fTriggerTracks->Delete();
804                 delete fTriggerTracks;
805         }
806
807         fTriggerTracks = 0x0;
808         fAssociatedTracks = 0x0;
809
810         // Tell the track cut object that the event is done.
811         fTrackCutsTrigger->EventIsDone(0);
812         fTrackCutsAssociated->EventIsDone(0);
813
814         PostData(1,fOutputList);
815
816 }
817
818 // -----------------------------------------------------------------------
819 void AliAnalysisTaskDiHadronPID::SelectCollisionCandidates(UInt_t offlineTriggerMask) {
820
821         // Overrides the method defined in AliAnalysisTaskSE. This is needed because
822         // the event selection is not done in the task, but in the AliAODEventCutsDiHadronPID class.
823
824         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
825         if (!fEventCuts) {cout << Form("%s -> ERROR: No AliAODEventCutsDiHadronPID class created for the analysis...",__func__) << endl; return;}
826
827         //fOfflineTriggerMask = offlineTriggerMask;
828         fEventCuts->SetTrigger(offlineTriggerMask);
829
830 }
831
832 // -----------------------------------------------------------------------
833 Bool_t AliAnalysisTaskDiHadronPID::LoadExtMismatchHistos() {
834
835         //
836         // Attempting to load a root file containing information needed
837         // to generate random TOF hits.
838         //
839
840         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
841
842         // Opening external TOF file.
843         if (fDebug > 0) cout<<"Trying to open TOFmismatchHistos.root ..."<<endl;
844         TFile* fin = 0x0;
845         fin = TFile::Open("alien:///alice/cern.ch/user/m/mveldhoe/rootfiles/TOFmismatchHistos.root");
846         if (!fin) {
847                 AliWarning("Couln't open TOFmismatchHistos, will not calculate mismatches...");
848                 fCalculateMismatch = kFALSE;
849                 return kFALSE;
850         }
851
852         // Check if the required histograms are present.
853         TH1F* tmp1 = (TH1F*)fin->Get("hNewT0Fill");
854         if (!tmp1) {
855                 AliWarning("Couln't find hNewT0Fill, will not calculate mismatches...");
856                 fCalculateMismatch = kFALSE;
857                 return kFALSE;  
858         }
859         TH2F* tmp2 = (TH2F*)fin->Get("hLvsEta");
860         if (!tmp2) {
861                 AliWarning("Couln't find hLvsEta, will not calculate mismatches...");
862                 fCalculateMismatch = kFALSE;
863                 return kFALSE;  
864         }       
865
866         // Make a deep copy of the files in the histogram.
867         fT0Fill = (TH1F*)tmp1->Clone("fT0Fill");
868         fLvsEta = (TH2F*)tmp2->Clone("fLvsEta");
869
870         // Close the external file.
871         AliInfo("Closing external file.");
872         fin->Close();
873
874         // Creating a TObjArray for LvsEta projections.
875         const Int_t nbinseta = fLvsEta->GetNbinsX();
876         fLvsEtaProjections = new TObjArray(nbinseta);
877         fLvsEtaProjections->SetOwner(kTRUE);
878
879         // Making the projections needed (excluding underflow/ overflow).
880         for (Int_t iEtaBin = 1; iEtaBin < (nbinseta + 1); iEtaBin++) {
881                 TH1F* tmp = (TH1F*)fLvsEta->ProjectionY(Form("LvsEtaProjection_%i",iEtaBin),iEtaBin,iEtaBin);
882                 tmp->SetDirectory(0);
883                 fLvsEtaProjections->AddAt(tmp,iEtaBin - 1);
884         }
885
886         return kTRUE;
887
888 }
889
890 // -----------------------------------------------------------------------
891 Double_t AliAnalysisTaskDiHadronPID::GenerateRandomHit(Double_t eta) {
892
893         //
894         // Returns a random TOF time.
895         //
896
897         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
898
899         // Default (error) value:
900         Double_t rndhittime = -1.e21;
901
902         // TOF mismatch flag is not turned on.
903         if (!fCalculateMismatch) {
904                 AliFatal("Called GenerateRandomHit() method, but flag fCalculateMismatch not set.");
905                 return rndhittime;
906         }
907
908         // TOF doesn't extend much further than 0.8.
909         if (TMath::Abs(eta) > 0.8) {
910                 if (fDebug) {AliInfo("Tried to get a random hit for a track with eta > 0.8.");}
911                 return rndhittime;
912         }
913
914         // Finding the bin of the eta.
915         TAxis* etaAxis = fLvsEta->GetXaxis();
916         Int_t etaBin = etaAxis->FindBin(eta);
917         if (etaBin == 0 || (etaBin == etaAxis->GetNbins() + 1)) {return rndhittime;}
918
919         const TH1F* lengthDistribution = (const TH1F*)fLvsEtaProjections->At(etaBin - 1);
920
921         if (!lengthDistribution) {
922                 AliFatal("length Distribution not found.");
923                 return rndhittime;
924         }
925
926         Double_t currentRndLength = lengthDistribution->GetRandom(); // in cm.
927
928         // Similar to Roberto's code.
929         Double_t currentRndTime = currentRndLength / (TMath::C() * 1.e2 / 1.e12);
930         Double_t t0fill = -1.26416e+04;
931         rndhittime = fT0Fill->GetRandom() - t0fill + currentRndTime;
932
933         return rndhittime;
934
935 }
936
937 // -----------------------------------------------------------------------
938 void AliAnalysisTaskDiHadronPID::PrintPoolManagerContents() {
939
940         //
941         // Prints out the current contents of the event pool manager.
942         //
943
944         // Determine the number of pools in the pool manager.
945         AliEventPool* poolin = fPoolMgr->GetEventPool(0,0);
946         Int_t NPoolsCentrality = 0;
947         while (poolin) {
948                 NPoolsCentrality++;
949                 poolin = fPoolMgr->GetEventPool(NPoolsCentrality,0);
950         } 
951
952         poolin = fPoolMgr->GetEventPool(0,0);
953         Int_t NPoolsVtxZ = 0;   
954         while (poolin) {
955                 NPoolsVtxZ++;
956                 poolin = fPoolMgr->GetEventPool(0,NPoolsVtxZ);
957         } 
958
959         // Loop over all Pools in the matrix of the pool manager.
960         cout<<" Pool manager contents: (Nevt,NTrack)"<<endl;
961         for (Int_t iCentrality = 0; iCentrality < NPoolsCentrality; iCentrality++) {
962                 cout<<Form("Centrality Bin: %2i --> ", iCentrality);
963
964                 for (Int_t iVtxZ = 0; iVtxZ < NPoolsVtxZ; iVtxZ++) {
965
966                         poolin = fPoolMgr->GetEventPool(iCentrality, iVtxZ);
967
968                         cout<<Form("(%2i,%4i) ",poolin->GetCurrentNEvents(), poolin->NTracksInPool());
969
970                 }
971
972                 cout<<endl;
973         }
974
975 }
976
977 // -----------------------------------------------------------------------
978 void AliAnalysisTaskDiHadronPID::Terminate(Option_t*) {;
979
980         //
981         // Called when task is done.
982         //
983
984         if (fDebug > 0) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
985
986         delete fT0Fill;
987         fT0Fill = 0x0;
988         delete fLvsEta;
989         fLvsEta = 0x0;
990         delete fLvsEtaProjections;
991         fLvsEtaProjections = 0x0;
992
993 }