]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskPIDqa.cxx
5e1d623f43aa1796d05ade72fc762c2e5e2490cb
[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     // the default scenario\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     // the "hybrid" scenario\r
439     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
440       TH2 *h=(TH2*)fListQAtpc->At(ispecie+AliPID::kSPECIESC);\r
441       if (!h) continue;\r
442       Double_t nSigma=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie, AliTPCPIDResponse::kdEdxHybrid);\r
443       h->Fill(mom,nSigma);\r
444     }\r
445     \r
446     // the "OROC" scenario\r
447     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
448       TH2 *h=(TH2*)fListQAtpc->At(ispecie+2*AliPID::kSPECIESC);\r
449       if (!h) continue;\r
450       Double_t nSigma=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie, AliTPCPIDResponse::kdEdxOROC);\r
451       //TSpline3* spline = fPIDResponse->GetTPCResponse().GetCurrentResponseFunction();\r
452       //std::cout<<ispecie<<" "<<nSigma<<" phi:"<<track->Phi()<<". "<<std::endl;\r
453       //if (spline) {cout<<spline->GetName()<<endl;}\r
454       //else {cout<<"NULL spline"<<endl;}\r
455       h->Fill(mom,nSigma);\r
456     }\r
457     \r
458     TH2 *h=(TH2*)fListQAtpc->At(3*AliPID::kSPECIESC);\r
459 \r
460     if (h) {\r
461       Double_t sig=track->GetTPCsignal();\r
462       h->Fill(mom,sig);\r
463     }\r
464   }\r
465 }\r
466 \r
467 //______________________________________________________________________________\r
468 void AliAnalysisTaskPIDqa::FillTRDqa()\r
469 {\r
470   //\r
471   // Fill PID qa histograms for the TRD\r
472   //\r
473   AliVEvent *event=InputEvent();\r
474   Int_t ntracks = event->GetNumberOfTracks();\r
475   for(Int_t itrack = 0; itrack <  ntracks; itrack++){\r
476     AliVTrack *track = (AliVTrack *)event->GetTrack(itrack);\r
477 \r
478     //\r
479     //basic track cuts\r
480     //\r
481     ULong_t status=track->GetStatus();\r
482     // not that nice. status bits not in virtual interface\r
483     // TPC refit + ITS refit + TPC pid + TRD out\r
484     if (!( (status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
485         !( (status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
486 //         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid  ) || //removes light nuclei. So it is out for the moment\r
487         !( (status & AliVTrack::kTRDout  ) == AliVTrack::kTRDout  )) continue;\r
488     \r
489     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
490     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
491     if (track->GetTPCNclsF()>0) {\r
492       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
493     }\r
494     \r
495     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
496 \r
497     Double_t likelihoods[AliPID::kSPECIES];\r
498     if(fPIDResponse->ComputeTRDProbability(track, AliPID::kSPECIES, likelihoods) != AliPIDResponse::kDetPidOk) continue;\r
499     Int_t ntracklets = 0;\r
500     Double_t momentum = -1.;\r
501     for(Int_t itl = 0; itl < 6; itl++)\r
502       if(track->GetTRDmomentum(itl) > 0.){\r
503         ntracklets++;\r
504         if(momentum < 0) momentum = track->GetTRDmomentum(itl);\r
505     } \r
506     for(Int_t ispecie = 0; ispecie < AliPID::kSPECIES; ispecie++){\r
507       TH2F *hLike = (TH2F *)fListQAtrd->At(ntracklets*AliPID::kSPECIES+ispecie);\r
508       if (hLike) hLike->Fill(momentum,likelihoods[ispecie]);\r
509     }\r
510   }\r
511 }\r
512 \r
513 //______________________________________________________________________________\r
514 void AliAnalysisTaskPIDqa::FillTOFqa()\r
515 {\r
516   //\r
517   // Fill TOF information\r
518   //\r
519   AliVEvent *event=InputEvent();\r
520 \r
521   Int_t ntracks=event->GetNumberOfTracks();\r
522   Int_t tracksAtTof = 0;\r
523   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
524     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
525 \r
526     //\r
527     //basic track cuts\r
528     //\r
529     ULong_t status=track->GetStatus();\r
530     // TPC refit + ITS refit +\r
531     // TOF out + TOFpid +\r
532     // kTIME\r
533     // (we don't use kTOFmismatch because it depends on TPC....)\r
534     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
535         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
536         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
537         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
538         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
539 \r
540     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
541     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
542     if (track->GetTPCNclsF()>0) {\r
543       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
544     }\r
545 \r
546     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
547 \r
548     tracksAtTof++;\r
549 \r
550     Double_t mom=track->P();\r
551 \r
552     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
553       TH2 *h=(TH2*)fListQAtof->At(ispecie);\r
554       if (!h) continue;\r
555       Double_t nSigma=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
556       h->Fill(mom,nSigma);\r
557     }\r
558 \r
559     TH2 *h=(TH2*)fListQAtof->FindObject("hSigP_TOF");\r
560     if (h) {\r
561       Double_t sig=track->GetTOFsignal()/1000.;\r
562       h->Fill(mom,sig);\r
563     }\r
564 \r
565     Int_t mask = fPIDResponse->GetTOFResponse().GetStartTimeMask(mom);\r
566     ((TH1F*)fListQAtof->FindObject("hStartTimeMask_TOF"))->Fill((Double_t)(mask+0.5));\r
567 \r
568     if (mom >= 0.75 && mom <= 1.25 ) {\r
569       Double_t nsigma= fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)AliPID::kPion);\r
570       if (mask == 0) {\r
571         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-Fill"))->Fill(nsigma);\r
572       } else if (mask == 1) {\r
573         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-TOF"))->Fill(nsigma);\r
574       } else if ( (mask == 2) || (mask == 4) || (mask == 6) ) {\r
575         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-T0"))->Fill(nsigma);\r
576       } else {\r
577         ((TH1F*)fListQAtof->FindObject("hNsigma_TOF_Pion_T0-Best"))->Fill(nsigma);\r
578       }\r
579     }\r
580 \r
581     Double_t res = (Double_t)fPIDResponse->GetTOFResponse().GetStartTimeRes(mom);\r
582     ((TH1F*)fListQAtof->FindObject("hStartTimeRes_TOF"))->Fill(res);\r
583 \r
584     AliESDEvent *esd = dynamic_cast<AliESDEvent *>(event);\r
585     if (esd) {\r
586       Double_t startTime = esd->GetT0TOF(0);\r
587       if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeAC_T0"))->Fill(startTime);\r
588       else {\r
589         startTime = esd->GetT0TOF(1);\r
590         if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeA_T0"))->Fill(startTime);\r
591         startTime = esd->GetT0TOF(2);\r
592         if (startTime < 90000) ((TH1F*)fListQAtof->FindObject("hStartTimeC_T0"))->Fill(startTime);\r
593       }\r
594     }\r
595   }\r
596   if (tracksAtTof > 0) {\r
597     ((TH1F* )fListQAtof->FindObject("hnTracksAt_TOF"))->Fill(tracksAtTof);\r
598     Int_t mask = fPIDResponse->GetTOFResponse().GetStartTimeMask(5.);\r
599     if (mask & 0x1) ((TH1F*)fListQAtof->FindObject("hT0MakerEff"))->Fill(tracksAtTof);\r
600   }\r
601 }\r
602 \r
603 \r
604 //______________________________________________________________________________\r
605 void AliAnalysisTaskPIDqa::FillEMCALqa()\r
606 {\r
607   //\r
608   // Fill PID qa histograms for the EMCAL\r
609   //\r
610 \r
611   AliVEvent *event=InputEvent();\r
612   \r
613   Int_t ntracks=event->GetNumberOfTracks();\r
614   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
615     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
616     \r
617     //\r
618     //basic track cuts\r
619     //\r
620     ULong_t status=track->GetStatus();\r
621     // not that nice. status bits not in virtual interface\r
622     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
623 \r
624     Double_t pt=track->Pt();\r
625    \r
626     //EMCAL nSigma (only for electrons at the moment)\r
627     TH2 *h=(TH2*)fListQAemcal->At(0);\r
628     if (!h) continue;\r
629     Double_t nSigma=fPIDResponse->NumberOfSigmasEMCAL(track, (AliPID::EParticleType)0);\r
630     h->Fill(pt,nSigma);\r
631     \r
632   }\r
633 \r
634    //EMCAL signal (E/p vs. pT) for electrons from V0\r
635   for(Int_t itrack = 0; itrack < fV0electrons->GetEntries(); itrack++){\r
636     AliVTrack *track=(AliVTrack*)fV0electrons->At(itrack);\r
637 \r
638     //\r
639     //basic track cuts\r
640     //\r
641     ULong_t status=track->GetStatus();\r
642     // not that nice. status bits not in virtual interface\r
643     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
644 \r
645     Double_t pt=track->Pt();\r
646 \r
647     TH2 *h=(TH2*)fListQAemcal->At(1);\r
648     if (h) {\r
649 \r
650       Int_t nMatchClus = track->GetEMCALcluster();\r
651       Double_t mom     = track->P();\r
652       Double_t eop     = -1.;\r
653 \r
654       if(nMatchClus > -1){\r
655     \r
656         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
657 \r
658         if(matchedClus){\r
659 \r
660           // matched cluster is EMCAL\r
661           if(matchedClus->IsEMCAL()){\r
662 \r
663             Double_t fClsE       = matchedClus->E();\r
664             eop                  = fClsE/mom;\r
665 \r
666             h->Fill(pt,eop);\r
667 \r
668           }\r
669         }\r
670       }\r
671     }\r
672   }\r
673 \r
674    //EMCAL signal (E/p vs. pT) for pions from V0\r
675   for(Int_t itrack = 0; itrack < fV0pions->GetEntries(); itrack++){\r
676     AliVTrack *track=(AliVTrack*)fV0pions->At(itrack);\r
677 \r
678     //\r
679     //basic track cuts\r
680     //\r
681     ULong_t status=track->GetStatus();\r
682     // not that nice. status bits not in virtual interface\r
683     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
684 \r
685     Double_t pt=track->Pt();\r
686 \r
687     TH2 *h=(TH2*)fListQAemcal->At(2);\r
688     if (h) {\r
689 \r
690       Int_t nMatchClus = track->GetEMCALcluster();\r
691       Double_t mom     = track->P();\r
692       Double_t eop     = -1.;\r
693 \r
694       if(nMatchClus > -1){\r
695     \r
696         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
697 \r
698         if(matchedClus){\r
699 \r
700           // matched cluster is EMCAL\r
701           if(matchedClus->IsEMCAL()){\r
702 \r
703             Double_t fClsE       = matchedClus->E();\r
704             eop                  = fClsE/mom;\r
705 \r
706             h->Fill(pt,eop);\r
707 \r
708           }\r
709         }\r
710       }\r
711     }\r
712   }\r
713 \r
714    //EMCAL signal (E/p vs. pT) for protons from V0\r
715   for(Int_t itrack = 0; itrack < fV0protons->GetEntries(); itrack++){\r
716     AliVTrack *track=(AliVTrack*)fV0protons->At(itrack);\r
717 \r
718     //\r
719     //basic track cuts\r
720     //\r
721     ULong_t status=track->GetStatus();\r
722     // not that nice. status bits not in virtual interface\r
723     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
724 \r
725     Double_t pt=track->Pt();\r
726 \r
727     TH2 *hP=(TH2*)fListQAemcal->At(3);\r
728     TH2 *hAP=(TH2*)fListQAemcal->At(4);\r
729     if (hP && hAP) {\r
730 \r
731       Int_t nMatchClus = track->GetEMCALcluster();\r
732       Double_t mom     = track->P();\r
733       Int_t charge     = track->Charge();             \r
734       Double_t eop     = -1.;\r
735 \r
736       if(nMatchClus > -1){\r
737     \r
738         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
739 \r
740         if(matchedClus){\r
741 \r
742           // matched cluster is EMCAL\r
743           if(matchedClus->IsEMCAL()){\r
744 \r
745             Double_t fClsE       = matchedClus->E();\r
746             eop                  = fClsE/mom;\r
747 \r
748             if(charge > 0)      hP->Fill(pt,eop);\r
749             else if(charge < 0) hAP->Fill(pt,eop);\r
750 \r
751           }\r
752         }\r
753       }\r
754     }\r
755   }\r
756 \r
757 }\r
758 \r
759 \r
760 //______________________________________________________________________________\r
761 void AliAnalysisTaskPIDqa::FillHMPIDqa()\r
762 {\r
763   //\r
764   // Fill PID qa histograms for the HMPID\r
765   //\r
766   \r
767   AliVEvent *event=InputEvent();\r
768   \r
769   Int_t ntracks=event->GetNumberOfTracks();\r
770   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
771     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
772     \r
773     //\r
774     //basic track cuts\r
775     //\r
776     ULong_t status=track->GetStatus();\r
777     // not that nice. status bits not in virtual interface\r
778     // TPC refit + ITS refit +\r
779     // TOF out + TOFpid +\r
780     // kTIME\r
781     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
782         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ) continue;\r
783 \r
784     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
785     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
786     if (track->GetTPCNclsF()>0) {\r
787       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
788     }\r
789 \r
790     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
791     \r
792     Double_t mom = track->P();\r
793     Double_t ckovAngle = track->GetHMPIDsignal();\r
794     \r
795     TH1F *hThetavsMom = (TH1F*)fListQAhmpid->At(0);;\r
796     \r
797     hThetavsMom->Fill(mom,ckovAngle);    \r
798   \r
799   }\r
800 }\r
801 //______________________________________________________________________________\r
802 void AliAnalysisTaskPIDqa::FillTOFHMPIDqa()\r
803 {\r
804   //\r
805   // Fill PID qa histograms for the HMPID\r
806   //\r
807   \r
808   AliVEvent *event=InputEvent();\r
809   \r
810   Int_t ntracks=event->GetNumberOfTracks();\r
811   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
812     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
813     \r
814     //\r
815     //basic track cuts\r
816     //\r
817     ULong_t status=track->GetStatus();\r
818     // not that nice. status bits not in virtual interface\r
819     // TPC refit + ITS refit +\r
820     // TOF out + TOFpid +\r
821     // kTIME\r
822     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
823         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
824         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
825         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
826         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
827 \r
828     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
829     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
830     if (track->GetTPCNclsF()>0) {\r
831       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
832     }\r
833 \r
834     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
835     \r
836     Double_t mom = track->P();\r
837     Double_t ckovAngle = track->GetHMPIDsignal();\r
838     \r
839     Double_t nSigmaTOF[3]; \r
840     TH1F *h[3];\r
841     \r
842     for (Int_t ispecie=2; ispecie<5; ++ispecie){\r
843       //TOF nSigma\r
844       nSigmaTOF[ispecie-2]=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
845       h[ispecie-2] = (TH1F*)fListQAtofhmpid->At(ispecie-2);}\r
846       \r
847     if(TMath::Abs(nSigmaTOF[0])<2)                                                              h[0]->Fill(mom,ckovAngle);\r
848     \r
849     if(TMath::Abs(nSigmaTOF[1])<2 && TMath::Abs(nSigmaTOF[0])>3)                                h[1]->Fill(mom,ckovAngle);\r
850 \r
851     if(TMath::Abs(nSigmaTOF[2])<2 && TMath::Abs(nSigmaTOF[1])>3 && TMath::Abs(nSigmaTOF[0])>3)  h[2]->Fill(mom,ckovAngle);\r
852       \r
853   }\r
854   \r
855 }\r
856 \r
857 //______________________________________________________________________________\r
858 void AliAnalysisTaskPIDqa::FillTPCTOFqa()\r
859 {\r
860   //\r
861   // Fill PID qa histograms for the TOF\r
862   //   Here also the TPC histograms after TOF selection are filled\r
863   //\r
864 \r
865   AliVEvent *event=InputEvent();\r
866 \r
867   Int_t ntracks=event->GetNumberOfTracks();\r
868   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
869     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
870 \r
871     //\r
872     //basic track cuts\r
873     //\r
874     ULong_t status=track->GetStatus();\r
875     // not that nice. status bits not in virtual interface\r
876     // TPC refit + ITS refit +\r
877     // TOF out + TOFpid +\r
878     // kTIME\r
879     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
880         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
881 //         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid ) || //removes light nuclei, so it is out for the moment\r
882         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
883         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
884         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
885 \r
886     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
887     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
888     if (track->GetTPCNclsF()>0) {\r
889       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
890     }\r
891 \r
892     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
893 \r
894 \r
895     Double_t mom=track->P();\r
896     Double_t momTPC=track->GetTPCmomentum();\r
897 \r
898     for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
899       //TOF nSigma\r
900       Double_t nSigmaTOF=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
901       Double_t nSigmaTPC=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie);\r
902 \r
903       //TPC after TOF cut\r
904       TH2 *h=(TH2*)fListQAtpctof->At(ispecie);\r
905       if (h && TMath::Abs(nSigmaTOF)<3.) h->Fill(momTPC,nSigmaTPC);\r
906 \r
907       //TOF after TPC cut\r
908       h=(TH2*)fListQAtpctof->At(ispecie+AliPID::kSPECIESC);\r
909       if (h && TMath::Abs(nSigmaTPC)<3.) h->Fill(mom,nSigmaTOF);\r
910 \r
911       //EMCAL after TOF and TPC cut\r
912       h=(TH2*)fListQAtpctof->At(ispecie+2*AliPID::kSPECIESC);\r
913       if (h && TMath::Abs(nSigmaTOF)<3. && TMath::Abs(nSigmaTPC)<3. ){\r
914 \r
915         Int_t nMatchClus = track->GetEMCALcluster();\r
916         Double_t pt      = track->Pt();\r
917         Double_t eop     = -1.;\r
918         \r
919         if(nMatchClus > -1){\r
920           \r
921           AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
922           \r
923           if(matchedClus){\r
924             \r
925             // matched cluster is EMCAL\r
926             if(matchedClus->IsEMCAL()){\r
927               \r
928               Double_t fClsE       = matchedClus->E();\r
929               eop                  = fClsE/mom;\r
930 \r
931               h->Fill(pt,eop);\r
932  \r
933               \r
934             }\r
935           }\r
936         }\r
937       }\r
938     }\r
939   }\r
940 }\r
941 \r
942 //______________________________________________________________________________\r
943 void AliAnalysisTaskPIDqa::SetupITSqa()\r
944 {\r
945   //\r
946   // Create the ITS qa objects\r
947   //\r
948   \r
949   TVectorD *vX=MakeLogBinning(200,.1,30);\r
950   \r
951   //ITS+TPC tracks\r
952   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
953     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_ITS_%s",AliPID::ParticleName(ispecie)),\r
954                               Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
955                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
956                               200,-10,10);\r
957     fListQAits->Add(hNsigmaP);\r
958   }\r
959   TH2F *hSig = new TH2F("hSigP_ITS",\r
960                         "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
961                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
962                         300,0,300);\r
963   fListQAits->Add(hSig);\r
964 \r
965   //ITS Standalone tracks\r
966   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
967     TH2F *hNsigmaPSA = new TH2F(Form("hNsigmaP_ITSSA_%s",AliPID::ParticleName(ispecie)),\r
968                                 Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
969                                 vX->GetNrows()-1,vX->GetMatrixArray(),\r
970                                 200,-10,10);\r
971     fListQAitsSA->Add(hNsigmaPSA);\r
972   }\r
973   TH2F *hSigSA = new TH2F("hSigP_ITSSA",\r
974                           "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
975                           vX->GetNrows()-1,vX->GetMatrixArray(),\r
976                           300,0,300);\r
977   fListQAitsSA->Add(hSigSA);\r
978   \r
979   //ITS Pure Standalone tracks\r
980   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
981     TH2F *hNsigmaPPureSA = new TH2F(Form("hNsigmaP_ITSPureSA_%s",AliPID::ParticleName(ispecie)),\r
982                                     Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
983                                     vX->GetNrows()-1,vX->GetMatrixArray(),\r
984                                     200,-10,10);\r
985     fListQAitsPureSA->Add(hNsigmaPPureSA);\r
986   }\r
987   TH2F *hSigPureSA = new TH2F("hSigP_ITSPureSA",\r
988                               "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
989                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
990                               300,0,300);\r
991   fListQAitsPureSA->Add(hSigPureSA);\r
992   \r
993   delete vX;  \r
994 }\r
995 \r
996 //______________________________________________________________________________\r
997 void AliAnalysisTaskPIDqa::SetupTPCqa()\r
998 {\r
999   //\r
1000   // Create the TPC qa objects\r
1001   //\r
1002   \r
1003   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1004   \r
1005   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1006     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_%s",AliPID::ParticleName(ispecie)),\r
1007                               Form("TPC n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1008                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1009                               200,-10,10);\r
1010     fListQAtpc->Add(hNsigmaP);\r
1011   }\r
1012 \r
1013   // the "hybrid" scenario\r
1014   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1015     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_%s_Hybrid",AliPID::ParticleName(ispecie)),\r
1016                               Form("TPC n#sigma %s vs. p (Hybrid gain scenario);p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1017                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1018                               200,-10,10);\r
1019     fListQAtpc->Add(hNsigmaP);\r
1020   }\r
1021    \r
1022   // the "OROC high" scenario\r
1023   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1024     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_%s_OROChigh",AliPID::ParticleName(ispecie)),\r
1025                               Form("TPC n#sigma %s vs. p (OROChigh gain scenario);p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1026                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1027                               200,-10,10);\r
1028     fListQAtpc->Add(hNsigmaP);\r
1029   }\r
1030   \r
1031   \r
1032   \r
1033   TH2F *hSig = new TH2F("hSigP_TPC",\r
1034                         "TPC signal vs. p;p [GeV]; TPC signal [arb. units]",\r
1035                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1036                         300,0,300);\r
1037   fListQAtpc->Add(hSig); //3*AliPID::kSPECIESC\r
1038 \r
1039   delete vX;  \r
1040 }\r
1041 \r
1042 //______________________________________________________________________________\r
1043 void AliAnalysisTaskPIDqa::SetupTRDqa()\r
1044 {\r
1045   //\r
1046   // Create the TRD qa objects\r
1047   //\r
1048   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1049   for(Int_t itl = 0; itl < 6; ++itl){\r
1050     for(Int_t ispecie = 0; ispecie < AliPID::kSPECIES; ispecie++){\r
1051       TH2F *hLikeP = new TH2F(Form("hLikeP_TRD_%dtls_%s", itl, AliPID::ParticleName(ispecie)),\r
1052                               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
1053                               vX->GetNrows()-1, vX->GetMatrixArray(),\r
1054                               100, 0., 1.);\r
1055       fListQAtrd->Add(hLikeP);\r
1056     }\r
1057   }\r
1058   delete vX;\r
1059 }\r
1060 \r
1061 //______________________________________________________________________________\r
1062 void AliAnalysisTaskPIDqa::SetupTOFqa()\r
1063 {\r
1064   //\r
1065   // Create the TOF qa objects\r
1066   //\r
1067   \r
1068   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1069 \r
1070   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1071     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_%s",AliPID::ParticleName(ispecie)),\r
1072                               Form("TOF n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1073                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1074                               200,-10,10);\r
1075     fListQAtof->Add(hNsigmaP);\r
1076   }\r
1077 \r
1078   // for Kaons PID we differentiate on Time Zero\r
1079   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
1080   fListQAtof->Add(hnSigT0Fill);\r
1081   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
1082   fListQAtof->Add(hnSigT0T0);\r
1083   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
1084   fListQAtof->Add(hnSigT0TOF);\r
1085   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
1086   fListQAtof->Add(hnSigT0Best);\r
1087 \r
1088   TH2F *hSig = new TH2F("hSigP_TOF",\r
1089                         "TOF signal vs. p;p [GeV]; TOF signal [ns]",\r
1090                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1091                         300,0,30);\r
1092 \r
1093   delete vX;\r
1094   \r
1095   fListQAtof->Add(hSig);\r
1096 \r
1097   TH1F *hStartTimeMaskTOF = new TH1F("hStartTimeMask_TOF","StartTime mask",8,0,8);\r
1098   fListQAtof->Add(hStartTimeMaskTOF);\r
1099   TH1F *hStartTimeResTOF = new TH1F("hStartTimeRes_TOF","StartTime resolution [ps]",100,0,500);\r
1100   fListQAtof->Add(hStartTimeResTOF);\r
1101 \r
1102   TH1F *hnTracksAtTOF = new TH1F("hnTracksAt_TOF","Matched tracks at TOF",20,0,20);\r
1103   fListQAtof->Add(hnTracksAtTOF);\r
1104   TH1F *hT0MakerEff = new TH1F("hT0MakerEff","Events with T0-TOF vs nTracks",20,0,20);\r
1105   fListQAtof->Add(hT0MakerEff);\r
1106 \r
1107   // this in principle should stay on a T0 PID QA, but are just the data prepared for TOF use\r
1108   TH1F *hStartTimeAT0 = new TH1F("hStartTimeA_T0","StartTime from T0A [ps]",1000,-1000,1000);\r
1109   fListQAtof->Add(hStartTimeAT0);\r
1110   TH1F *hStartTimeCT0 = new TH1F("hStartTimeC_T0","StartTime from T0C [ps]",1000,-1000,1000);\r
1111   fListQAtof->Add(hStartTimeCT0);\r
1112   TH1F *hStartTimeACT0 = new TH1F("hStartTimeAC_T0","StartTime from T0AC [ps]",1000,-1000,1000);;\r
1113   fListQAtof->Add(hStartTimeACT0);\r
1114 }\r
1115 \r
1116 //______________________________________________________________________________\r
1117 void AliAnalysisTaskPIDqa::SetupEMCALqa()\r
1118 {\r
1119   //\r
1120   // Create the EMCAL qa objects\r
1121   //\r
1122 \r
1123   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1124   \r
1125   TH2F *hNsigmaPt = new TH2F(Form("hNsigmaPt_EMCAL_%s",AliPID::ParticleName(0)),\r
1126                              Form("EMCAL n#sigma %s vs. p_{T};p_{T} [GeV]; n#sigma",AliPID::ParticleName(0)),\r
1127                              vX->GetNrows()-1,vX->GetMatrixArray(),\r
1128                              200,-10,10);\r
1129   fListQAemcal->Add(hNsigmaPt);  \r
1130   \r
1131 \r
1132   TH2F *hSigPtEle = new TH2F("hSigPt_EMCAL_Ele",\r
1133                         "EMCAL signal (E/p) vs. p_{T} for electrons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1134                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1135                         200,0,2);\r
1136   fListQAemcal->Add(hSigPtEle);\r
1137 \r
1138   TH2F *hSigPtPions = new TH2F("hSigPt_EMCAL_Pions",\r
1139                         "EMCAL signal (E/p) vs. p_{T} for pions;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1140                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1141                         200,0,2);\r
1142   fListQAemcal->Add(hSigPtPions);\r
1143 \r
1144   TH2F *hSigPtProtons = new TH2F("hSigPt_EMCAL_Protons",\r
1145                         "EMCAL signal (E/p) vs. p_{T} for protons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1146                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1147                         200,0,2);\r
1148   fListQAemcal->Add(hSigPtProtons);\r
1149 \r
1150   TH2F *hSigPtAntiProtons = new TH2F("hSigPt_EMCAL_Antiprotons",\r
1151                         "EMCAL signal (E/p) vs. p_{T} for antiprotons;p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
1152                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
1153                         200,0,2);\r
1154   fListQAemcal->Add(hSigPtAntiProtons);\r
1155 \r
1156   delete vX;  \r
1157 }\r
1158 \r
1159 //______________________________________________________________________________\r
1160 void AliAnalysisTaskPIDqa::SetupHMPIDqa()\r
1161 {\r
1162   //\r
1163   // Create the HMPID qa objects\r
1164   //\r
1165   \r
1166   TH2F *hCkovAnglevsMom   = new TH2F("hCkovAnglevsMom",  "Cherenkov angle vs momnetum",500,0,5.,500,0,1);\r
1167   fListQAhmpid->Add(hCkovAnglevsMom);\r
1168   \r
1169 }\r
1170 \r
1171 //______________________________________________________________________________\r
1172 void AliAnalysisTaskPIDqa::SetupTOFHMPIDqa()\r
1173 {\r
1174   //\r
1175   // Create the HMPID qa objects\r
1176   //\r
1177   \r
1178   TH2F *hCkovAnglevsMomPion   = new TH2F("hCkovAnglevsMom_pion",  "Cherenkov angle vs momnetum for pions",500,0,5.,500,0,1);\r
1179   fListQAtofhmpid->Add(hCkovAnglevsMomPion);\r
1180   \r
1181   TH2F *hCkovAnglevsMomKaon   = new TH2F("hCkovAnglevsMom_kaon",  "Cherenkov angle vs momnetum for kaons",500,0,5.,500,0,1);\r
1182   fListQAtofhmpid->Add(hCkovAnglevsMomKaon);\r
1183   \r
1184   TH2F *hCkovAnglevsMomProton = new TH2F("hCkovAnglevsMom_proton","Cherenkov angle vs momnetum for protons",500,0,5.,500,0,1);\r
1185   fListQAtofhmpid->Add(hCkovAnglevsMomProton);\r
1186   \r
1187   \r
1188 }  \r
1189 \r
1190 //______________________________________________________________________________\r
1191 void AliAnalysisTaskPIDqa::SetupTPCTOFqa()\r
1192 {\r
1193   //\r
1194   // Create the qa objects for TPC + TOF combination\r
1195   //\r
1196   \r
1197   TVectorD *vX=MakeLogBinning(200,.1,30);\r
1198 \r
1199   //TPC signals after TOF cut\r
1200   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1201     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_TOF_%s",AliPID::ParticleName(ispecie)),\r
1202                               Form("TPC n#sigma %s vs. p (after TOF 3#sigma cut);p_{TPC} [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1203                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1204                               200,-10,10);\r
1205     fListQAtpctof->Add(hNsigmaP);\r
1206   }\r
1207 \r
1208   //TOF signals after TPC cut\r
1209   for (Int_t ispecie=0; ispecie<AliPID::kSPECIESC; ++ispecie){\r
1210     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
1211                               Form("TOF n#sigma %s vs. p (after TPC n#sigma cut);p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
1212                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
1213                               200,-10,10);\r
1214     fListQAtpctof->Add(hNsigmaP);\r
1215   }\r
1216 \r
1217   //EMCAL signal after TOF and TPC cut\r
1218   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
1219     TH2F *heopPt = new TH2F(Form("heopPt_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
1220                             Form("EMCAL signal (E/p) %s vs. p_{T};p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",AliPID::ParticleName(ispecie)),\r
1221                             vX->GetNrows()-1,vX->GetMatrixArray(),\r
1222                             200,0,2);\r
1223     fListQAtpctof->Add(heopPt);\r
1224   }\r
1225 \r
1226   delete vX;\r
1227 }\r
1228 //______________________________________________________________________________\r
1229 void AliAnalysisTaskPIDqa::SetupV0qa()\r
1230 {\r
1231   //\r
1232   // Create the qa objects for V0 Kine cuts\r
1233   //\r
1234   \r
1235   TH2F *hArmenteros  = new TH2F("hArmenteros",  "Armenteros plot",200,-1.,1.,200,0.,0.4);\r
1236   fListQAV0->Add(hArmenteros);\r
1237  \r
1238 }\r
1239 \r
1240 //______________________________________________________________________________\r
1241 TVectorD* AliAnalysisTaskPIDqa::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
1242 {\r
1243   //\r
1244   // Make logarithmic binning\r
1245   // the user has to delete the array afterwards!!!\r
1246   //\r
1247   \r
1248   //check limits\r
1249   if (xmin<1e-20 || xmax<1e-20){\r
1250     AliError("For Log binning xmin and xmax must be > 1e-20. Using linear binning instead!");\r
1251     return MakeLinBinning(nbinsX, xmin, xmax);\r
1252   }\r
1253   if (xmax<xmin){\r
1254     Double_t tmp=xmin;\r
1255     xmin=xmax;\r
1256     xmax=tmp;\r
1257   }\r
1258   TVectorD *binLim=new TVectorD(nbinsX+1);\r
1259   Double_t first=xmin;\r
1260   Double_t last=xmax;\r
1261   Double_t expMax=TMath::Log(last/first);\r
1262   for (Int_t i=0; i<nbinsX+1; ++i){\r
1263     (*binLim)[i]=first*TMath::Exp(expMax/nbinsX*(Double_t)i);\r
1264   }\r
1265   return binLim;\r
1266 }\r
1267 \r
1268 //______________________________________________________________________________\r
1269 TVectorD* AliAnalysisTaskPIDqa::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
1270 {\r
1271   //\r
1272   // Make linear binning\r
1273   // the user has to delete the array afterwards!!!\r
1274   //\r
1275   if (xmax<xmin){\r
1276     Double_t tmp=xmin;\r
1277     xmin=xmax;\r
1278     xmax=tmp;\r
1279   }\r
1280   TVectorD *binLim=new TVectorD(nbinsX+1);\r
1281   Double_t first=xmin;\r
1282   Double_t last=xmax;\r
1283   Double_t binWidth=(last-first)/nbinsX;\r
1284   for (Int_t i=0; i<nbinsX+1; ++i){\r
1285     (*binLim)[i]=first+binWidth*(Double_t)i;\r
1286   }\r
1287   return binLim;\r
1288 }\r
1289 \r
1290 //_____________________________________________________________________________\r
1291 TVectorD* AliAnalysisTaskPIDqa::MakeArbitraryBinning(const char* bins)\r
1292 {\r
1293   //\r
1294   // Make arbitrary binning, bins separated by a ','\r
1295   //\r
1296   TString limits(bins);\r
1297   if (limits.IsNull()){\r
1298     AliError("Bin Limit string is empty, cannot add the variable");\r
1299     return 0x0;\r
1300   }\r
1301   \r
1302   TObjArray *arr=limits.Tokenize(",");\r
1303   Int_t nLimits=arr->GetEntries();\r
1304   if (nLimits<2){\r
1305     AliError("Need at leas 2 bin limits, cannot add the variable");\r
1306     delete arr;\r
1307     return 0x0;\r
1308   }\r
1309   \r
1310   TVectorD *binLimits=new TVectorD(nLimits);\r
1311   for (Int_t iLim=0; iLim<nLimits; ++iLim){\r
1312     (*binLimits)[iLim]=(static_cast<TObjString*>(arr->At(iLim)))->GetString().Atof();\r
1313   }\r
1314   \r
1315   delete arr;\r
1316   return binLimits;\r
1317 }\r
1318 \r