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 pures 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 pure. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 // AliAnalysisTaskPIDflowQA:
21 //origin: Marek Chojnacki, Marek.Chojnacki@cern.ch
22 //modified: Mikolaj Krzewicki, Mikolaj.Krzewicki@cern.ch
26 #include "AliAnalysisTaskPIDflowQA.h"
27 #include "AliAnalysisManager.h"
28 #include "AliESDEvent.h"
30 #include "AliMCEvent.h"
35 #include "AliVEvent.h"
37 #include "AliCDBManager.h"
38 #include "AliFlowEventCuts.h"
39 #include "AliFlowTrackCuts.h"
40 #include "AliVEventHandler.h"
41 #include "AliInputEventHandler.h"
45 ClassImp( AliAnalysisTaskPIDflowQA)
47 //________________________________________________________________________
48 AliAnalysisTaskPIDflowQA:: AliAnalysisTaskPIDflowQA():
49 AliAnalysisTaskSE("AliAnalysisTaskPIDflowQA"),
55 fUseDebugFile(kFALSE),
74 fTOFbetaAfterElectronsCuts(NULL),
75 fTOFbetaAfterPionCuts(NULL),
76 fTOFbetaAfterKaonCuts(NULL),
77 fTOFbetaAfterProtonCuts(NULL),
78 fTPCsignalAfterPionCuts(NULL),
79 fTPCsignalAfterKaonCuts(NULL),
80 fTPCsignalAfterProtonCuts(NULL),
81 fTOFbetaAfterElectronsCuts1(NULL),
82 fTOFbetaAfterPionCuts1(NULL),
83 fTOFbetaAfterKaonCuts1(NULL),
84 fTOFbetaAfterProtonCuts1(NULL),
85 fTPCsignalAfterPionCuts1(NULL),
86 fTPCsignalAfterKaonCuts1(NULL),
87 fTPCsignalAfterProtonCuts1(NULL),
89 fTOFbetaPiafter(NULL),
92 fTPCsignalPiafter(NULL),
93 fTPCsignalKafter(NULL),
94 fTPCsignalPafter(NULL),
95 fTOFyieldSelEmcE(NULL),
96 fTOFyieldSelPimcE(NULL),
97 fTOFyieldSelKmcE(NULL),
98 fTOFyieldSelPmcE(NULL),
99 fTOFyieldSelEmcM(NULL),
100 fTOFyieldSelPimcM(NULL),
101 fTOFyieldSelKmcM(NULL),
102 fTOFyieldSelPmcM(NULL),
103 fTOFyieldSelEmcPi(NULL),
104 fTOFyieldSelPimcPi(NULL),
105 fTOFyieldSelKmcPi(NULL),
106 fTOFyieldSelPmcPi(NULL),
107 fTOFyieldSelEmcK(NULL),
108 fTOFyieldSelPimcK(NULL),
109 fTOFyieldSelKmcK(NULL),
110 fTOFyieldSelPmcK(NULL),
111 fTOFyieldSelEmcP(NULL),
112 fTOFyieldSelPimcP(NULL),
113 fTOFyieldSelKmcP(NULL),
114 fTOFyieldSelPmcP(NULL),
115 fTPCyieldSelEmcE(NULL),
116 fTPCyieldSelPimcE(NULL),
117 fTPCyieldSelKmcE(NULL),
118 fTPCyieldSelPmcE(NULL),
119 fTPCyieldSelEmcM(NULL),
120 fTPCyieldSelPimcM(NULL),
121 fTPCyieldSelKmcM(NULL),
122 fTPCyieldSelPmcM(NULL),
123 fTPCyieldSelEmcPi(NULL),
124 fTPCyieldSelPimcPi(NULL),
125 fTPCyieldSelKmcPi(NULL),
126 fTPCyieldSelPmcPi(NULL),
127 fTPCyieldSelEmcK(NULL),
128 fTPCyieldSelPimcK(NULL),
129 fTPCyieldSelKmcK(NULL),
130 fTPCyieldSelPmcK(NULL),
131 fTPCyieldSelEmcP(NULL),
132 fTPCyieldSelPimcP(NULL),
133 fTPCyieldSelKmcP(NULL),
134 fTPCyieldSelPmcP(NULL),
135 fTPCdedxAfterTOFpidPions(NULL),
136 fTPCdedxAfterTOFpidKaons(NULL),
137 fTPCdedxAfterTOFpidProtons(NULL),
140 fTPCvsGlobalMult(NULL),
141 fStandardGlobalCuts(NULL),
142 fStandardTPCCuts(NULL),
143 fCutsTOFElectrons(NULL),
146 fCutsTOFProtons(NULL),
147 fCutsTPCElectrons(NULL),
150 fCutsTPCProtons(NULL),
156 //________________________________________________________________________
157 AliAnalysisTaskPIDflowQA:: AliAnalysisTaskPIDflowQA(const char *name):
158 AliAnalysisTaskSE(name),
164 fUseDebugFile(kFALSE),
170 fTPCsignalPimc(NULL),
183 fTOFbetaAfterElectronsCuts(NULL),
184 fTOFbetaAfterPionCuts(NULL),
185 fTOFbetaAfterKaonCuts(NULL),
186 fTOFbetaAfterProtonCuts(NULL),
187 fTPCsignalAfterPionCuts(NULL),
188 fTPCsignalAfterKaonCuts(NULL),
189 fTPCsignalAfterProtonCuts(NULL),
190 fTOFbetaAfterElectronsCuts1(NULL),
191 fTOFbetaAfterPionCuts1(NULL),
192 fTOFbetaAfterKaonCuts1(NULL),
193 fTOFbetaAfterProtonCuts1(NULL),
194 fTPCsignalAfterPionCuts1(NULL),
195 fTPCsignalAfterKaonCuts1(NULL),
196 fTPCsignalAfterProtonCuts1(NULL),
197 fTOFbetaEafter(NULL),
198 fTOFbetaPiafter(NULL),
199 fTOFbetaKafter(NULL),
200 fTOFbetaPafter(NULL),
201 fTPCsignalPiafter(NULL),
202 fTPCsignalKafter(NULL),
203 fTPCsignalPafter(NULL),
204 fTOFyieldSelEmcE(NULL),
205 fTOFyieldSelPimcE(NULL),
206 fTOFyieldSelKmcE(NULL),
207 fTOFyieldSelPmcE(NULL),
208 fTOFyieldSelEmcM(NULL),
209 fTOFyieldSelPimcM(NULL),
210 fTOFyieldSelKmcM(NULL),
211 fTOFyieldSelPmcM(NULL),
212 fTOFyieldSelEmcPi(NULL),
213 fTOFyieldSelPimcPi(NULL),
214 fTOFyieldSelKmcPi(NULL),
215 fTOFyieldSelPmcPi(NULL),
216 fTOFyieldSelEmcK(NULL),
217 fTOFyieldSelPimcK(NULL),
218 fTOFyieldSelKmcK(NULL),
219 fTOFyieldSelPmcK(NULL),
220 fTOFyieldSelEmcP(NULL),
221 fTOFyieldSelPimcP(NULL),
222 fTOFyieldSelKmcP(NULL),
223 fTOFyieldSelPmcP(NULL),
224 fTPCyieldSelEmcE(NULL),
225 fTPCyieldSelPimcE(NULL),
226 fTPCyieldSelKmcE(NULL),
227 fTPCyieldSelPmcE(NULL),
228 fTPCyieldSelEmcM(NULL),
229 fTPCyieldSelPimcM(NULL),
230 fTPCyieldSelKmcM(NULL),
231 fTPCyieldSelPmcM(NULL),
232 fTPCyieldSelEmcPi(NULL),
233 fTPCyieldSelPimcPi(NULL),
234 fTPCyieldSelKmcPi(NULL),
235 fTPCyieldSelPmcPi(NULL),
236 fTPCyieldSelEmcK(NULL),
237 fTPCyieldSelPimcK(NULL),
238 fTPCyieldSelKmcK(NULL),
239 fTPCyieldSelPmcK(NULL),
240 fTPCyieldSelEmcP(NULL),
241 fTPCyieldSelPimcP(NULL),
242 fTPCyieldSelKmcP(NULL),
243 fTPCyieldSelPmcP(NULL),
244 fTPCdedxAfterTOFpidPions(NULL),
245 fTPCdedxAfterTOFpidKaons(NULL),
246 fTPCdedxAfterTOFpidProtons(NULL),
249 fTPCvsGlobalMult(NULL),
250 fStandardGlobalCuts(NULL),
251 fStandardTPCCuts(NULL),
252 fCutsTOFElectrons(NULL),
255 fCutsTOFProtons(NULL),
256 fCutsTPCElectrons(NULL),
259 fCutsTPCProtons(NULL),
263 fESDpid=new AliESDpid();
266 fESDpid->GetTPCResponse().SetBetheBlochParameters(0.0283086,
272 //fESDpid->GetTPCResponse().SetBetheBlochParameters(1.28949/50.,
274 // TMath::Exp(-3.21763e+01),
278 DefineOutput(1, TList::Class());
281 //________________________________________________________________________
282 void AliAnalysisTaskPIDflowQA::UserCreateOutputObjects()
284 //UserCreateOutputObject
285 if (fOutputList) fOutputList->Delete();
287 fOutputList=new TList();
288 fOutputList->SetOwner(kTRUE);
292 const Int_t npredec=50;
293 Double_t tabx[ndec*npredec+1];
294 for (Int_t i=0; i<ndec; i++)
296 for (Int_t j=0; j<npredec; j++)
298 tabx[npredec*i+j]=TMath::Power(10,((Double_t)i)+((Double_t)startvalue)+((Double_t)j)/((Double_t)npredec));
301 tabx[ndec*npredec]=TMath::Power(10,ndec+startvalue);
304 Double_t binsPtDummy[kPtBins+1];
306 for(int i=1; i<=kPtBins+1; i++)
308 if(binsPtDummy[i-1]+0.05<1.01)
309 binsPtDummy[i]=binsPtDummy[i-1]+0.05;
311 binsPtDummy[i]=binsPtDummy[i-1]+0.1;
315 Double_t binsPDummy[kPBins+1];
317 for(int i=1; i<=kPBins+1; i++)
319 if(binsPDummy[i-1]+0.05<1.01)
320 binsPDummy[i]=binsPDummy[i-1]+0.05;
322 binsPDummy[i]=binsPDummy[i-1]+0.1;
325 fTPCsignal=new TH2F("fTPCsignal",";p [GeV/c];dEdx",kPBins,binsPDummy,500,0,500);
326 fOutputList->Add(fTPCsignal);
327 fTPCdedxAfterTOFpidPions=new TH2F("fTPCsignalAfterTOFpions",";p [GeV/c];dEdx",kPBins,binsPDummy,500,0,500);
328 fOutputList->Add(fTPCdedxAfterTOFpidPions);
329 fTPCdedxAfterTOFpidKaons=new TH2F("fTPCsignalAfterTOFkaons",";p [GeV/c];dEdx",kPBins,binsPDummy,500,0,500);
330 fOutputList->Add(fTPCdedxAfterTOFpidKaons);
331 fTPCdedxAfterTOFpidProtons=new TH2F("fTPCsignalAfterTOFprotons",";p [GeV/c];dEdx",kPBins,binsPDummy,500,0,500);
332 fOutputList->Add(fTPCdedxAfterTOFpidProtons);
333 fTPCsignalAfterPionCuts=new TH2F("fTPCsignalAfterPionCuts",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
334 fTPCsignalAfterKaonCuts=new TH2F("fTPCsignalAfterKaonCuts",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
335 fTPCsignalAfterProtonCuts=new TH2F("fTPCsignalAfterProtonCuts",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
336 fOutputList->Add(fTPCsignalAfterPionCuts);
337 fOutputList->Add(fTPCsignalAfterKaonCuts);
338 fOutputList->Add(fTPCsignalAfterProtonCuts);
339 fTPCsignalAfterPionCuts1=new TH2F("fTPCsignalAfterPionCuts1",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
340 fTPCsignalAfterKaonCuts1=new TH2F("fTPCsignalAfterKaonCuts1",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
341 fTPCsignalAfterProtonCuts1=new TH2F("fTPCsignalAfterProtonCuts1",";p [GeV/c];dE/dx",kPBins,binsPDummy,500, 0., 500.);//
342 fOutputList->Add(fTPCsignalAfterPionCuts1);
343 fOutputList->Add(fTPCsignalAfterKaonCuts1);
344 fOutputList->Add(fTPCsignalAfterProtonCuts1);
346 fTPCsignalPi=new TH2F("fTPCsignalPi",";p [GeV/c];signal",kPBins,binsPDummy,300,-2,2);//TPC PID signal as function of p for pi+
347 fTPCsignalK=new TH2F("fTPCsignalK",";p [GeV/c];signal",kPBins,binsPDummy,300,-2,2);//TPC PID signal as function of p for K+
348 fTPCsignalP=new TH2F("fTPCsignalP",";p [GeV/c];signal",kPBins,binsPDummy,300,-2,2);//TPC PID signal as function of p for p
349 fOutputList->Add(fTPCsignalPi);
350 fOutputList->Add(fTPCsignalK);
351 fOutputList->Add(fTPCsignalP);
352 fTPCsignalPiafter=new TH2F("fTPCsignalPiafter",";p [GeV/c];(dE/dx-dE/dx_{#pi})/dE/dx_{#pi}",kPBins,binsPDummy,300, -2., 2.);//
353 fTPCsignalKafter=new TH2F("fTPCsignalKafter",";p [GeV/c];(dE/dx-dE/dx_{K})/dE/dx_{K}",kPBins,binsPDummy,300, -2., 2.);//
354 fTPCsignalPafter=new TH2F("fTPCsignalPafter",";p [GeV/c];(dE/dx-dE/dx_{p})/dE/dx_{p}",kPBins,binsPDummy,300, -2., 2.);//
355 fOutputList->Add(fTPCsignalPiafter);
356 fOutputList->Add(fTPCsignalKafter);
357 fOutputList->Add(fTPCsignalPafter);
361 fTPCsignalPimc=new TH2F("fTPCsignalPionsMC",";p [GeV/c];signal",kPBins,binsPDummy,600,-2,2);//TPC PID signal as function of pt for pi+
362 fTPCsignalKmc=new TH2F("fTPCsignalKaonsMC",";p [GeV/c];signal",kPBins,binsPDummy,600,-2,2);//TPC PID signal as function of pt for K+
363 fTPCsignalPmc=new TH2F("fTPCsignalProtonsMC",";p [GeV/c];signal",kPBins,binsPDummy,600,-2,2);//TPC PID signal as function of pt for p
364 fOutputList->Add(fTPCsignalPimc);
365 fOutputList->Add(fTPCsignalKmc);
366 fOutputList->Add(fTPCsignalPmc);
369 fTOFbeta=new TH2F("fTOFbeta",";p[GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
370 fOutputList->Add(fTOFbeta);
371 fTOFbetaE=new TH2F("fTOFbetaE",";p [GeV/c];#beta-#beta_{#pi}",kPBins,binsPDummy,500, -0.25, 0.25);//
372 fTOFbetaPi=new TH2F("fTOFbetaPi",";p [GeV/c];#beta-#beta_{#pi}",kPBins,binsPDummy,500, -0.25, 0.25);//
373 fTOFbetaK=new TH2F("fTOFbetaK",";p [GeV/c];#beta-#beta_{K}",kPBins,binsPDummy,500, -0.25, 0.25);//
374 fTOFbetaP=new TH2F("fTOFbetaP",";p [GeV/c];#beta-#beta_{p}",kPBins,binsPDummy,500, -0.25, 0.25);//
375 fOutputList->Add(fTOFbetaE);
376 fOutputList->Add(fTOFbetaPi);
377 fOutputList->Add(fTOFbetaK);
378 fOutputList->Add(fTOFbetaP);
380 fTOFinvbeta=new TH2F("fTOFinvbeta",";p[GeV/c];1/#beta",kPBins,binsPDummy,1000, 0.90, 2.5);//
381 fOutputList->Add(fTOFinvbeta);
382 fTOFinvbetaE=new TH2F("fTOFinvbetaE",";p [GeV/c];1/#beta-1/#beta_{#pi}",kPBins,binsPDummy,600, -0.3, 0.3);//
383 fTOFinvbetaPi=new TH2F("fTOFinvbetaPi",";p [GeV/c];1/#beta-1/#beta_{#pi}",kPBins,binsPDummy,600, -0.3, 0.3);//
384 fTOFinvbetaK=new TH2F("fTOFinvbetaK",";p [GeV/c];1/#beta-1/#beta_{K}",kPBins,binsPDummy,600, -0.3, 0.3);//
385 fTOFinvbetaP=new TH2F("fTOFinvbetaP",";p [GeV/c];1/#beta-1/#beta_{p}",kPBins,binsPDummy,600, -0.3, 0.3);//
386 fOutputList->Add(fTOFinvbetaE);
387 fOutputList->Add(fTOFinvbetaPi);
388 fOutputList->Add(fTOFinvbetaK);
389 fOutputList->Add(fTOFinvbetaP);
391 fTOFbetaAfterElectronsCuts=new TH2F("fTOFbetaAfterElectronsCuts",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
392 fTOFbetaAfterPionCuts=new TH2F("fTOFbetaAfterPionCuts",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
393 fTOFbetaAfterKaonCuts=new TH2F("fTOFbetaAfterKaonCuts",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
394 fTOFbetaAfterProtonCuts=new TH2F("fTOFbetaAfterProtonCuts",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
395 fOutputList->Add(fTOFbetaAfterElectronsCuts);
396 fOutputList->Add(fTOFbetaAfterPionCuts);
397 fOutputList->Add(fTOFbetaAfterKaonCuts);
398 fOutputList->Add(fTOFbetaAfterProtonCuts);
399 fTOFbetaAfterElectronsCuts1=new TH2F("fTOFbetaAfterElectronsCuts1",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
400 fTOFbetaAfterPionCuts1=new TH2F("fTOFbetaAfterPionCuts1",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
401 fTOFbetaAfterKaonCuts1=new TH2F("fTOFbetaAfterKaonCuts1",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
402 fTOFbetaAfterProtonCuts1=new TH2F("fTOFbetaAfterProtonCuts1",";p [GeV/c];#beta",kPBins,binsPDummy,1000, 0.4, 1.1);//
403 fOutputList->Add(fTOFbetaAfterElectronsCuts1);
404 fOutputList->Add(fTOFbetaAfterPionCuts1);
405 fOutputList->Add(fTOFbetaAfterKaonCuts1);
406 fOutputList->Add(fTOFbetaAfterProtonCuts1);
408 fTOFbetaEafter=new TH2F("fTOFbetaEafter",";p [GeV/c];#beta-#beta_{#pi}",kPBins,binsPDummy,500, -0.25, 0.25 );//
409 fTOFbetaPiafter=new TH2F("fTOFbetaPiafter",";p [GeV/c];#beta-#beta_{#pi}",kPBins,binsPDummy,500, -0.25, 0.25 );//
410 fTOFbetaKafter=new TH2F("fTOFbetaKafter",";p [GeV/c];#beta-#beta_{K}",kPBins,binsPDummy,500, -0.25, 0.25 );//
411 fTOFbetaPafter=new TH2F("fTOFbetaPafter",";p [GeV/c];#beta-#beta_{p}",kPBins,binsPDummy,500, -0.25, 0.25 );//
412 fOutputList->Add(fTOFbetaEafter);
413 fOutputList->Add(fTOFbetaPiafter);
414 fOutputList->Add(fTOFbetaKafter);
415 fOutputList->Add(fTOFbetaPafter);
419 fTOFyieldSelEmcE = new TH1F("fTOFyieldSelEmcE",";p [Gev/c];",kPBins,binsPDummy);
420 fTOFyieldSelPimcE = new TH1F("fTOFyieldSelPimcE",";p [Gev/c];",kPBins,binsPDummy);
421 fTOFyieldSelKmcE = new TH1F("fTOFyieldSelKmcE",";p [Gev/c];",kPBins,binsPDummy);
422 fTOFyieldSelPmcE = new TH1F("fTOFyieldSelPmcE",";p [Gev/c];",kPBins,binsPDummy);
423 fOutputList->Add(fTOFyieldSelEmcE);
424 fOutputList->Add(fTOFyieldSelPimcE);
425 fOutputList->Add(fTOFyieldSelKmcE);
426 fOutputList->Add(fTOFyieldSelPmcE);
427 fTOFyieldSelEmcM = new TH1F("fTOFyieldSelEmcM",";p [Gev/c];",kPBins,binsPDummy);
428 fTOFyieldSelPimcM = new TH1F("fTOFyieldSelPimcM",";p [Gev/c];",kPBins,binsPDummy);
429 fTOFyieldSelKmcM = new TH1F("fTOFyieldSelKmcM",";p [Gev/c];",kPBins,binsPDummy);
430 fTOFyieldSelPmcM = new TH1F("fTOFyieldSelPmcM",";p [Gev/c];",kPBins,binsPDummy);
431 fOutputList->Add(fTOFyieldSelEmcM);
432 fOutputList->Add(fTOFyieldSelPimcM);
433 fOutputList->Add(fTOFyieldSelKmcM);
434 fOutputList->Add(fTOFyieldSelPmcM);
435 fTOFyieldSelEmcPi = new TH1F("fTOFyieldSelEmcPi",";p [Gev/c];",kPBins,binsPDummy);
436 fTOFyieldSelPimcPi = new TH1F("fTOFyieldSelPimcPi",";p [Gev/c];",kPBins,binsPDummy);
437 fTOFyieldSelKmcPi = new TH1F("fTOFyieldSelKmcPi",";p [Gev/c];",kPBins,binsPDummy);
438 fTOFyieldSelPmcPi = new TH1F("fTOFyieldSelPmcPi",";p [Gev/c];",kPBins,binsPDummy);
439 fOutputList->Add(fTOFyieldSelEmcPi);
440 fOutputList->Add(fTOFyieldSelPimcPi);
441 fOutputList->Add(fTOFyieldSelKmcPi);
442 fOutputList->Add(fTOFyieldSelPmcPi);
443 fTOFyieldSelEmcK = new TH1F("fTOFyieldSelEmcK",";p [Gev/c];",kPBins,binsPDummy);
444 fTOFyieldSelPimcK = new TH1F("fTOFyieldSelPimcK",";p [Gev/c];",kPBins,binsPDummy);
445 fTOFyieldSelKmcK = new TH1F("fTOFyieldSelKmcK",";p [Gev/c];",kPBins,binsPDummy);
446 fTOFyieldSelPmcK = new TH1F("fTOFyieldSelPmcK",";p [Gev/c];",kPBins,binsPDummy);
447 fOutputList->Add(fTOFyieldSelEmcK);
448 fOutputList->Add(fTOFyieldSelPimcK);
449 fOutputList->Add(fTOFyieldSelKmcK);
450 fOutputList->Add(fTOFyieldSelPmcK);
451 fTOFyieldSelEmcP = new TH1F("fTOFyieldSelEmcP",";p [Gev/c];",kPBins,binsPDummy);
452 fTOFyieldSelPimcP = new TH1F("fTOFyieldSelPimcP",";p [Gev/c];",kPBins,binsPDummy);
453 fTOFyieldSelKmcP = new TH1F("fTOFyieldSelKmcP",";p [Gev/c];",kPBins,binsPDummy);
454 fTOFyieldSelPmcP = new TH1F("fTOFyieldSelPmcP",";p [Gev/c];",kPBins,binsPDummy);
455 fOutputList->Add(fTOFyieldSelEmcP);
456 fOutputList->Add(fTOFyieldSelPimcP);
457 fOutputList->Add(fTOFyieldSelKmcP);
458 fOutputList->Add(fTOFyieldSelPmcP);
459 fTPCyieldSelEmcE = new TH1F("fTPCyieldSelEmcE",";p [Gev/c];",kPBins,binsPDummy);
460 fTPCyieldSelPimcE = new TH1F("fTPCyieldSelPimcE",";p [Gev/c];",kPBins,binsPDummy);
461 fTPCyieldSelKmcE = new TH1F("fTPCyieldSelKmcE",";p [Gev/c];",kPBins,binsPDummy);
462 fTPCyieldSelPmcE = new TH1F("fTPCyieldSelPmcE",";p [Gev/c];",kPBins,binsPDummy);
463 fOutputList->Add(fTPCyieldSelEmcE);
464 fOutputList->Add(fTPCyieldSelPimcE);
465 fOutputList->Add(fTPCyieldSelKmcE);
466 fOutputList->Add(fTPCyieldSelPmcE);
467 fTPCyieldSelEmcM = new TH1F("fTPCyieldSelEmcM",";p [Gev/c];",kPBins,binsPDummy);
468 fTPCyieldSelPimcM = new TH1F("fTPCyieldSelPimcM",";p [Gev/c];",kPBins,binsPDummy);
469 fTPCyieldSelKmcM = new TH1F("fTPCyieldSelKmcM",";p [Gev/c];",kPBins,binsPDummy);
470 fTPCyieldSelPmcM = new TH1F("fTPCyieldSelPmcM",";p [Gev/c];",kPBins,binsPDummy);
471 fOutputList->Add(fTPCyieldSelEmcM);
472 fOutputList->Add(fTPCyieldSelPimcM);
473 fOutputList->Add(fTPCyieldSelKmcM);
474 fOutputList->Add(fTPCyieldSelPmcM);
475 fTPCyieldSelEmcPi = new TH1F("fTPCyieldSelEmcPi",";p [Gev/c];",kPBins,binsPDummy);
476 fTPCyieldSelPimcPi = new TH1F("fTPCyieldSelPimcPi",";p [Gev/c];",kPBins,binsPDummy);
477 fTPCyieldSelKmcPi = new TH1F("fTPCyieldSelKmcPi",";p [Gev/c];",kPBins,binsPDummy);
478 fTPCyieldSelPmcPi = new TH1F("fTPCyieldSelPmcPi",";p [Gev/c];",kPBins,binsPDummy);
479 fOutputList->Add(fTPCyieldSelEmcPi);
480 fOutputList->Add(fTPCyieldSelPimcPi);
481 fOutputList->Add(fTPCyieldSelKmcPi);
482 fOutputList->Add(fTPCyieldSelPmcPi);
483 fTPCyieldSelEmcK = new TH1F("fTPCyieldSelEmcK",";p [Gev/c];",kPBins,binsPDummy);
484 fTPCyieldSelPimcK = new TH1F("fTPCyieldSelPimcK",";p [Gev/c];",kPBins,binsPDummy);
485 fTPCyieldSelKmcK = new TH1F("fTPCyieldSelKmcK",";p [Gev/c];",kPBins,binsPDummy);
486 fTPCyieldSelPmcK = new TH1F("fTPCyieldSelPmcK",";p [Gev/c];",kPBins,binsPDummy);
487 fOutputList->Add(fTPCyieldSelEmcK);
488 fOutputList->Add(fTPCyieldSelPimcK);
489 fOutputList->Add(fTPCyieldSelKmcK);
490 fOutputList->Add(fTPCyieldSelPmcK);
491 fTPCyieldSelEmcP = new TH1F("fTPCyieldSelEmcP",";p [Gev/c];",kPBins,binsPDummy);
492 fTPCyieldSelPimcP = new TH1F("fTPCyieldSelPimcP",";p [Gev/c];",kPBins,binsPDummy);
493 fTPCyieldSelKmcP = new TH1F("fTPCyieldSelKmcP",";p [Gev/c];",kPBins,binsPDummy);
494 fTPCyieldSelPmcP = new TH1F("fTPCyieldSelPmcP",";p [Gev/c];",kPBins,binsPDummy);
495 fOutputList->Add(fTPCyieldSelEmcP);
496 fOutputList->Add(fTPCyieldSelPimcP);
497 fOutputList->Add(fTPCyieldSelKmcP);
498 fOutputList->Add(fTPCyieldSelPmcP);
501 fPvsPt=new TH2F("fPvsPt","p vs p_{t};p [GeV/c];p_{t} [GeV/c]",kPBins,binsPDummy,kPtBins,binsPtDummy);
502 fOutputList->Add(fPvsPt);
504 fMeanPvsP = new TProfile("fMeanPvsP","Mean P vs P;p [Gev/c];<p> [GeV/c]",kPBins,binsPDummy);
505 fOutputList->Add(fMeanPvsP);
507 fTPCvsGlobalMult = new TH2F("fTPCvsGlobalMult","TPC only vs Global track multiplicity;global;TPC only",500,0,2500,500,0,3500);
508 fOutputList->Add(fTPCvsGlobalMult);
510 fStandardGlobalCuts = AliFlowTrackCuts::GetStandardGlobalTrackCuts2010();
511 fStandardTPCCuts = AliFlowTrackCuts::GetStandardTPCStandaloneTrackCuts2010();
513 fCutsTOFElectrons = new AliFlowTrackCuts("tof electron cuts");
514 fCutsTOFElectrons->SetPID(AliPID::kElectron, AliFlowTrackCuts::kTOFbeta);
515 fCutsTOFElectrons->SetQA();
516 fCutsTOFPions = new AliFlowTrackCuts("tof pion cuts");
517 fCutsTOFPions->SetPID(AliPID::kPion, AliFlowTrackCuts::kTOFbeta);
518 fCutsTOFPions->SetQA();
519 fCutsTOFKaons = new AliFlowTrackCuts("tof kaon cuts");
520 fCutsTOFKaons->SetPID(AliPID::kKaon, AliFlowTrackCuts::kTOFbeta);
521 fCutsTOFKaons->SetQA();
522 fCutsTOFProtons = new AliFlowTrackCuts("tof proton cuts");
523 fCutsTOFProtons->SetPID(AliPID::kProton, AliFlowTrackCuts::kTOFbeta);
524 fCutsTOFProtons->SetQA();
525 fCutsTPCElectrons = new AliFlowTrackCuts("tpc electron cuts");
526 fCutsTPCElectrons->SetPID(AliPID::kElectron, AliFlowTrackCuts::kTPCdedx);
527 fCutsTPCElectrons->SetQA();
528 fCutsTPCPions = new AliFlowTrackCuts("tpc pion cuts");
529 fCutsTPCPions->SetPID(AliPID::kPion, AliFlowTrackCuts::kTPCdedx);
530 fCutsTPCPions->SetQA();
531 fCutsTPCKaons = new AliFlowTrackCuts("tpc kaon cuts");
532 fCutsTPCKaons->SetPID(AliPID::kKaon, AliFlowTrackCuts::kTPCdedx);
533 fCutsTPCKaons->SetQA();
534 fCutsTPCProtons = new AliFlowTrackCuts("tpc proton cuts");
535 fCutsTPCProtons->SetPID(AliPID::kProton, AliFlowTrackCuts::kTPCdedx);
536 fCutsTPCProtons->SetQA();
538 //fOutputList->Add(fESDpid);
541 h = static_cast<TH1*>(fCutsTOFPions->GetQA()->At(0));
542 h->SetName(Form("pion direct %s",h->GetName()));
543 h = static_cast<TH1*>(fCutsTOFKaons->GetQA()->At(0));
544 h->SetName(Form("kaon direct %s",h->GetName()));
545 h = static_cast<TH1*>(fCutsTOFProtons->GetQA()->At(0));
546 h->SetName(Form("proton direct %s",h->GetName()));
547 fOutputList->Add(fCutsTOFPions->GetQA()->At(0));
548 fOutputList->Add(fCutsTOFKaons->GetQA()->At(0));
549 fOutputList->Add(fCutsTOFProtons->GetQA()->At(0));
551 h = static_cast<TH1*>(fCutsTOFPions->GetQA()->At(1));
552 h->SetName(Form("pion direct %s",h->GetName()));
553 h = static_cast<TH1*>(fCutsTOFKaons->GetQA()->At(1));
554 h->SetName(Form("kaon direct %s",h->GetName()));
555 h = static_cast<TH1*>(fCutsTOFProtons->GetQA()->At(1));
556 h->SetName(Form("proton direct %s",h->GetName()));
557 fOutputList->Add(fCutsTOFPions->GetQA()->At(1));
558 fOutputList->Add(fCutsTOFKaons->GetQA()->At(1));
559 fOutputList->Add(fCutsTOFProtons->GetQA()->At(1));
561 h = static_cast<TH1*>(fCutsTOFPions->GetQA()->At(2));
562 h->SetName(Form("pion direct %s",h->GetName()));
563 h = static_cast<TH1*>(fCutsTOFKaons->GetQA()->At(2));
564 h->SetName(Form("kaon direct %s",h->GetName()));
565 h = static_cast<TH1*>(fCutsTOFProtons->GetQA()->At(2));
566 h->SetName(Form("proton direct %s",h->GetName()));
567 fOutputList->Add(fCutsTOFPions->GetQA()->At(2));
568 fOutputList->Add(fCutsTOFKaons->GetQA()->At(2));
569 fOutputList->Add(fCutsTOFProtons->GetQA()->At(2));
571 h = static_cast<TH1*>(fCutsTOFPions->GetQA()->At(3));
572 h->SetName(Form("pion direct %s",h->GetName()));
573 h = static_cast<TH1*>(fCutsTOFKaons->GetQA()->At(3));
574 h->SetName(Form("kaon direct %s",h->GetName()));
575 h = static_cast<TH1*>(fCutsTOFProtons->GetQA()->At(3));
576 h->SetName(Form("proton direct %s",h->GetName()));
577 fOutputList->Add(fCutsTOFPions->GetQA()->At(3));
578 fOutputList->Add(fCutsTOFKaons->GetQA()->At(3));
579 fOutputList->Add(fCutsTOFProtons->GetQA()->At(3));
581 if (fUseDebugFile) fFile = fopen("debug.txt","w");
583 PostData(1, fOutputList);
586 //________________________________________________________________________
587 void AliAnalysisTaskPIDflowQA::UserExec(Option_t *)
589 fESD = dynamic_cast<AliESDEvent*> (InputEvent());
592 //do the calibration bit
593 fESDpid->SetTOFResponse(fESD,AliESDpid::kTOF_T0); // to use T0-TOF
594 fESDpid->MakePID(fESD,kFALSE);
596 if(!fCuts || !fEventCuts)
598 Printf("No CUTS Defined.........\n");
599 PostData(1, fOutputList);
603 if (!(fEventCuts->IsSelected(fESD)))
609 AliMCEvent* mcEvent=NULL;
612 mcEvent = (AliMCEvent*) MCEvent();
613 Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
614 stack = mcEvent->Stack();
617 Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
618 Int_t nTracks=fESD->GetNumberOfTracks();
620 AliESDtrack *trackESD=0;
622 for(int tr1=0; tr1<nTracks; tr1++)
624 trackESD=fESD->GetTrack(tr1);
625 if (!trackESD) continue;
627 Double_t p=trackESD->GetP();
628 Double_t pt=trackESD->Pt();
630 if(!(fCuts->IsSelected(trackESD))) continue;
633 if(fMC) label=trackESD->GetLabel();
638 TParticle* particle2 = stack->Particle(TMath::Abs(label));
639 pdgcode=particle2->GetPdgCode();
643 fMeanPvsP->Fill(p,p);
645 pidTPC(trackESD,pdgcode);
646 pidTOF(trackESD,pdgcode);
649 //check the correlation between the global and TPConly number of tracks
650 fStandardGlobalCuts->SetEvent(fESD);
651 fStandardTPCCuts->SetEvent(fESD);
652 Int_t multGlobal = fStandardGlobalCuts->Count();
653 Int_t multTPC = fStandardTPCCuts->Count();
654 fTPCvsGlobalMult->Fill(multGlobal,multTPC);
658 const AliESDVertex* pvtx = fESD->GetPrimaryVertex();
659 const AliESDVertex* tpcvtx = fESD->GetPrimaryVertexTPC();
660 const AliESDVertex* spdvtx = fESD->GetPrimaryVertexSPD();
661 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
662 AliVEventHandler* handler = mgr->GetInputEventHandler();
663 TTree* tree = handler->GetTree();
664 TFile* file = tree->GetCurrentFile();
665 if (multTPC>(23+1.216*multGlobal) || multTPC<(-20+1.087*multGlobal))
667 fprintf(fFile, "%i %i %s %i\n",multTPC,multGlobal,file->GetName(),fESD->GetEventNumberInFile());
668 fprintf(fFile, " primary vertex: x: %.2f, y: %.2f, z: %.2f, n: %i\n", pvtx->GetX(), pvtx->GetY(), pvtx->GetZ(), pvtx->GetNContributors());
669 fprintf(fFile, " SPD vertex: x: %.2f, y: %.2f, z: %.2f, n: %i\n", spdvtx->GetX(), spdvtx->GetY(), spdvtx->GetZ(), spdvtx->GetNContributors());
670 fprintf(fFile, " TPC vertex: x: %.2f, y: %.2f, z: %.2f, n: %i\n", tpcvtx->GetX(), tpcvtx->GetY(), tpcvtx->GetZ(), tpcvtx->GetNContributors());
675 //________________________________________________________________________
676 void AliAnalysisTaskPIDflowQA::Terminate(Option_t *)
684 Printf("AliAnalysisTaskPIDflowQA: end of Terminate");
688 //________________________________________________________________________
689 void AliAnalysisTaskPIDflowQA::pidTPC(AliESDtrack* t, Int_t pdgcode)
692 const AliExternalTrackParam* innerParam = t->GetInnerParam();
693 if (!innerParam) return;
694 Double_t pinTPCglobal=innerParam->GetP();
695 Double_t tpcSignal =t ->GetTPCsignal();
696 Float_t p=innerParam->P();
697 Float_t pt=innerParam->Pt();
698 Float_t sigPion = fESDpid->GetTPCResponse().GetExpectedSignal(pinTPCglobal, AliPID::kPion);
699 Float_t sigKaon = fESDpid->GetTPCResponse().GetExpectedSignal(pinTPCglobal, AliPID::kKaon);
700 Float_t sigProton = fESDpid->GetTPCResponse().GetExpectedSignal(pinTPCglobal, AliPID::kProton);
701 if(!(sigPion>0.0&&sigKaon>0.0&&sigProton>0.0))
704 fTPCsignal->Fill(pinTPCglobal,tpcSignal);
706 fTPCsignalPi->Fill(p,(tpcSignal-sigPion)/sigPion);
707 fTPCsignalK->Fill(p,(tpcSignal-sigKaon)/sigKaon);
708 fTPCsignalP->Fill(p,(tpcSignal-sigProton)/sigProton);
710 if (fCutsTPCPions->IsSelected(t))
712 fTPCsignalAfterPionCuts->Fill(p,tpcSignal);
713 fTPCsignalPiafter->Fill(p,(tpcSignal-sigPion)/sigPion);
716 if (TMath::Abs(pdgcode)==11)
718 fTPCyieldSelPimcE->Fill(pt);
720 if (TMath::Abs(pdgcode)==13)
722 fTPCyieldSelPimcM->Fill(pt);
724 if (TMath::Abs(pdgcode)==211)
726 fTPCyieldSelPimcPi->Fill(pt);
728 else if(TMath::Abs(pdgcode)==321)
730 fTPCyieldSelPimcK->Fill(pt);
732 else if (TMath::Abs(pdgcode)==2212)
734 fTPCyieldSelPimcP->Fill(pt);
738 if (fCutsTPCKaons->IsSelected(t))
740 fTPCsignalAfterKaonCuts->Fill(p,tpcSignal);
741 fTPCsignalKafter->Fill(p,(tpcSignal-sigKaon)/sigKaon);
744 if (TMath::Abs(pdgcode)==11)
746 fTPCyieldSelKmcE->Fill(pt);
748 if (TMath::Abs(pdgcode)==13)
750 fTPCyieldSelKmcM->Fill(pt);
752 if (TMath::Abs(pdgcode)==211)
754 fTPCyieldSelKmcPi->Fill(pt);
756 else if(TMath::Abs(pdgcode)==321)
758 fTPCyieldSelKmcK->Fill(pt);
760 else if (TMath::Abs(pdgcode)==2212)
762 fTPCyieldSelKmcP->Fill(pt);
766 if (fCutsTPCProtons->IsSelected(t))
768 fTPCsignalAfterProtonCuts->Fill(p,tpcSignal);
769 fTPCsignalPafter->Fill(p,(tpcSignal-sigProton)/sigProton);
772 if (TMath::Abs(pdgcode)==11)
774 fTPCyieldSelPmcE->Fill(pt);
776 if (TMath::Abs(pdgcode)==13)
778 fTPCyieldSelPmcM->Fill(pt);
780 if (TMath::Abs(pdgcode)==211)
782 fTPCyieldSelPmcPi->Fill(pt);
784 else if(TMath::Abs(pdgcode)==321)
786 fTPCyieldSelPmcK->Fill(pt);
788 else if (TMath::Abs(pdgcode)==2212)
790 fTPCyieldSelPmcP->Fill(pt);
797 if(TMath::Abs(pdgcode)==211)
798 fTPCsignalPimc->Fill(p,(tpcSignal-sigPion)/sigPion);
799 else if(TMath::Abs(pdgcode)==321)
800 fTPCsignalKmc->Fill(p,(tpcSignal-sigKaon)/sigKaon);
801 else if (TMath::Abs(pdgcode)==2212)
802 fTPCsignalPmc->Fill(p,(tpcSignal-sigProton)/sigProton);
806 //______________________________________________________________________________
807 void AliAnalysisTaskPIDflowQA::pidTOF(AliESDtrack* track, Int_t pdgcode)
810 Bool_t goodtrack = (track) &&
811 (track->GetStatus() & AliESDtrack::kTOFpid) &&
812 (track->GetTOFsignal() > 12000) &&
813 (track->GetTOFsignal() < 100000) &&
814 (track->GetIntegratedLength() > 365) &&
815 !(track->GetStatus() & AliESDtrack::kTOFmismatch);
817 if (!goodtrack) return;
819 const AliExternalTrackParam* innerParam = track->GetInnerParam();
820 if (!innerParam) return;
821 Double_t pinTPCglobal=innerParam->GetP();
822 Double_t tpcSignal =track->GetTPCsignal();
823 const Float_t c = 2.99792457999999984e-02;
824 Float_t p = track->GetP();
825 Float_t pt = track->Pt();
826 Float_t l = track->GetIntegratedLength();
827 Float_t trackT0 = fESDpid->GetTOFResponse().GetStartTime(p);
828 Float_t timeTOF = track->GetTOFsignal()- trackT0;
829 Double_t integratedTimes[5] = {-1.0,-1.0,-1.0,-1.0,-1.0};
830 track->GetIntegratedTimes(integratedTimes);
833 Float_t beta = l/timeTOF/c;
834 Float_t betaHypothesis[5] = {0.0,0.0,0.0,0.0,0.0};
835 Float_t betadiff[5] = {0.0,0.0,0.0,0.0,0.0};
836 for (Int_t i=0;i<5;i++)
838 betaHypothesis[i] = l/integratedTimes[i]/c;
839 betadiff[i] = beta-betaHypothesis[i];
843 Float_t invbeta = 1/beta;
844 Float_t invbetaHypothesis[5] = {0.0,0.0,0.0,0.0,0.0};
845 Float_t invbetadiff[5] = {0.0,0.0,0.0,0.0,0.0};
846 for (Int_t i=0;i<5;i++)
848 invbetaHypothesis[i] = 1/betaHypothesis[i];
849 invbetadiff[i] = invbeta-invbetaHypothesis[i];
852 /////////////simple cuts
853 Bool_t isPion = ( (betadiff[2]<0.015) && (betadiff[2]>-0.015) &&
854 (betadiff[3]>0.025) &&
855 (betadiff[4]>0.03) );
857 Bool_t isKaon = ( (betadiff[3]<0.015) && (betadiff[3]>-0.015) &&
858 (betadiff[2]<-0.03) &&
859 (betadiff[4]>0.03) );
861 Bool_t isProton = ( (betadiff[4]<0.015) && (betadiff[4]>-0.015) &&
862 (betadiff[3]<-0.025) &&
863 (betadiff[2]<-0.025) );
865 if (isPion) fTOFbetaAfterPionCuts1->Fill(p,beta);
866 if (isKaon) fTOFbetaAfterKaonCuts1->Fill(p,beta);
867 if (isProton) fTOFbetaAfterProtonCuts1->Fill(p,beta);
870 fTOFbeta->Fill(p,beta);
871 fTOFbetaE->Fill(p,betadiff[0]);
872 fTOFbetaPi->Fill(p,betadiff[2]);
873 fTOFbetaK->Fill(p,betadiff[3]);
874 fTOFbetaP->Fill(p,betadiff[4]);
876 fTOFinvbeta->Fill(p,invbeta);
877 fTOFinvbetaE->Fill(p,invbetadiff[0]);
878 fTOFinvbetaPi->Fill(p,invbetadiff[2]);
879 fTOFinvbetaK->Fill(p,invbetadiff[3]);
880 fTOFinvbetaP->Fill(p,invbetadiff[4]);
882 if (fCutsTOFElectrons->PassesTOFbetaCut(track))
884 fTOFbetaAfterElectronsCuts->Fill(p,beta);
885 fTOFbetaEafter->Fill(p,beta-betaHypothesis[0]);
888 if (TMath::Abs(pdgcode)==11)
890 fTOFyieldSelEmcE->Fill(pt);
892 if (TMath::Abs(pdgcode)==13)
894 fTOFyieldSelEmcM->Fill(pt);
896 if (TMath::Abs(pdgcode)==211)
898 fTOFyieldSelEmcPi->Fill(pt);
900 else if(TMath::Abs(pdgcode)==321)
902 fTOFyieldSelEmcK->Fill(pt);
904 else if (TMath::Abs(pdgcode)==2212)
906 fTOFyieldSelEmcP->Fill(pt);
910 if (fCutsTOFPions->PassesTOFbetaCut(track))
912 fTPCdedxAfterTOFpidPions->Fill(pinTPCglobal,tpcSignal);
913 fTOFbetaAfterPionCuts->Fill(p,beta);
914 fTOFbetaPiafter->Fill(p,beta-betaHypothesis[2]);
917 if (TMath::Abs(pdgcode)==11)
919 fTOFyieldSelPimcE->Fill(pt);
921 if (TMath::Abs(pdgcode)==13)
923 fTOFyieldSelPimcM->Fill(pt);
925 if (TMath::Abs(pdgcode)==211)
927 fTOFyieldSelPimcPi->Fill(pt);
929 else if(TMath::Abs(pdgcode)==321)
931 fTOFyieldSelPimcK->Fill(pt);
933 else if (TMath::Abs(pdgcode)==2212)
935 fTOFyieldSelPimcP->Fill(pt);
939 if (fCutsTOFKaons->PassesTOFbetaCut(track))
941 fTPCdedxAfterTOFpidKaons->Fill(pinTPCglobal,tpcSignal);
942 fTOFbetaAfterKaonCuts->Fill(p,beta);
943 fTOFbetaKafter->Fill(p,beta-betaHypothesis[3]);
946 if (TMath::Abs(pdgcode)==11)
948 fTOFyieldSelKmcE->Fill(pt);
950 if (TMath::Abs(pdgcode)==13)
952 fTOFyieldSelKmcM->Fill(pt);
954 if (TMath::Abs(pdgcode)==211)
956 fTOFyieldSelKmcPi->Fill(pt);
958 else if(TMath::Abs(pdgcode)==321)
960 fTOFyieldSelKmcK->Fill(pt);
962 else if (TMath::Abs(pdgcode)==2212)
964 fTOFyieldSelKmcP->Fill(pt);
968 if (fCutsTOFProtons->PassesTOFbetaCut(track))
970 fTPCdedxAfterTOFpidProtons->Fill(pinTPCglobal,tpcSignal);
971 fTOFbetaAfterProtonCuts->Fill(p,beta);
972 fTOFbetaPafter->Fill(p,beta-betaHypothesis[4]);
975 if (TMath::Abs(pdgcode)==11)
977 fTOFyieldSelPmcE->Fill(pt);
979 if (TMath::Abs(pdgcode)==13)
981 fTOFyieldSelPmcM->Fill(pt);
983 if (TMath::Abs(pdgcode)==211)
985 fTOFyieldSelPmcPi->Fill(pt);
987 else if(TMath::Abs(pdgcode)==321)
989 fTOFyieldSelPmcK->Fill(pt);
991 else if (TMath::Abs(pdgcode)==2212)
993 fTOFyieldSelPmcP->Fill(pt);
1000 //______________________________________________________________________________
1001 Float_t AliAnalysisTaskPIDflowQA::Beta(Float_t m, Float_t p)
1003 //get theoretical beta
1004 return TMath::Sqrt(1. / (1. + m * m / (p * p)));