]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskPIDqa.cxx
o add lists for ITSsa and ITSpureSA
[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 "AliAnalysisTaskPIDqa.h"\r
42 \r
43 \r
44 ClassImp(AliAnalysisTaskPIDqa)\r
45 \r
46 //______________________________________________________________________________\r
47 AliAnalysisTaskPIDqa::AliAnalysisTaskPIDqa():\r
48 AliAnalysisTaskSE(),\r
49 fPIDResponse(0x0),\r
50 fListQA(0x0),\r
51 fListQAits(0x0),\r
52 fListQAitsSA(0x0),\r
53 fListQAitsPureSA(0x0),\r
54 fListQAtpc(0x0),\r
55 fListQAtrd(0x0),\r
56 fListQAtof(0x0),\r
57 fListQAemcal(0x0),\r
58 fListQAtpctof(0x0)\r
59 {\r
60   //\r
61   // Dummy constructor\r
62   //\r
63 }\r
64 \r
65 //______________________________________________________________________________\r
66 AliAnalysisTaskPIDqa::AliAnalysisTaskPIDqa(const char* name):\r
67 AliAnalysisTaskSE(name),\r
68 fPIDResponse(0x0),\r
69 fListQA(0x0),\r
70 fListQAits(0x0),\r
71 fListQAitsSA(0x0),\r
72 fListQAitsPureSA(0x0),\r
73 fListQAtpc(0x0),\r
74 fListQAtrd(0x0),\r
75 fListQAtof(0x0),\r
76 fListQAemcal(0x0),\r
77 fListQAtpctof(0x0)\r
78 {\r
79   //\r
80   // Default constructor\r
81   //\r
82   DefineInput(0,TChain::Class());\r
83   DefineOutput(1,TList::Class());\r
84 }\r
85 \r
86 //______________________________________________________________________________\r
87 AliAnalysisTaskPIDqa::~AliAnalysisTaskPIDqa()\r
88 {\r
89   //\r
90   // Destructor\r
91   //\r
92 \r
93 }\r
94 \r
95 //______________________________________________________________________________\r
96 void AliAnalysisTaskPIDqa::UserCreateOutputObjects()\r
97 {\r
98   //\r
99   // Create the output QA objects\r
100   //\r
101 \r
102   AliLog::SetClassDebugLevel("AliAnalysisTaskPIDqa",10);\r
103 \r
104   //input hander\r
105   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();\r
106   AliInputEventHandler *inputHandler=dynamic_cast<AliInputEventHandler*>(man->GetInputEventHandler());\r
107   if (!inputHandler) AliFatal("Input handler needed");\r
108 \r
109   //pid response object\r
110   fPIDResponse=inputHandler->GetPIDResponse();\r
111   if (!fPIDResponse) AliError("PIDResponse object was not created");\r
112   \r
113   //\r
114   fListQA=new TList;\r
115   fListQA->SetOwner();\r
116   \r
117   fListQAits=new TList;\r
118   fListQAits->SetOwner();\r
119   fListQAits->SetName("ITS");\r
120 \r
121   fListQAitsSA=new TList;\r
122   fListQAitsSA->SetOwner();\r
123   fListQAitsSA->SetName("ITS");\r
124 \r
125   fListQAitsPureSA=new TList;\r
126   fListQAitsPureSA->SetOwner();\r
127   fListQAitsPureSA->SetName("ITS");\r
128 \r
129   fListQAtpc=new TList;\r
130   fListQAtpc->SetOwner();\r
131   fListQAtpc->SetName("TPC");\r
132   \r
133   fListQAtrd=new TList;\r
134   fListQAtrd->SetOwner();\r
135   fListQAtrd->SetName("TRD");\r
136   \r
137   fListQAtof=new TList;\r
138   fListQAtof->SetOwner();\r
139   fListQAtof->SetName("TOF");\r
140   \r
141   fListQAemcal=new TList;\r
142   fListQAemcal->SetOwner();\r
143   fListQAemcal->SetName("EMCAL");\r
144 \r
145   fListQAtpctof=new TList;\r
146   fListQAtpctof->SetOwner();\r
147   fListQAtpctof->SetName("TPC_TOF");\r
148 \r
149   fListQA->Add(fListQAits);\r
150   fListQA->Add(fListQAitsSA);\r
151   fListQA->Add(fListQAitsPureSA);\r
152   fListQA->Add(fListQAtpc);\r
153   fListQA->Add(fListQAtrd);\r
154   fListQA->Add(fListQAtof);\r
155   fListQA->Add(fListQAemcal);\r
156   fListQA->Add(fListQAtpctof);\r
157 \r
158   SetupITSqa();\r
159   SetupTPCqa();\r
160   SetupTRDqa();\r
161   SetupTOFqa();\r
162   SetupEMCALqa();\r
163   SetupTPCTOFqa();\r
164 \r
165   PostData(1,fListQA);\r
166 }\r
167 \r
168 \r
169 //______________________________________________________________________________\r
170 void AliAnalysisTaskPIDqa::UserExec(Option_t */*option*/)\r
171 {\r
172   //\r
173   // Setup the PID response functions and fill the QA histograms\r
174   //\r
175 \r
176   AliVEvent *event=InputEvent();\r
177   if (!event||!fPIDResponse) return;\r
178 \r
179   \r
180   FillITSqa();\r
181   FillTPCqa();\r
182   FillTRDqa();\r
183   FillTOFqa();\r
184   FillEMCALqa();\r
185   FillTPCTOFqa();\r
186 \r
187   PostData(1,fListQA);\r
188 }\r
189 \r
190 //______________________________________________________________________________\r
191 void AliAnalysisTaskPIDqa::FillITSqa()\r
192 {\r
193   //\r
194   // Fill PID qa histograms for the ITS\r
195   //\r
196 \r
197   AliVEvent *event=InputEvent();\r
198   \r
199   Int_t ntracks=event->GetNumberOfTracks();\r
200   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
201     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
202     ULong_t status=track->GetStatus();\r
203     // not that nice. status bits not in virtual interface\r
204     // ITS refit + ITS pid\r
205     if (!( ( (status & AliVTrack::kITSrefit)==AliVTrack::kITSrefit ) &&\r
206            ( (status & AliVTrack::kITSpid  )==AliVTrack::kITSpid   ) )) continue;\r
207     Double_t mom=track->P();\r
208     \r
209     for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
210       TH2 *h=(TH2*)fListQAits->At(ispecie);\r
211       if (!h) continue;\r
212       Double_t nSigma=fPIDResponse->NumberOfSigmasITS(track, (AliPID::EParticleType)ispecie);\r
213       h->Fill(mom,nSigma);\r
214     }\r
215     \r
216     TH2 *h=(TH2*)fListQAits->At(AliPID::kSPECIES);\r
217     if (h) {\r
218       Double_t sig=track->GetITSsignal();\r
219       h->Fill(mom,sig);\r
220     }\r
221     \r
222   }\r
223 }\r
224 \r
225 //______________________________________________________________________________\r
226 void AliAnalysisTaskPIDqa::FillTPCqa()\r
227 {\r
228   //\r
229   // Fill PID qa histograms for the TPC\r
230   //\r
231   \r
232   AliVEvent *event=InputEvent();\r
233   \r
234   Int_t ntracks=event->GetNumberOfTracks();\r
235   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
236     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
237     \r
238     //\r
239     //basic track cuts\r
240     //\r
241     ULong_t status=track->GetStatus();\r
242     // not that nice. status bits not in virtual interface\r
243     // TPC refit + ITS refit + TPC pid\r
244     if (!( (status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
245         !( (status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
246         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid  ) ) continue;\r
247     \r
248     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
249     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
250     if (track->GetTPCNclsF()>0) {\r
251       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
252     }\r
253     \r
254     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
255     \r
256     Double_t mom=track->GetTPCmomentum();\r
257     \r
258     for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
259       TH2 *h=(TH2*)fListQAtpc->At(ispecie);\r
260       if (!h) continue;\r
261       Double_t nSigma=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie);\r
262       h->Fill(mom,nSigma);\r
263     }\r
264     \r
265     TH2 *h=(TH2*)fListQAtpc->At(AliPID::kSPECIES);\r
266     if (h) {\r
267       Double_t sig=track->GetTPCsignal();\r
268       h->Fill(mom,sig);\r
269     }\r
270   }\r
271 }\r
272 \r
273 //______________________________________________________________________________\r
274 void AliAnalysisTaskPIDqa::FillTRDqa()\r
275 {\r
276   //\r
277   // Fill PID qa histograms for the TRD\r
278   //\r
279 \r
280 }\r
281 \r
282 //______________________________________________________________________________\r
283 void AliAnalysisTaskPIDqa::FillTOFqa()\r
284 {\r
285   //\r
286   // Fill PID qa histograms for the TOF\r
287   //\r
288 \r
289   AliVEvent *event=InputEvent();\r
290   \r
291   Int_t ntracks=event->GetNumberOfTracks();\r
292   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
293     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
294     \r
295     //\r
296     //basic track cuts\r
297     //\r
298     ULong_t status=track->GetStatus();\r
299     // not that nice. status bits not in virtual interface\r
300     // TPC refit + ITS refit +\r
301     // TOF out + TOFpid +\r
302     // kTIME\r
303     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
304         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
305         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
306         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
307         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
308     \r
309     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
310     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
311     if (track->GetTPCNclsF()>0) {\r
312       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
313     }\r
314     \r
315     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
316     \r
317     \r
318     Double_t mom=track->P();\r
319     \r
320     for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
321       //TOF nSigma\r
322       TH2 *h=(TH2*)fListQAtof->At(ispecie);\r
323       if (!h) continue;\r
324       Double_t nSigma=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
325       h->Fill(mom,nSigma);\r
326     }\r
327     \r
328     TH2 *h=(TH2*)fListQAtof->At(AliPID::kSPECIES);\r
329     if (h) {\r
330       Double_t sig=track->GetTOFsignal();\r
331       h->Fill(mom,sig);\r
332     }\r
333     \r
334   }\r
335 }\r
336 \r
337 //______________________________________________________________________________\r
338 void AliAnalysisTaskPIDqa::FillEMCALqa()\r
339 {\r
340   //\r
341   // Fill PID qa histograms for the EMCAL\r
342   //\r
343 \r
344   AliVEvent *event=InputEvent();\r
345   \r
346   Int_t ntracks=event->GetNumberOfTracks();\r
347   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
348     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
349     \r
350     //\r
351     //basic track cuts\r
352     //\r
353     ULong_t status=track->GetStatus();\r
354     // not that nice. status bits not in virtual interface\r
355     // TPC refit + ITS refit +\r
356     // TOF out + TOFpid +\r
357     // kTIME\r
358     if (!( (status & AliVTrack::kEMCALmatch) == AliVTrack::kEMCALmatch) ) continue;\r
359 \r
360     Double_t pt=track->Pt();\r
361    \r
362     //EMCAL nSigma (only for electrons at the moment)\r
363     TH2 *h=(TH2*)fListQAemcal->At(0);\r
364     if (!h) continue;\r
365     Double_t nSigma=fPIDResponse->NumberOfSigmasEMCAL(track, (AliPID::EParticleType)0);\r
366     h->Fill(pt,nSigma);\r
367     \r
368     //EMCAL signal (E/p vs. pT)\r
369     h=(TH2*)fListQAemcal->At(1);\r
370     if (h) {\r
371 \r
372       Int_t nMatchClus = track->GetEMCALcluster();\r
373       Double_t mom     = track->P();\r
374       Double_t eop     = -1.;\r
375 \r
376       if(nMatchClus > -1){\r
377     \r
378         AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
379     \r
380         if(matchedClus){\r
381 \r
382           // matched cluster is EMCAL\r
383           if(matchedClus->IsEMCAL()){\r
384       \r
385             Double_t fClsE       = matchedClus->E();\r
386             eop                  = fClsE/mom;\r
387 \r
388             h->Fill(pt,eop);\r
389             \r
390           }\r
391         }\r
392       }\r
393       else{\r
394         Printf("status status = AliVTrack::kEMCALmatch, BUT no matched cluster!");\r
395       }\r
396     }\r
397   }\r
398 }\r
399 \r
400 //______________________________________________________________________________\r
401 void AliAnalysisTaskPIDqa::FillTPCTOFqa()\r
402 {\r
403   //\r
404   // Fill PID qa histograms for the TOF\r
405   //   Here also the TPC histograms after TOF selection are filled\r
406   //\r
407 \r
408   AliVEvent *event=InputEvent();\r
409 \r
410   Int_t ntracks=event->GetNumberOfTracks();\r
411   for(Int_t itrack = 0; itrack < ntracks; itrack++){\r
412     AliVTrack *track=(AliVTrack*)event->GetTrack(itrack);\r
413 \r
414     //\r
415     //basic track cuts\r
416     //\r
417     ULong_t status=track->GetStatus();\r
418     // not that nice. status bits not in virtual interface\r
419     // TPC refit + ITS refit +\r
420     // TOF out + TOFpid +\r
421     // kTIME\r
422     if (!((status & AliVTrack::kTPCrefit) == AliVTrack::kTPCrefit) ||\r
423         !((status & AliVTrack::kITSrefit) == AliVTrack::kITSrefit) ||\r
424         !( (status & AliVTrack::kTPCpid  ) == AliVTrack::kTPCpid ) ||\r
425         !((status & AliVTrack::kTOFout  ) == AliVTrack::kTOFout  ) ||\r
426         !((status & AliVTrack::kTOFpid  ) == AliVTrack::kTOFpid  ) ||\r
427         !((status & AliVTrack::kTIME    ) == AliVTrack::kTIME    ) ) continue;\r
428 \r
429     Float_t nCrossedRowsTPC = track->GetTPCClusterInfo(2,1);\r
430     Float_t  ratioCrossedRowsOverFindableClustersTPC = 1.0;\r
431     if (track->GetTPCNclsF()>0) {\r
432       ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC/track->GetTPCNclsF();\r
433     }\r
434 \r
435     if ( nCrossedRowsTPC<70 || ratioCrossedRowsOverFindableClustersTPC<.8 ) continue;\r
436 \r
437 \r
438     Double_t mom=track->P();\r
439     Double_t momTPC=track->GetTPCmomentum();\r
440 \r
441     for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
442       //TOF nSigma\r
443       Double_t nSigmaTOF=fPIDResponse->NumberOfSigmasTOF(track, (AliPID::EParticleType)ispecie);\r
444       Double_t nSigmaTPC=fPIDResponse->NumberOfSigmasTPC(track, (AliPID::EParticleType)ispecie);\r
445 \r
446       //TPC after TOF cut\r
447       TH2 *h=(TH2*)fListQAtpctof->At(ispecie);\r
448       if (h && TMath::Abs(nSigmaTOF)<3.) h->Fill(momTPC,nSigmaTPC);\r
449 \r
450       //TOF after TPC cut\r
451       h=(TH2*)fListQAtpctof->At(ispecie+AliPID::kSPECIES);\r
452       if (h && TMath::Abs(nSigmaTPC)<3.) h->Fill(mom,nSigmaTOF);\r
453 \r
454       //EMCAL after TOF and TPC cut\r
455       h=(TH2*)fListQAtpctof->At(ispecie+2*AliPID::kSPECIES);\r
456       if (h && TMath::Abs(nSigmaTOF)<3. && TMath::Abs(nSigmaTPC)<3. ){\r
457 \r
458         Int_t nMatchClus = track->GetEMCALcluster();\r
459         Double_t pt      = track->Pt();\r
460         Double_t eop     = -1.;\r
461         \r
462         if(nMatchClus > -1){\r
463           \r
464           AliVCluster *matchedClus = (AliVCluster*)event->GetCaloCluster(nMatchClus);\r
465           \r
466           if(matchedClus){\r
467             \r
468             // matched cluster is EMCAL\r
469             if(matchedClus->IsEMCAL()){\r
470               \r
471               Double_t fClsE       = matchedClus->E();\r
472               eop                  = fClsE/mom;\r
473 \r
474               h->Fill(pt,eop);\r
475  \r
476               \r
477             }\r
478           }\r
479         }\r
480       }\r
481     }\r
482   }\r
483 }\r
484 \r
485 //______________________________________________________________________________\r
486 void AliAnalysisTaskPIDqa::SetupITSqa()\r
487 {\r
488   //\r
489   // Create the ITS qa objects\r
490   //\r
491   \r
492   TVectorD *vX=MakeLogBinning(200,.1,30);\r
493   \r
494   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
495     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_ITS_%s",AliPID::ParticleName(ispecie)),\r
496                               Form("ITS n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
497                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
498                               200,-10,10);\r
499     fListQAits->Add(hNsigmaP);\r
500   }\r
501   \r
502   \r
503   TH2F *hSig = new TH2F("hSigP_ITS",\r
504                         "ITS signal vs. p;p [GeV]; ITS signal [arb. units]",\r
505                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
506                         300,0,300);\r
507   fListQAits->Add(hSig);\r
508 \r
509   delete vX;  \r
510 }\r
511 \r
512 //______________________________________________________________________________\r
513 void AliAnalysisTaskPIDqa::SetupTPCqa()\r
514 {\r
515   //\r
516   // Create the TPC qa objects\r
517   //\r
518   \r
519   TVectorD *vX=MakeLogBinning(200,.1,30);\r
520   \r
521   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
522     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_%s",AliPID::ParticleName(ispecie)),\r
523                               Form("TPC n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
524                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
525                               200,-10,10);\r
526     fListQAtpc->Add(hNsigmaP);\r
527   }\r
528   \r
529   \r
530   TH2F *hSig = new TH2F("hSigP_TPC",\r
531                         "TPC signal vs. p;p [GeV]; TPC signal [arb. units]",\r
532                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
533                         300,0,300);\r
534   fListQAtpc->Add(hSig);\r
535 \r
536   delete vX;  \r
537 }\r
538 \r
539 //______________________________________________________________________________\r
540 void AliAnalysisTaskPIDqa::SetupTRDqa()\r
541 {\r
542   //\r
543   // Create the TRD qa objects\r
544   //\r
545   \r
546 }\r
547 \r
548 //______________________________________________________________________________\r
549 void AliAnalysisTaskPIDqa::SetupTOFqa()\r
550 {\r
551   //\r
552   // Create the TOF qa objects\r
553   //\r
554   \r
555   TVectorD *vX=MakeLogBinning(200,.1,30);\r
556   \r
557   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
558     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_%s",AliPID::ParticleName(ispecie)),\r
559                               Form("TOF n#sigma %s vs. p;p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
560                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
561                               200,-10,10);\r
562     fListQAtof->Add(hNsigmaP);\r
563   }\r
564   \r
565   \r
566   TH2F *hSig = new TH2F("hSigP_TOF",\r
567                         "TOF signal vs. p;p [GeV]; TOF signal [arb. units]",\r
568                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
569                         300,0,300);\r
570   fListQAtof->Add(hSig);\r
571 \r
572   delete vX;  \r
573 }\r
574 \r
575 //______________________________________________________________________________\r
576 void AliAnalysisTaskPIDqa::SetupEMCALqa()\r
577 {\r
578   //\r
579   // Create the EMCAL qa objects\r
580   //\r
581 \r
582   TVectorD *vX=MakeLogBinning(200,.1,30);\r
583   \r
584   TH2F *hNsigmaPt = new TH2F(Form("hNsigmaPt_EMCAL_%s",AliPID::ParticleName(0)),\r
585                              Form("EMCAL n#sigma %s vs. p_{T};p_{T} [GeV]; n#sigma",AliPID::ParticleName(0)),\r
586                              vX->GetNrows()-1,vX->GetMatrixArray(),\r
587                              200,-10,10);\r
588   fListQAemcal->Add(hNsigmaPt);  \r
589   \r
590   TH2F *hSigPt = new TH2F("hSigPt_EMCAL",\r
591                         "EMCAL signal (E/p) vs. p_{T};p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",\r
592                         vX->GetNrows()-1,vX->GetMatrixArray(),\r
593                         200,0,2);\r
594   fListQAemcal->Add(hSigPt);\r
595 \r
596   delete vX;  \r
597 }\r
598 \r
599 //______________________________________________________________________________\r
600 void AliAnalysisTaskPIDqa::SetupTPCTOFqa()\r
601 {\r
602   //\r
603   // Create the qa objects for TPC + TOF combination\r
604   //\r
605   \r
606   TVectorD *vX=MakeLogBinning(200,.1,30);\r
607 \r
608   //TPC signals after TOF cut\r
609   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
610     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TPC_TOF_%s",AliPID::ParticleName(ispecie)),\r
611                               Form("TPC n#sigma %s vs. p (after TOF 3#sigma cut);p_{TPC} [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
612                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
613                               200,-10,10);\r
614     fListQAtpctof->Add(hNsigmaP);\r
615   }\r
616 \r
617   //TOF signals after TPC cut\r
618   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
619     TH2F *hNsigmaP = new TH2F(Form("hNsigmaP_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
620                               Form("TOF n#sigma %s vs. p (after TPC n#sigma cut);p [GeV]; n#sigma",AliPID::ParticleName(ispecie)),\r
621                               vX->GetNrows()-1,vX->GetMatrixArray(),\r
622                               200,-10,10);\r
623     fListQAtpctof->Add(hNsigmaP);\r
624   }\r
625 \r
626   //EMCAL signal after TOF and TPC cut\r
627   for (Int_t ispecie=0; ispecie<AliPID::kSPECIES; ++ispecie){\r
628     TH2F *heopPt = new TH2F(Form("heopPt_TOF_TPC_%s",AliPID::ParticleName(ispecie)),\r
629                             Form("EMCAL signal (E/p) %s vs. p_{T};p_{T} [GeV]; EMCAL signal (E/p) [arb. units]",AliPID::ParticleName(ispecie)),\r
630                             vX->GetNrows()-1,vX->GetMatrixArray(),\r
631                             200,0,2);\r
632     fListQAtpctof->Add(heopPt);\r
633   }\r
634 \r
635   delete vX;\r
636 }\r
637 \r
638 //______________________________________________________________________________\r
639 TVectorD* AliAnalysisTaskPIDqa::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
640 {\r
641   //\r
642   // Make logarithmic binning\r
643   // the user has to delete the array afterwards!!!\r
644   //\r
645   \r
646   //check limits\r
647   if (xmin<1e-20 || xmax<1e-20){\r
648     AliError("For Log binning xmin and xmax must be > 1e-20. Using linear binning instead!");\r
649     return MakeLinBinning(nbinsX, xmin, xmax);\r
650   }\r
651   if (xmax<xmin){\r
652     Double_t tmp=xmin;\r
653     xmin=xmax;\r
654     xmax=tmp;\r
655   }\r
656   TVectorD *binLim=new TVectorD(nbinsX+1);\r
657   Double_t first=xmin;\r
658   Double_t last=xmax;\r
659   Double_t expMax=TMath::Log(last/first);\r
660   for (Int_t i=0; i<nbinsX+1; ++i){\r
661     (*binLim)[i]=first*TMath::Exp(expMax/nbinsX*(Double_t)i);\r
662   }\r
663   return binLim;\r
664 }\r
665 \r
666 //______________________________________________________________________________\r
667 TVectorD* AliAnalysisTaskPIDqa::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)\r
668 {\r
669   //\r
670   // Make linear binning\r
671   // the user has to delete the array afterwards!!!\r
672   //\r
673   if (xmax<xmin){\r
674     Double_t tmp=xmin;\r
675     xmin=xmax;\r
676     xmax=tmp;\r
677   }\r
678   TVectorD *binLim=new TVectorD(nbinsX+1);\r
679   Double_t first=xmin;\r
680   Double_t last=xmax;\r
681   Double_t binWidth=(last-first)/nbinsX;\r
682   for (Int_t i=0; i<nbinsX+1; ++i){\r
683     (*binLim)[i]=first+binWidth*(Double_t)i;\r
684   }\r
685   return binLim;\r
686 }\r
687 \r
688 //_____________________________________________________________________________\r
689 TVectorD* AliAnalysisTaskPIDqa::MakeArbitraryBinning(const char* bins)\r
690 {\r
691   //\r
692   // Make arbitrary binning, bins separated by a ','\r
693   //\r
694   TString limits(bins);\r
695   if (limits.IsNull()){\r
696     AliError("Bin Limit string is empty, cannot add the variable");\r
697     return 0x0;\r
698   }\r
699   \r
700   TObjArray *arr=limits.Tokenize(",");\r
701   Int_t nLimits=arr->GetEntries();\r
702   if (nLimits<2){\r
703     AliError("Need at leas 2 bin limits, cannot add the variable");\r
704     delete arr;\r
705     return 0x0;\r
706   }\r
707   \r
708   TVectorD *binLimits=new TVectorD(nLimits);\r
709   for (Int_t iLim=0; iLim<nLimits; ++iLim){\r
710     (*binLimits)[iLim]=(static_cast<TObjString*>(arr->At(iLim)))->GetString().Atof();\r
711   }\r
712   \r
713   delete arr;\r
714   return binLimits;\r
715 }\r