]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskPIDqa.cxx
- Reshuffling of the particle codes in AliPID. Now the light nuclei are between the
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisTaskPIDqa.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 \r
16 /* $Id: AliAnalysisTaskPIDqa.cxx 43811 2010-09-23 14:13:31Z wiechula $ */\r
17 #include <TList.h>\r
18 #include <TVectorD.h>\r
19 #include <TObjArray.h>\r
20 #include <TH2.h>\r
21 #include <TFile.h>\r
22 #include <TPRegexp.h>\r
23 #include <TChain.h>\r
24 #include <TF1.h>\r
25 #include <TSpline.h>\r
26 \r
27 #include <AliAnalysisManager.h>\r
28 #include <AliInputEventHandler.h>\r
29 #include <AliVEventHandler.h>\r
30 #include <AliVEvent.h>\r
31 #include <AliVParticle.h>\r
32 #include <AliVTrack.h>\r
33 #include <AliLog.h>\r
34 #include <AliPID.h>\r
35 #include <AliPIDResponse.h>\r
36 #include <AliITSPIDResponse.h>\r
37 #include <AliTPCPIDResponse.h>\r
38 #include <AliTRDPIDResponse.h>\r
39 #include <AliTOFPIDResponse.h>\r
40 \r
41 #include <AliESDEvent.h>\r
42 #include <AliAODEvent.h>\r
43 #include <AliESDv0.h>\r
44 #include <AliAODv0.h>\r
45 #include <AliESDv0KineCuts.h>\r
46 \r
47 #include "AliAnalysisTaskPIDqa.h"\r
48 \r
49 \r
50 ClassImp(AliAnalysisTaskPIDqa)\r
51 \r
52 //______________________________________________________________________________\r
53 AliAnalysisTaskPIDqa::AliAnalysisTaskPIDqa():\r
54 AliAnalysisTaskSE(),\r
55 fPIDResponse(0x0),\r
56 fV0cuts(0x0),\r
57 fV0electrons(0x0),\r
58 fV0pions(0x0),\r
59 fV0kaons(0x0),\r
60 fV0protons(0x0),\r
61 fListQA(0x0),\r
62 fListQAits(0x0),\r
63 fListQAitsSA(0x0),\r
64 fListQAitsPureSA(0x0),\r
65 fListQAtpc(0x0),\r
66 fListQAtrd(0x0),\r
67 fListQAtof(0x0),\r
68 fListQAemcal(0x0),\r
69 fListQAhmpid(0x0),\r
70 fListQAtofhmpid(0x0),\r
71 fListQAtpctof(0x0),\r
72 fListQAV0(0x0)\r
73 {\r
74   //\r
75   // Dummy constructor\r
76   //\r
77 }\r
78 \r
79 //______________________________________________________________________________\r
80 AliAnalysisTaskPIDqa::AliAnalysisTaskPIDqa(const char* name):\r
81 AliAnalysisTaskSE(name),\r
82 fPIDResponse(0x0),\r
83 fV0cuts(0x0),\r
84 fV0electrons(0x0),\r
85 fV0pions(0x0),\r
86 fV0kaons(0x0),\r
87 fV0protons(0x0),\r
88 fListQA(0x0),\r
89 fListQAits(0x0),\r
90 fListQAitsSA(0x0),\r
91 fListQAitsPureSA(0x0),\r
92 fListQAtpc(0x0),\r
93 fListQAtrd(0x0),\r
94 fListQAtof(0x0),\r
95 fListQAemcal(0x0),\r
96 fListQAhmpid(0x0),\r
97 fListQAtofhmpid(0x0),\r
98 fListQAtpctof(0x0),\r
99 fListQAV0(0x0)\r
100 {\r
101   //\r
102   // Default constructor\r
103   //\r
104   DefineInput(0,TChain::Class());\r
105   DefineOutput(1,TList::Class());\r
106 }\r
107 \r
108 //______________________________________________________________________________\r
109 AliAnalysisTaskPIDqa::~AliAnalysisTaskPIDqa()\r
110 {\r
111   //\r
112   // Destructor\r
113   //\r
114 \r
115   delete fV0cuts;\r
116   delete fV0electrons;\r
117   delete fV0pions;\r
118   delete fV0kaons;\r
119   delete fV0protons;\r
120 \r
121   if (!AliAnalysisManager::GetAnalysisManager()->IsProofMode()) delete fListQA;\r
122 }\r
123 \r
124 //______________________________________________________________________________\r
125 void AliAnalysisTaskPIDqa::UserCreateOutputObjects()\r
126 {\r
127   //\r
128   // Create the output QA objects\r
129   //\r
130 \r
131   AliLog::SetClassDebugLevel("AliAnalysisTaskPIDqa",10);\r
132 \r
133   //input hander\r
134   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();\r
135   AliInputEventHandler *inputHandler=dynamic_cast<AliInputEventHandler*>(man->GetInputEventHandler());\r
136   if (!inputHandler) AliFatal("Input handler needed");\r
137 \r
138   //pid response object\r
139   fPIDResponse=inputHandler->GetPIDResponse();\r
140   if (!fPIDResponse) AliError("PIDResponse object was not created");\r
141   \r
142   // V0 Kine cuts \r
143   fV0cuts = new AliESDv0KineCuts;\r
144  \r
145   // V0 PID Obj arrays\r
146   fV0electrons = new TObjArray;\r
147   fV0pions     = new TObjArray;\r
148   fV0kaons     = new TObjArray;\r
149   fV0protons   = new TObjArray;\r
150 \r
151   //\r
152   fListQA=new TList;\r
153   fListQA->SetOwner();\r
154   \r
155   fListQAits=new TList;\r
156   fListQAits->SetOwner();\r
157   fListQAits->SetName("ITS");\r
158 \r
159   fListQAitsSA=new TList;\r
160   fListQAitsSA->SetOwner();\r
161   fListQAitsSA->SetName("ITS_SA");\r
162 \r
163   fListQAitsPureSA=new TList;\r
164   fListQAitsPureSA->SetOwner();\r
165   fListQAitsPureSA->SetName("ITS_PureSA");\r
166 \r
167   fListQAtpc=new TList;\r
168   fListQAtpc->SetOwner();\r
169   fListQAtpc->SetName("TPC");\r
170   \r
171   fListQAtrd=new TList;\r
172   fListQAtrd->SetOwner();\r
173   fListQAtrd->SetName("TRD");\r
174   \r
175   fListQAtof=new TList;\r
176   fListQAtof->SetOwner();\r
177   fListQAtof->SetName("TOF");\r
178   \r
179   fListQAemcal=new TList;\r
180   fListQAemcal->SetOwner();\r
181   fListQAemcal->SetName("EMCAL");\r
182   \r
183   fListQAhmpid=new TList;\r
184   fListQAhmpid->SetOwner();\r
185   fListQAhmpid->SetName("HMPID");\r
186   \r
187   fListQAtpctof=new TList;\r
188   fListQAtpctof->SetOwner();\r
189   fListQAtpctof->SetName("TPC_TOF");\r
190 \r
191   fListQAtofhmpid=new TList;\r
192   fListQAtofhmpid->SetOwner();\r
193   fListQAtofhmpid->SetName("TOF_HMPID");\r
194   \r
195   fListQAV0=new TList;\r
196   fListQAV0->SetOwner();\r
197   fListQAV0->SetName("V0decay");\r
198   \r
199   fListQA->Add(fListQAits);\r
200   fListQA->Add(fListQAitsSA);\r
201   fListQA->Add(fListQAitsPureSA);\r
202   fListQA->Add(fListQAtpc);\r
203   fListQA->Add(fListQAtrd);\r
204   fListQA->Add(fListQAtof);\r
205   fListQA->Add(fListQAemcal);\r
206   fListQA->Add(fListQAhmpid);\r
207   fListQA->Add(fListQAtpctof);\r
208   fListQA->Add(fListQAtofhmpid);\r
209   fListQA->Add(fListQAV0);\r
210 \r
211   SetupITSqa();\r
212   SetupTPCqa();\r
213   SetupTRDqa();\r
214   SetupTOFqa();\r
215   SetupEMCALqa();\r
216   SetupHMPIDqa();\r
217   SetupTPCTOFqa();\r
218   SetupTOFHMPIDqa();\r
219   SetupV0qa();\r
220   \r
221   PostData(1,fListQA);\r
222 }\r
223 \r
224 \r
225 //______________________________________________________________________________\r
226 void AliAnalysisTaskPIDqa::UserExec(Option_t */*option*/)\r
227 {\r
228   //\r
229   // Setup the PID response functions and fill the QA histograms\r
230   //\r
231 \r
232   AliVEvent *event=InputEvent();\r
233   if (!event||!fPIDResponse) return;\r
234 \r
235   // Start with the V0 task (only possible for ESDs?)\r
236   FillV0PIDlist();\r
237   \r
238   FillITSqa();\r
239   FillTPCqa();\r
240   FillTRDqa();\r
241   FillTOFqa();\r
242   FillEMCALqa();\r
243   FillHMPIDqa();\r
244   \r
245   //combined detector QA\r
246   FillTPCTOFqa();\r
247   FillTOFHMPIDqa();\r
248   \r
249   // Clear the V0 PID arrays\r
250   ClearV0PIDlist();\r
251 \r
252 \r
253   \r
254   PostData(1,fListQA);\r
255 }\r
256 \r
257 //______________________________________________________________________________\r
258 void  AliAnalysisTaskPIDqa::FillV0PIDlist(){\r
259 \r
260   //\r
261   // Fill the PID object arrays holding the pointers to identified particle tracks\r
262   //\r
263 \r
264   // Dynamic cast to ESD events (DO NOTHING for AOD events)\r
265   AliESDEvent *event = dynamic_cast<AliESDEvent *>(InputEvent());\r
266   if ( !event )  return;\r
267   \r
268   if(TString(event->GetBeamType())=="Pb-Pb" || TString(event->GetBeamType())=="A-A"){\r
269     fV0cuts->SetMode(AliESDv0KineCuts::kPurity,AliESDv0KineCuts::kPbPb); \r
270   }\r
271   else{\r
272     fV0cuts->SetMode(AliESDv0KineCuts::kPurity,AliESDv0KineCuts::kPP); \r
273   }\r
274 \r
275   // V0 selection\r
276   // set event\r
277   fV0cuts->SetEvent(event);\r
278 \r
279   // loop over V0 particles\r
280   for(Int_t iv0=0; iv0<event->GetNumberOfV0s();iv0++){\r
281 \r
282     AliESDv0 *v0 = (AliESDv0 *) event->GetV0(iv0);\r
283  \r
284     if(!v0) continue;\r
285     if(v0->GetOnFlyStatus()) continue; \r
286   \r
287     // Get the particle selection \r
288     Bool_t foundV0 = kFALSE;\r
289     Int_t pdgV0, pdgP, pdgN;\r
290 \r
291     foundV0 = fV0cuts->ProcessV0(v0, pdgV0, pdgP, pdgN);\r
292     if(!foundV0) continue;\r
293     \r
294     Int_t iTrackP = v0->GetPindex();  // positive track\r
295     Int_t iTrackN = v0->GetNindex();  // negative track\r
296 \r
297     // v0 Armenteros plot (QA)\r
298     Float_t armVar[2] = {0.0,0.0};\r
299     fV0cuts->Armenteros(v0, armVar);\r
300 \r
301     TH2 *h=(TH2*)fListQAV0->At(0);\r
302     if (!h) continue;\r
303     h->Fill(armVar[0],armVar[1]);\r
304 \r
305     // fill the Object arrays\r
306     // positive particles\r
307     if( pdgP == -11){\r
308       fV0electrons->Add((AliVTrack*)event->GetTrack(iTrackP));\r
309     }\r
310     else if( pdgP == 211){\r
311       fV0pions->Add((AliVTrack*)event->GetTrack(iTrackP));\r
312     }\r
313     else if( pdgP == 321){\r
314       fV0kaons->Add((AliVTrack*)event->GetTrack(iTrackP));\r
315     }\r
316     else if( pdgP == 2212){\r
317       fV0protons->Add((AliVTrack*)event->GetTrack(iTrackP));\r
318     }\r
319 \r
320     // negative particles\r
321     if( pdgN == 11){\r
322       fV0electrons->Add((AliVTrack*)event->GetTrack(iTrackN));\r
323     }\r
324     else if( pdgN == -211){\r
325       fV0pions->Add((AliVTrack*)event->GetTrack(iTrackN));\r
326     }\r
327     else if( pdgN == -321){\r
328       fV0kaons->Add((AliVTrack*)event->GetTrack(iTrackN));\r
329     }\r
330     else if( pdgN == -2212){\r
331       fV0protons->Add((AliVTrack*)event->GetTrack(iTrackN));\r
332     }\r
333   \r
334 \r
335   }\r
336 }\r
337 //______________________________________________________________________________\r
338 void  AliAnalysisTaskPIDqa::ClearV0PIDlist(){\r
339 \r
340   //\r
341   // Clear the PID object arrays\r
342   //\r
343 \r
344   fV0electrons->Clear();\r
345   fV0pions->Clear();\r
346   fV0kaons->Clear();\r
347   fV0protons->Clear();\r
348 \r
349 }\r
350 //______________________________________________________________________________\r
351 void AliAnalysisTaskPIDqa::FillITSqa()\r
352 {\r
353   //\r
354   // Fill PID qa histograms for the ITS\r
355   //\r
356 \r
357   AliVEvent *event=InputEvent();\r
358   \r
359   Int_t ntracks=event->GetNumberOfTracks();\r
360   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
361     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
362     ULong_t status=track->GetStatus();\r
363     // not that nice. status bits not in virtual interface\r
364     // ITS refit + ITS pid selection\r
365     if (!( ( (status & AliVTrack::kITSrefit)==AliVTrack::kITSrefit ) ||\r
366            ! ( (status & AliVTrack::kITSpid  )==AliVTrack::kITSpid   ) )) continue;\r
367     Double_t mom=track->P();\r
368     \r
369     TList *theList = 0x0;\r
370     if(( (status & AliVTrack::kTPCin)==AliVTrack::kTPCin )){\r
371       //ITS+TPC tracks\r
372       theList=fListQAits;\r
373     }else{\r
374       if(!( (status & AliVTrack::kITSpureSA)==AliVTrack::kITSpureSA )){ \r
375         //ITS Standalone tracks\r
376         theList=fListQAitsSA;\r
377       }else{\r
378         //ITS Pure Standalone tracks\r
379         theList=fListQAitsPureSA;\r
380       }\r
381     }\r
382     \r
383     \r
384     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
385       TH2 *h=(TH2*)theList->At(ispecie);\r
386       if (!h) continue;\r
387       Double_t nSigma=fPIDResponse->NumberOfSigmasITS(track, (AliPID::EParticleType)ispecie);\r
388       h->Fill(mom,nSigma);\r
389     }\r
390     TH2 *h=(TH2*)theList->At(AliPID::kSPECIESC);\r
391     if (h) {\r
392       Double_t sig=track->GetITSsignal();\r
393       h->Fill(mom,sig);\r
394     }\r
395   }\r
396 }\r
397 \r
398 //______________________________________________________________________________\r
399 void AliAnalysisTaskPIDqa::FillTPCqa()\r
400 {\r
401   //\r
402   // Fill PID qa histograms for the TPC\r
403   //\r
404   \r
405   AliVEvent *event=InputEvent();\r
406   \r
407   Int_t ntracks=event->GetNumberOfTracks();\r
408   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
409     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
410     \r
411     //\r
412     //basic track cuts\r
413     //\r
414     ULong_t status=track->GetStatus();\r
415     // not that nice. status bits not in virtual interface\r
416     // TPC refit + ITS refit + TPC pid\r
417     if (!( (status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
418         !( (status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ) continue;\r
419 \r
420     // The TPC pid cut removes the light nuclei (>5 sigma from proton line)\r
421     //||        !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid  )\r
422     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
423     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
424     if (track->GetTPCNclsF()>0) {\r
425       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
426     }\r
427     \r
428     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
429     \r
430     Double_t mom=track->GetTPCmomentum();\r
431     \r
432     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
433       TH2 *h=(TH2*)fListQAtpc->At(ispecie);\r
434       if (!h) continue;\r
435       Double_t nSigma=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie);\r
436       h->Fill(mom,nSigma);\r
437     }\r
438     \r
439     TH2 *h=(TH2*)fListQAtpc->At(AliPID::kSPECIESC);\r
440     if (h) {\r
441       Double_t sig=track->GetTPCsignal();\r
442       h->Fill(mom,sig);\r
443     }\r
444   }\r
445 }\r
446 \r
447 //______________________________________________________________________________\r
448 void AliAnalysisTaskPIDqa::FillTRDqa()\r
449 {\r
450   //\r
451   // Fill PID qa histograms for the TRD\r
452   //\r
453   AliVEvent *event=InputEvent();\r
454   Int_t ntracks = event->GetNumberOfTracks();\r
455   for(Int_t itrack = 0; itrack <  ntracks; itrack++){\r
456     AliVTrack *track = (AliVTrack *)event->GetTrack(itrack);\r
457 \r
458     //\r
459     //basic track cuts\r
460     //\r
461     ULong_t status=track->GetStatus();\r
462     // not that nice. status bits not in virtual interface\r
463     // TPC refit + ITS refit + TPC pid + TRD out\r
464     if (!( (status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
465         !( (status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
466 //         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid  ) || //removes light nuclei. So it is out for the moment\r
467         !( (status & AliVTrack::kTRDout  ) == AliVTrack::kTRDout  )) continue;\r
468     \r
469     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
470     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
471     if (track->GetTPCNclsF()>0) {\r
472       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
473     }\r
474     \r
475     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
476 \r
477     Double_t likelihoods[AliPID::kSPECIES];\r
478     if(fPIDResponse->ComputeTRDProbability(track, AliPID::kSPECIES, likelihoods) != AliPIDResponse::kDetPidOk) continue;\r
479     Int_t ntracklets = 0;\r
480     Double_t momentum = -1.;\r
481     for(Int_t itl = 0; itl < 6; itl++)\r
482       if(track->GetTRDmomentum(itl) > 0.){\r
483         ntracklets++;\r
484         if(momentum < 0) momentum = track->GetTRDmomentum(itl);\r
485     } \r
486     for(Int_t ispecie = 0; ispecie < AliPID::kSPECIES; ispecie++){\r
487       TH2F *hLike = (TH2F *)fListQAtrd->At(ntracklets*AliPID::kSPECIES+ispecie);\r
488       if (hLike) hLike->Fill(momentum,likelihoods[ispecie]);\r
489     }\r
490   }\r
491 }\r
492 \r
493 //______________________________________________________________________________\r
494 void AliAnalysisTaskPIDqa::FillTOFqa()\r
495 {\r
496   //\r
497   // Fill TOF information\r
498   //\r
499   AliVEvent *event=InputEvent();\r
500 \r
501   Int_t ntracks=event->GetNumberOfTracks();\r
502   Int_t tracksAtTof = 0;\r
503   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
504     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
505 \r
506     //\r
507     //basic track cuts\r
508     //\r
509     ULong_t status=track->GetStatus();\r
510     // TPC refit + ITS refit +\r
511     // TOF out + TOFpid +\r
512     // kTIME\r
513     // (we don't use kTOFmismatch because it depends on TPC....)\r
514     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
515         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
516         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
517         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
518         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
519 \r
520     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
521     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
522     if (track->GetTPCNclsF()>0) {\r
523       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
524     }\r
525 \r
526     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
527 \r
528     tracksAtTof++;\r
529 \r
530     Double_t mom=track->P();\r
531 \r
532     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
533       TH2 *h=(TH2*)fListQAtof->At(ispecie);\r
534       if (!h) continue;\r
535       Double_t nSigma=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
536       h->Fill(mom,nSigma);\r
537     }\r
538 \r
539     TH2 *h=(TH2*)fListQAtof->FindObject("hSigP_TOF");\r
540     if (h) {\r
541       Double_t sig=track->GetTOFsignal()/1000.;\r
542       h->Fill(mom,sig);\r
543     }\r
544 \r
545     Int_t mask = fPIDResponse->GetTOFResponse().GetStartTimeMask(mom);\r
546     ((TH1F*)fListQAtof->FindObject("hStartTimeMask_TOF"))->Fill((Double_t)(mask+0.5));\r
547 \r
548     if (mom >= 0.75 && mom <= 1.25 ) {\r
549       Double_t nsigma= fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)AliPID::kPion);\r
550       if (mask == 0) {\r
551         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-Fill"))->Fill(nsigma);\r
552       } else if (mask == 1) {\r
553         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-TOF"))->Fill(nsigma);\r
554       } else if ( (mask == 2) || (mask == 4) || (mask == 6) ) {\r
555         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-T0"))->Fill(nsigma);\r
556       } else {\r
557         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-Best"))->Fill(nsigma);\r
558       }\r
559     }\r
560 \r
561     Double_t res = (Double_t)fPIDResponse->GetTOFResponse().GetStartTimeRes(mom);\r
562     ((TH1F*)fListQAtof->FindObject("hStartTimeRes_TOF"))->Fill(res);\r
563 \r
564     AliESDEvent *esd = dynamic_cast<AliESDEvent *>(event);\r
565     if (esd) {\r
566       Double_t startTime = esd->GetT0TOF(0);\r
567       if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeAC_T0"))->Fill(startTime);\r
568       else {\r
569         startTime = esd->GetT0TOF(1);\r
570         if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeA_T0"))->Fill(startTime);\r
571         startTime = esd->GetT0TOF(2);\r
572         if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeC_T0"))->Fill(startTime);\r
573       }\r
574     }\r
575   }\r
576   if (tracksAtTof > 0) {\r
577     ((TH1F* )fListQAtof->FindObject("hnTracksAt_TOF"))->Fill(tracksAtTof);\r
578     Int_t mask = fPIDResponse->GetTOFResponse().GetStartTimeMask(5.);\r
579     if (mask & 0x1) ((TH1F*)fListQAtof->FindObject("hT0MakerEff"))->Fill(tracksAtTof);\r
580   }\r
581 }\r
582 \r
583 \r
584 //______________________________________________________________________________\r
585 void AliAnalysisTaskPIDqa::FillEMCALqa()\r
586 {\r
587   //\r
588   // Fill PID qa histograms for the EMCAL\r
589   //\r
590 \r
591   AliVEvent *event=InputEvent();\r
592   \r
593   Int_t ntracks=event->GetNumberOfTracks();\r
594   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
595     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
596     \r
597     //\r
598     //basic track cuts\r
599     //\r
600     ULong_t status=track->GetStatus();\r
601     // not that nice. status bits not in virtual interface\r
602     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
603 \r
604     Double_t pt=track->Pt();\r
605    \r
606     //EMCAL nSigma (only for electrons at the moment)\r
607     TH2 *h=(TH2*)fListQAemcal->At(0);\r
608     if (!h) continue;\r
609     Double_t nSigma=fPIDResponse->NumberOfSigmasEMCAL(track, (AliPID::EParticleType)0);\r
610     h->Fill(pt,nSigma);\r
611     \r
612   }\r
613 \r
614    //EMCAL signal (E/p vs. pT) for electrons from V0\r
615   for(Int_t itrack = 0; itrack < fV0electrons->GetEntries(); itrack++){\r
616     AliVTrack *track=(AliVTrack*)fV0electrons->At(itrack);\r
617 \r
618     //\r
619     //basic track cuts\r
620     //\r
621     ULong_t status=track->GetStatus();\r
622     // not that nice. status bits not in virtual interface\r
623     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
624 \r
625     Double_t pt=track->Pt();\r
626 \r
627     TH2 *h=(TH2*)fListQAemcal->At(1);\r
628     if (h) {\r
629 \r
630       Int_t nMatchClus = track->GetEMCALcluster();\r
631       Double_t mom     = track->P();\r
632       Double_t eop     = -1.;\r
633 \r
634       if(nMatchClus > -1){\r
635     \r
636         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
637 \r
638         if(matchedClus){\r
639 \r
640           // matched cluster is EMCAL\r
641           if(matchedClus->IsEMCAL()){\r
642 \r
643             Double_t fClsE       = matchedClus->E();\r
644             eop                  = fClsE/mom;\r
645 \r
646             h->Fill(pt,eop);\r
647 \r
648           }\r
649         }\r
650       }\r
651     }\r
652   }\r
653 \r
654    //EMCAL signal (E/p vs. pT) for pions from V0\r
655   for(Int_t itrack = 0; itrack < fV0pions->GetEntries(); itrack++){\r
656     AliVTrack *track=(AliVTrack*)fV0pions->At(itrack);\r
657 \r
658     //\r
659     //basic track cuts\r
660     //\r
661     ULong_t status=track->GetStatus();\r
662     // not that nice. status bits not in virtual interface\r
663     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
664 \r
665     Double_t pt=track->Pt();\r
666 \r
667     TH2 *h=(TH2*)fListQAemcal->At(2);\r
668     if (h) {\r
669 \r
670       Int_t nMatchClus = track->GetEMCALcluster();\r
671       Double_t mom     = track->P();\r
672       Double_t eop     = -1.;\r
673 \r
674       if(nMatchClus > -1){\r
675     \r
676         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
677 \r
678         if(matchedClus){\r
679 \r
680           // matched cluster is EMCAL\r
681           if(matchedClus->IsEMCAL()){\r
682 \r
683             Double_t fClsE       = matchedClus->E();\r
684             eop                  = fClsE/mom;\r
685 \r
686             h->Fill(pt,eop);\r
687 \r
688           }\r
689         }\r
690       }\r
691     }\r
692   }\r
693 \r
694    //EMCAL signal (E/p vs. pT) for protons from V0\r
695   for(Int_t itrack = 0; itrack < fV0protons->GetEntries(); itrack++){\r
696     AliVTrack *track=(AliVTrack*)fV0protons->At(itrack);\r
697 \r
698     //\r
699     //basic track cuts\r
700     //\r
701     ULong_t status=track->GetStatus();\r
702     // not that nice. status bits not in virtual interface\r
703     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
704 \r
705     Double_t pt=track->Pt();\r
706 \r
707     TH2 *hP=(TH2*)fListQAemcal->At(3);\r
708     TH2 *hAP=(TH2*)fListQAemcal->At(4);\r
709     if (hP && hAP) {\r
710 \r
711       Int_t nMatchClus = track->GetEMCALcluster();\r
712       Double_t mom     = track->P();\r
713       Int_t charge     = track->Charge();             \r
714       Double_t eop     = -1.;\r
715 \r
716       if(nMatchClus > -1){\r
717     \r
718         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
719 \r
720         if(matchedClus){\r
721 \r
722           // matched cluster is EMCAL\r
723           if(matchedClus->IsEMCAL()){\r
724 \r
725             Double_t fClsE       = matchedClus->E();\r
726             eop                  = fClsE/mom;\r
727 \r
728             if(charge > 0)      hP->Fill(pt,eop);\r
729             else if(charge < 0) hAP->Fill(pt,eop);\r
730 \r
731           }\r
732         }\r
733       }\r
734     }\r
735   }\r
736 \r
737 }\r
738 \r
739 \r
740 //______________________________________________________________________________\r
741 void AliAnalysisTaskPIDqa::FillHMPIDqa()\r
742 {\r
743   //\r
744   // Fill PID qa histograms for the HMPID\r
745   //\r
746   \r
747   AliVEvent *event=InputEvent();\r
748   \r
749   Int_t ntracks=event->GetNumberOfTracks();\r
750   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
751     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
752     \r
753     //\r
754     //basic track cuts\r
755     //\r
756     ULong_t status=track->GetStatus();\r
757     // not that nice. status bits not in virtual interface\r
758     // TPC refit + ITS refit +\r
759     // TOF out + TOFpid +\r
760     // kTIME\r
761     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
762         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ) continue;\r
763 \r
764     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
765     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
766     if (track->GetTPCNclsF()>0) {\r
767       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
768     }\r
769 \r
770     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
771     \r
772     Double_t mom = track->P();\r
773     Double_t ckovAngle = track->GetHMPIDsignal();\r
774     \r
775     TH1F *hThetavsMom = (TH1F*)fListQAhmpid->At(0);;\r
776     \r
777     hThetavsMom->Fill(mom,ckovAngle);    \r
778   \r
779   }\r
780 }\r
781 //______________________________________________________________________________\r
782 void AliAnalysisTaskPIDqa::FillTOFHMPIDqa()\r
783 {\r
784   //\r
785   // Fill PID qa histograms for the HMPID\r
786   //\r
787   \r
788   AliVEvent *event=InputEvent();\r
789   \r
790   Int_t ntracks=event->GetNumberOfTracks();\r
791   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
792     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
793     \r
794     //\r
795     //basic track cuts\r
796     //\r
797     ULong_t status=track->GetStatus();\r
798     // not that nice. status bits not in virtual interface\r
799     // TPC refit + ITS refit +\r
800     // TOF out + TOFpid +\r
801     // kTIME\r
802     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
803         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
804         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
805         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
806         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
807 \r
808     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
809     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
810     if (track->GetTPCNclsF()>0) {\r
811       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
812     }\r
813 \r
814     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
815     \r
816     Double_t mom = track->P();\r
817     Double_t ckovAngle = track->GetHMPIDsignal();\r
818     \r
819     Double_t nSigmaTOF[3]; \r
820     TH1F *h[3];\r
821     \r
822     for (Int_t ispecie=2; ispecie<5; ++ispecie){\r
823       //TOF nSigma\r
824       nSigmaTOF[ispecie-2]=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
825       h[ispecie-2] = (TH1F*)fListQAtofhmpid->At(ispecie-2);}\r
826       \r
827     if(TMath::Abs(nSigmaTOF[0])<2)                                                              h[0]->Fill(mom,ckovAngle);\r
828     \r
829     if(TMath::Abs(nSigmaTOF[1])<2 && TMath::Abs(nSigmaTOF[0])>3)                                h[1]->Fill(mom,ckovAngle);\r
830 \r
831     if(TMath::Abs(nSigmaTOF[2])<2 && TMath::Abs(nSigmaTOF[1])>3 && TMath::Abs(nSigmaTOF[0])>3)  h[2]->Fill(mom,ckovAngle);\r
832       \r
833   }\r
834   \r
835 }\r
836 \r
837 //______________________________________________________________________________\r
838 void AliAnalysisTaskPIDqa::FillTPCTOFqa()\r
839 {\r
840   //\r
841   // Fill PID qa histograms for the TOF\r
842   //   Here also the TPC histograms after TOF selection are filled\r
843   //\r
844 \r
845   AliVEvent *event=InputEvent();\r
846 \r
847   Int_t ntracks=event->GetNumberOfTracks();\r
848   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
849     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
850 \r
851     //\r
852     //basic track cuts\r
853     //\r
854     ULong_t status=track->GetStatus();\r
855     // not that nice. status bits not in virtual interface\r
856     // TPC refit + ITS refit +\r
857     // TOF out + TOFpid +\r
858     // kTIME\r
859     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
860         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
861 //         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid ) || //removes light nuclei, so it is out for the moment\r
862         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
863         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
864         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
865 \r
866     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
867     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
868     if (track->GetTPCNclsF()>0) {\r
869       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
870     }\r
871 \r
872     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
873 \r
874 \r
875     Double_t mom=track->P();\r
876     Double_t momTPC=track->GetTPCmomentum();\r
877 \r
878     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
879       //TOF nSigma\r
880       Double_t nSigmaTOF=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
881       Double_t nSigmaTPC=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie);\r
882 \r
883       //TPC after TOF cut\r
884       TH2 *h=(TH2*)fListQAtpctof->At(ispecie);\r
885       if (h && TMath::Abs(nSigmaTOF)<3.) h->Fill(momTPC,nSigmaTPC);\r
886 \r
887       //TOF after TPC cut\r
888       h=(TH2*)fListQAtpctof->At(ispecie+AliPID::kSPECIESC);\r
889       if (h && TMath::Abs(nSigmaTPC)<3.) h->Fill(mom,nSigmaTOF);\r
890 \r
891       //EMCAL after TOF and TPC cut\r
892       h=(TH2*)fListQAtpctof->At(ispecie+2*AliPID::kSPECIESC);\r
893       if (h && TMath::Abs(nSigmaTOF)<3. && TMath::Abs(nSigmaTPC)<3. ){\r
894 \r
895         Int_t nMatchClus = track->GetEMCALcluster();\r
896         Double_t pt      = track->Pt();\r
897         Double_t eop     = -1.;\r
898         \r
899         if(nMatchClus > -1){\r
900           \r
901           AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
902           \r
903           if(matchedClus){\r
904             \r
905             // matched cluster is EMCAL\r
906             if(matchedClus->IsEMCAL()){\r
907               \r
908               Double_t fClsE       = matchedClus->E();\r
909               eop                  = fClsE/mom;\r
910 \r
911               h->Fill(pt,eop);\r
912  \r
913               \r
914             }\r
915           }\r
916         }\r
917       }\r
918     }\r
919   }\r
920 }\r
921 \r
922 //______________________________________________________________________________\r
923 void AliAnalysisTaskPIDqa::SetupITSqa()\r
924 {\r
925   //\r
926   // Create the ITS qa objects\r
927   //\r
928   \r
929   TVectorD *vX=MakeLogBinning(200,.1,30);\r
930   \r
931   //ITS+TPC tracks\r
932   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
933     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_ITS_%s",AliPID::ParticleName(ispecie)),\r
934                               Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
935                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
936                               200,-10,10);\r
937     fListQAits->Add(hNsigmaP);\r
938   }\r
939   TH2F *hSig = new TH2F("hSigP_ITS",\r
940                         "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
941                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
942                         300,0,300);\r
943   fListQAits->Add(hSig);\r
944 \r
945   //ITS Standalone tracks\r
946   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
947     TH2F *hNsigmaPSA = new TH2F(Form("hNsigmaP_ITSSA_%s",AliPID::ParticleName(ispecie)),\r
948                                 Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
949                                 vX->GetNrows()-1,vX->GetMatrixArray(),\r
950                                 200,-10,10);\r
951     fListQAitsSA->Add(hNsigmaPSA);\r
952   }\r
953   TH2F *hSigSA = new TH2F("hSigP_ITSSA",\r
954                           "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
955                           vX->GetNrows()-1,vX->GetMatrixArray(),\r
956                           300,0,300);\r
957   fListQAitsSA->Add(hSigSA);\r
958   \r
959   //ITS Pure Standalone tracks\r
960   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
961     TH2F *hNsigmaPPureSA = new TH2F(Form("hNsigmaP_ITSPureSA_%s",AliPID::ParticleName(ispecie)),\r
962                                     Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
963                                     vX->GetNrows()-1,vX->GetMatrixArray(),\r
964                                     200,-10,10);\r
965     fListQAitsPureSA->Add(hNsigmaPPureSA);\r
966   }\r
967   TH2F *hSigPureSA = new TH2F("hSigP_ITSPureSA",\r
968                               "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
969                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
970                               300,0,300);\r
971   fListQAitsPureSA->Add(hSigPureSA);\r
972   \r
973   delete vX;  \r
974 }\r
975 \r
976 //______________________________________________________________________________\r
977 void AliAnalysisTaskPIDqa::SetupTPCqa()\r
978 {\r
979   //\r
980   // Create the TPC qa objects\r
981   //\r
982   \r
983   TVectorD *vX=MakeLogBinning(200,.1,30);\r
984   \r
985   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
986     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_%s",AliPID::ParticleName(ispecie)),\r
987                               Form("TPC n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
988                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
989                               200,-10,10);\r
990     fListQAtpc->Add(hNsigmaP);\r
991   }\r
992   \r
993   \r
994   TH2F *hSig = new TH2F("hSigP_TPC",\r
995                         "TPC signal vs. p;p [GeV]; TPC signal [arb. units]",\r
996                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
997                         300,0,300);\r
998   fListQAtpc->Add(hSig);\r
999 \r
1000   delete vX;  \r
1001 }\r
1002 \r
1003 //______________________________________________________________________________\r
1004 void AliAnalysisTaskPIDqa::SetupTRDqa()\r
1005 {\r
1006   //\r
1007   // Create the TRD qa objects\r
1008   //\r
1009   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1010   for(Int_t itl = 0; itl < 6; ++itl){\r
1011     for(Int_t ispecie = 0; ispecie < AliPID::kSPECIES; ispecie++){\r
1012       TH2F *hLikeP = new TH2F(Form("hLikeP_TRD_%dtls_%s", itl, AliPID::ParticleName(ispecie)),\r
1013                               Form("TRD Likelihood to be %s %s for tracks having %d %s; p (GeV/c); TRD %s Likelihood", ispecie == 0 ? "an" : "a", AliPID::ParticleName(ispecie), itl+1, itl == 0 ? "tracklet" : "tracklets", AliPID::ParticleName(ispecie)),\r
1014                               vX->GetNrows()-1, vX->GetMatrixArray(),\r
1015                               100, 0., 1.);\r
1016       fListQAtrd->Add(hLikeP);\r
1017     }\r
1018   }\r
1019   delete vX;\r
1020 }\r
1021 \r
1022 //______________________________________________________________________________\r
1023 void AliAnalysisTaskPIDqa::SetupTOFqa()\r
1024 {\r
1025   //\r
1026   // Create the TOF qa objects\r
1027   //\r
1028   \r
1029   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1030 \r
1031   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1032     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_%s",AliPID::ParticleName(ispecie)),\r
1033                               Form("TOF n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1034                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1035                               200,-10,10);\r
1036     fListQAtof->Add(hNsigmaP);\r
1037   }\r
1038 \r
1039   // for Kaons PID we differentiate on Time Zero\r
1040   TH1F *hnSigT0Fill = new TH1F("hNsigma_TOF_Pion_T0-Fill","TOF n#sigma (Pion) T0-FILL [0.75-1.25. GeV/c]",200,-10,10);\r
1041   fListQAtof->Add(hnSigT0Fill);\r
1042   TH1F *hnSigT0T0 = new TH1F("hNsigma_TOF_Pion_T0-T0","TOF n#sigma (Pion) T0-T0 [0.75-1.25 GeV/c]",200,-10,10);\r
1043   fListQAtof->Add(hnSigT0T0);\r
1044   TH1F *hnSigT0TOF = new TH1F("hNsigma_TOF_Pion_T0-TOF","TOF n#sigma (Pion) T0-TOF [0.75-1.25 GeV/c]",200,-10,10);\r
1045   fListQAtof->Add(hnSigT0TOF);\r
1046   TH1F *hnSigT0Best = new TH1F("hNsigma_TOF_Pion_T0-Best","TOF n#sigma (Pion) T0-Best [0.75-1.25 GeV/c]",200,-10,10);\r
1047   fListQAtof->Add(hnSigT0Best);\r
1048 \r
1049   TH2F *hSig = new TH2F("hSigP_TOF",\r
1050                         "TOF signal vs. p;p [GeV]; TOF signal [ns]",\r
1051                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1052                         300,0,30);\r
1053 \r
1054   delete vX;\r
1055   \r
1056   fListQAtof->Add(hSig);\r
1057 \r
1058   TH1F *hStartTimeMaskTOF = new TH1F("hStartTimeMask_TOF","StartTime mask",8,0,8);\r
1059   fListQAtof->Add(hStartTimeMaskTOF);\r
1060   TH1F *hStartTimeResTOF = new TH1F("hStartTimeRes_TOF","StartTime resolution [ps]",100,0,500);\r
1061   fListQAtof->Add(hStartTimeResTOF);\r
1062 \r
1063   TH1F *hnTracksAtTOF = new TH1F("hnTracksAt_TOF","Matched tracks at TOF",20,0,20);\r
1064   fListQAtof->Add(hnTracksAtTOF);\r
1065   TH1F *hT0MakerEff = new TH1F("hT0MakerEff","Events with T0-TOF vs nTracks",20,0,20);\r
1066   fListQAtof->Add(hT0MakerEff);\r
1067 \r
1068   // this in principle should stay on a T0 PID QA, but are just the data prepared for TOF use\r
1069   TH1F *hStartTimeAT0 = new TH1F("hStartTimeA_T0","StartTime from T0A [ps]",1000,-1000,1000);\r
1070   fListQAtof->Add(hStartTimeAT0);\r
1071   TH1F *hStartTimeCT0 = new TH1F("hStartTimeC_T0","StartTime from T0C [ps]",1000,-1000,1000);\r
1072   fListQAtof->Add(hStartTimeCT0);\r
1073   TH1F *hStartTimeACT0 = new TH1F("hStartTimeAC_T0","StartTime from T0AC [ps]",1000,-1000,1000);;\r
1074   fListQAtof->Add(hStartTimeACT0);\r
1075 }\r
1076 \r
1077 //______________________________________________________________________________\r
1078 void AliAnalysisTaskPIDqa::SetupEMCALqa()\r
1079 {\r
1080   //\r
1081   // Create the EMCAL qa objects\r
1082   //\r
1083 \r
1084   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1085   \r
1086   TH2F *hNsigmaPt = new TH2F(Form("hNsigmaPt_EMCAL_%s",AliPID::ParticleName(0)),\r
1087                              Form("EMCAL n#sigma %s vs. p_{T};p_{T} [GeV]; n#sigma",AliPID::ParticleName(0)),\r
1088                              vX->GetNrows()-1,vX->GetMatrixArray(),\r
1089                              200,-10,10);\r
1090   fListQAemcal->Add(hNsigmaPt);  \r
1091   \r
1092 \r
1093   TH2F *hSigPtEle = new TH2F("hSigPt_EMCAL_Ele",\r
1094                         "EMCAL signal (E/p) vs. p_{T} for electrons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1095                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1096                         200,0,2);\r
1097   fListQAemcal->Add(hSigPtEle);\r
1098 \r
1099   TH2F *hSigPtPions = new TH2F("hSigPt_EMCAL_Pions",\r
1100                         "EMCAL signal (E/p) vs. p_{T} for pions;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1101                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1102                         200,0,2);\r
1103   fListQAemcal->Add(hSigPtPions);\r
1104 \r
1105   TH2F *hSigPtProtons = new TH2F("hSigPt_EMCAL_Protons",\r
1106                         "EMCAL signal (E/p) vs. p_{T} for protons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1107                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1108                         200,0,2);\r
1109   fListQAemcal->Add(hSigPtProtons);\r
1110 \r
1111   TH2F *hSigPtAntiProtons = new TH2F("hSigPt_EMCAL_Antiprotons",\r
1112                         "EMCAL signal (E/p) vs. p_{T} for antiprotons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1113                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1114                         200,0,2);\r
1115   fListQAemcal->Add(hSigPtAntiProtons);\r
1116 \r
1117   delete vX;  \r
1118 }\r
1119 \r
1120 //______________________________________________________________________________\r
1121 void AliAnalysisTaskPIDqa::SetupHMPIDqa()\r
1122 {\r
1123   //\r
1124   // Create the HMPID qa objects\r
1125   //\r
1126   \r
1127   TH2F *hCkovAnglevsMom   = new TH2F("hCkovAnglevsMom",  "Cherenkov angle vs momnetum",500,0,5.,500,0,1);\r
1128   fListQAhmpid->Add(hCkovAnglevsMom);\r
1129   \r
1130 }\r
1131 \r
1132 //______________________________________________________________________________\r
1133 void AliAnalysisTaskPIDqa::SetupTOFHMPIDqa()\r
1134 {\r
1135   //\r
1136   // Create the HMPID qa objects\r
1137   //\r
1138   \r
1139   TH2F *hCkovAnglevsMomPion   = new TH2F("hCkovAnglevsMom_pion",  "Cherenkov angle vs momnetum for pions",500,0,5.,500,0,1);\r
1140   fListQAtofhmpid->Add(hCkovAnglevsMomPion);\r
1141   \r
1142   TH2F *hCkovAnglevsMomKaon   = new TH2F("hCkovAnglevsMom_kaon",  "Cherenkov angle vs momnetum for kaons",500,0,5.,500,0,1);\r
1143   fListQAtofhmpid->Add(hCkovAnglevsMomKaon);\r
1144   \r
1145   TH2F *hCkovAnglevsMomProton = new TH2F("hCkovAnglevsMom_proton","Cherenkov angle vs momnetum for protons",500,0,5.,500,0,1);\r
1146   fListQAtofhmpid->Add(hCkovAnglevsMomProton);\r
1147   \r
1148   \r
1149 }  \r
1150 \r
1151 //______________________________________________________________________________\r
1152 void AliAnalysisTaskPIDqa::SetupTPCTOFqa()\r
1153 {\r
1154   //\r
1155   // Create the qa objects for TPC + TOF combination\r
1156   //\r
1157   \r
1158   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1159 \r
1160   //TPC signals after TOF cut\r
1161   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1162     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_TOF_%s",AliPID::ParticleName(ispecie)),\r
1163                               Form("TPC n#sigma %s vs. p (after TOF 3#sigma cut);p_{TPC} [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1164                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1165                               200,-10,10);\r
1166     fListQAtpctof->Add(hNsigmaP);\r
1167   }\r
1168 \r
1169   //TOF signals after TPC cut\r
1170   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1171     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
1172                               Form("TOF n#sigma %s vs. p (after TPC n#sigma cut);p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1173                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1174                               200,-10,10);\r
1175     fListQAtpctof->Add(hNsigmaP);\r
1176   }\r
1177 \r
1178   //EMCAL signal after TOF and TPC cut\r
1179   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
1180     TH2F *heopPt = new TH2F(Form("heopPt_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
1181                             Form("EMCAL signal (E/p) %s vs. p_{T};p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",AliPID::ParticleName(ispecie)),\r
1182                             vX->GetNrows()-1,vX->GetMatrixArray(),\r
1183                             200,0,2);\r
1184     fListQAtpctof->Add(heopPt);\r
1185   }\r
1186 \r
1187   delete vX;\r
1188 }\r
1189 //______________________________________________________________________________\r
1190 void AliAnalysisTaskPIDqa::SetupV0qa()\r
1191 {\r
1192   //\r
1193   // Create the qa objects for V0 Kine cuts\r
1194   //\r
1195   \r
1196   TH2F *hArmenteros  = new TH2F("hArmenteros",  "Armenteros plot",200,-1.,1.,200,0.,0.4);\r
1197   fListQAV0->Add(hArmenteros);\r
1198  \r
1199 }\r
1200 \r
1201 //______________________________________________________________________________\r
1202 TVectorD* AliAnalysisTaskPIDqa::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
1203 {\r
1204   //\r
1205   // Make logarithmic binning\r
1206   // the user has to delete the array afterwards!!!\r
1207   //\r
1208   \r
1209   //check limits\r
1210   if (xmin<1e-20 || xmax<1e-20){\r
1211     AliError("For Log binning xmin and xmax must be > 1e-20. Using linear binning instead!");\r
1212     return MakeLinBinning(nbinsX, xmin, xmax);\r
1213   }\r
1214   if (xmax<xmin){\r
1215     Double_t tmp=xmin;\r
1216     xmin=xmax;\r
1217     xmax=tmp;\r
1218   }\r
1219   TVectorD *binLim=new TVectorD(nbinsX+1);\r
1220   Double_t first=xmin;\r
1221   Double_t last=xmax;\r
1222   Double_t expMax=TMath::Log(last/first);\r
1223   for (Int_t i=0; i<nbinsX+1; ++i){\r
1224     (*binLim)[i]=first*TMath::Exp(expMax/nbinsX*(Double_t)i);\r
1225   }\r
1226   return binLim;\r
1227 }\r
1228 \r
1229 //______________________________________________________________________________\r
1230 TVectorD* AliAnalysisTaskPIDqa::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
1231 {\r
1232   //\r
1233   // Make linear binning\r
1234   // the user has to delete the array afterwards!!!\r
1235   //\r
1236   if (xmax<xmin){\r
1237     Double_t tmp=xmin;\r
1238     xmin=xmax;\r
1239     xmax=tmp;\r
1240   }\r
1241   TVectorD *binLim=new TVectorD(nbinsX+1);\r
1242   Double_t first=xmin;\r
1243   Double_t last=xmax;\r
1244   Double_t binWidth=(last-first)/nbinsX;\r
1245   for (Int_t i=0; i<nbinsX+1; ++i){\r
1246     (*binLim)[i]=first+binWidth*(Double_t)i;\r
1247   }\r
1248   return binLim;\r
1249 }\r
1250 \r
1251 //_____________________________________________________________________________\r
1252 TVectorD* AliAnalysisTaskPIDqa::MakeArbitraryBinning(const char* bins)\r
1253 {\r
1254   //\r
1255   // Make arbitrary binning, bins separated by a ','\r
1256   //\r
1257   TString limits(bins);\r
1258   if (limits.IsNull()){\r
1259     AliError("Bin Limit string is empty, cannot add the variable");\r
1260     return 0x0;\r
1261   }\r
1262   \r
1263   TObjArray *arr=limits.Tokenize(",");\r
1264   Int_t nLimits=arr->GetEntries();\r
1265   if (nLimits<2){\r
1266     AliError("Need at leas 2 bin limits, cannot add the variable");\r
1267     delete arr;\r
1268     return 0x0;\r
1269   }\r
1270   \r
1271   TVectorD *binLimits=new TVectorD(nLimits);\r
1272   for (Int_t iLim=0; iLim<nLimits; ++iLim){\r
1273     (*binLimits)[iLim]=(static_cast<TObjString*>(arr->At(iLim)))->GetString().Atof();\r
1274   }\r
1275   \r
1276   delete arr;\r
1277   return binLimits;\r
1278 }\r
1279 \r