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