]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TOF/AliAnalysisTaskTOFqaID.cxx
Reverting previous commit and going back to rev.64119Subversion Repositories:
[u/mrichter/AliRoot.git] / PWGPP / TOF / AliAnalysisTaskTOFqaID.cxx
1 /*  created by fbellini@cern.ch on 29/04/2013 */
2 /*  last modified by fbellini   on 19/08/2013 */
3
4
5 #ifndef ALIANALYSISTASKTOFQAID_CXX
6 #define ALIANALYSISTASKTOFQAID_CXX
7
8 #include "TChain.h"
9 #include "TTree.h"
10 #include "TH1F.h"
11 #include "TH2F.h"
12 #include "TProfile.h"
13 #include "TCanvas.h"
14 #include "AliAnalysisTaskSE.h"
15 #include "AliAnalysisManager.h"
16 #include "AliVEvent.h"
17 #include "AliVTrack.h"
18 #include "AliESDEvent.h"
19 #include "AliMCEvent.h"
20 #include "AliMCParticle.h"
21 #include "AliESDInputHandler.h"
22 #include "AliMCEventHandler.h"
23 #include "AliESDpid.h"
24 #include "AliCDBManager.h"
25 #include "AliTOFcalib.h"
26 #include "AliTOFT0maker.h"
27 #include "AliTOFT0v1.h"
28 #include "AliAnalysisTaskTOFqaID.h"
29 #include "AliAnalysisFilter.h"
30 #include "AliESDtrackCuts.h"
31 #include "AliLog.h"
32 #include "AliTOFRawStream.h"
33 #include "AliTOFGeometry.h"
34
35 ClassImp(AliAnalysisTaskTOFqaID)
36
37 //________________________________________________________________________
38 AliAnalysisTaskTOFqaID::AliAnalysisTaskTOFqaID() :
39   fRunNumber(0), 
40   fESD(0x0), 
41   fMCevent(0x0),
42   fTrackFilter(0x0), 
43   fVertex(0x0),
44   fESDpid(new AliESDpid()),
45   fTOFT0v1(new AliTOFT0v1(fESDpid)),
46   fTOFHeader(0x0),
47   fEnableAdvancedCheck(kFALSE),
48   fExpTimeBinWidth(24.4),
49   fExpTimeRangeMin(-25010.),
50   fExpTimeRangeMax(25010.),
51   fExpTimeSmallRangeMin(-5002.),
52   fExpTimeSmallRangeMax(5002.),
53   fnExpTimeBins(1),
54   fnExpTimeSmallBins(1),  
55   fMyTimeZeroTOF(1e20),
56   fMyTimeZeroTOFsigma(1e20),
57   fMyTimeZeroTOFtracks(-1),
58   fIsMC(kFALSE),
59   fSelectedPdg(0),
60   fP(1e10),
61   fPt(1e10),
62   fEta(1e10),
63   fPhi(1e10),
64   fTPCOuterPhi(1e10),
65   fL(1e10),
66   fMatchingMomCut(1e10),
67   fTof(1e10),
68   fHlist(0x0),
69   fHlistTimeZero(0x0),
70   fHlistPID(0x0),
71   fHlistTRD(0x0),
72   fHlistTrigger(0x0)
73 {
74   // Default constructor
75    
76    for (Int_t j=0;j<3;j++ ) {
77      if (j<3) {
78        fT0[j]=0.0;
79        fNTOFtracks[j]=0;
80      }
81      fSigmaSpecie[j]=0.0;
82      fTrkExpTimes[j]=0.0;
83      fThExpTimes[j]=0.0;
84    }
85 }
86 //________________________________________________________________________
87 AliAnalysisTaskTOFqaID::AliAnalysisTaskTOFqaID(const char *name) : 
88   AliAnalysisTaskSE(name), 
89   fRunNumber(0), 
90   fESD(0x0), 
91   fMCevent(0x0),
92   fTrackFilter(0x0),
93   fVertex(0x0),
94   fESDpid(new AliESDpid()),
95   fTOFT0v1(new AliTOFT0v1(fESDpid)),
96   fTOFHeader(0x0),
97   fEnableAdvancedCheck(kFALSE),
98   fExpTimeBinWidth(24.4),
99   fExpTimeRangeMin(-25010.),
100   fExpTimeRangeMax(25010.),
101   fExpTimeSmallRangeMin(-5002.),
102   fExpTimeSmallRangeMax(5002.),
103   fnExpTimeBins(1),
104   fnExpTimeSmallBins(1),
105   fMyTimeZeroTOF(1e20),
106   fMyTimeZeroTOFsigma(1e20),
107   fMyTimeZeroTOFtracks(-1),
108   fIsMC(kFALSE),
109   fSelectedPdg(0),
110   fP(1e10),
111   fPt(1e10),
112   fEta(1e10),
113   fPhi(1e10),
114   fTPCOuterPhi(1e10),
115   fL(1e10),
116   fMatchingMomCut(1.0),
117   fTof(1e10),
118   fHlist(0x0),
119   fHlistTimeZero(0x0),
120   fHlistPID(0x0),
121   fHlistTRD(0x0),
122   fHlistTrigger(0x0)
123  {
124   // Constructor
125   // Define input and output slots here
126    Info("AliAnalysisTaskTOFqaID","Calling Constructor");
127    
128    for (Int_t j=0;j<5;j++ ) {
129      if (j<3){ 
130        fT0[j]=0.0;
131        fNTOFtracks[j]=0;
132      }
133      fSigmaSpecie[j]=0.0;
134      fTrkExpTimes[j]=0.0;
135      fThExpTimes[j]=0.0;
136    }
137    // Input slot #0 works with a TChain
138    DefineInput(0, TChain::Class());
139    
140    // Output slot #0 writes into a TH1 container
141    // Output slot #1 writes into a user defined  container
142    DefineOutput(1, TList::Class());
143    DefineOutput(2, TList::Class());
144    DefineOutput(3, TList::Class());
145    DefineOutput(4, TList::Class());
146    DefineOutput(5, TList::Class());
147    
148  }
149
150 //________________________________________________________________________
151 AliAnalysisTaskTOFqaID::AliAnalysisTaskTOFqaID(const AliAnalysisTaskTOFqaID& copy) 
152 : AliAnalysisTaskSE(), 
153   fRunNumber(copy.fRunNumber), 
154   fESD(copy.fESD), 
155   fMCevent(copy.fMCevent),
156   fTrackFilter(copy.fTrackFilter), 
157   fVertex(copy.fVertex),
158   fESDpid(copy.fESDpid),
159   fTOFT0v1(copy.fTOFT0v1),
160   fTOFHeader(copy.fTOFHeader),
161   fEnableAdvancedCheck(copy.fEnableAdvancedCheck),
162   fExpTimeBinWidth(copy.fExpTimeBinWidth),
163   fExpTimeRangeMin(copy.fExpTimeRangeMin),
164   fExpTimeRangeMax(copy.fExpTimeRangeMax),
165   fExpTimeSmallRangeMin(copy.fExpTimeSmallRangeMin),
166   fExpTimeSmallRangeMax(copy.fExpTimeSmallRangeMax),
167   fnExpTimeBins(copy.fnExpTimeBins),
168   fnExpTimeSmallBins(copy.fnExpTimeSmallBins),
169   fMyTimeZeroTOF(copy.fMyTimeZeroTOF),
170   fMyTimeZeroTOFsigma(copy.fMyTimeZeroTOFsigma),
171   fMyTimeZeroTOFtracks(copy.fMyTimeZeroTOFtracks),
172   fIsMC(copy.fIsMC),
173   fSelectedPdg(copy.fSelectedPdg),
174   fP(copy.fP),
175   fPt(copy.fPt),
176   fEta(copy.fEta),
177   fPhi(copy.fPhi),
178   fTPCOuterPhi(copy.fTPCOuterPhi),
179   fL(copy.fL),
180   fMatchingMomCut(copy.fMatchingMomCut),
181   fTof(copy.fTof),
182   fHlist(copy.fHlist),
183   fHlistTimeZero(copy.fHlistTimeZero),
184   fHlistPID(copy.fHlistPID),
185   fHlistTRD(copy.fHlistTRD),
186   fHlistTrigger(copy.fHlistTrigger)
187 {
188   // Copy constructor
189    for (Int_t j=0;j<5;j++ ) {
190      if (j<3) { 
191        fT0[j]=copy.fT0[j];
192        fNTOFtracks[j]=copy.fNTOFtracks[j];
193      }
194      fSigmaSpecie[j]=copy.fSigmaSpecie[j];
195      fTrkExpTimes[j]=copy.fTrkExpTimes[j];
196      fThExpTimes[j]=copy.fThExpTimes[j];
197    }
198   
199
200 }
201
202 //___________________________________________________________________________
203 AliAnalysisTaskTOFqaID& AliAnalysisTaskTOFqaID::operator=(const AliAnalysisTaskTOFqaID& copy) 
204 {
205   //
206   // Assignment operator
207   //
208   if (this!=&copy) {
209     AliAnalysisTaskSE::operator=(copy) ;
210     fRunNumber=copy.fRunNumber; 
211     fESD=copy.fESD;
212     fMCevent=copy.fMCevent;
213     fTrackFilter=copy.fTrackFilter;
214     fVertex=copy.fVertex;
215     fESDpid=copy.fESDpid;
216     fTOFT0v1=copy.fTOFT0v1;
217     fTOFHeader=copy.fTOFHeader;
218     fEnableAdvancedCheck=copy.fEnableAdvancedCheck;
219     fExpTimeBinWidth=copy.fExpTimeBinWidth;
220     fExpTimeRangeMin=copy.fExpTimeRangeMin;
221     fExpTimeRangeMax=copy.fExpTimeRangeMax;
222     fExpTimeSmallRangeMin=copy.fExpTimeSmallRangeMin;
223     fExpTimeSmallRangeMax=copy.fExpTimeSmallRangeMax;
224     fnExpTimeBins=copy.fnExpTimeBins;
225     fnExpTimeSmallBins=copy.fnExpTimeSmallBins;
226     fMyTimeZeroTOF=copy.fMyTimeZeroTOF;
227     fMyTimeZeroTOFsigma=copy.fMyTimeZeroTOFsigma;
228     fMyTimeZeroTOFtracks=copy.fMyTimeZeroTOFtracks;
229     for (Int_t j=0;j<5;j++ ) {
230       if (j<3) {
231         fT0[j]=copy.fT0[j];
232         fNTOFtracks[j]=copy.fNTOFtracks[j];
233       }
234       fSigmaSpecie[j]=copy.fSigmaSpecie[j];
235       fTrkExpTimes[j]=copy.fTrkExpTimes[j];
236       fThExpTimes[j]=copy.fThExpTimes[j];
237     }
238     fIsMC=copy.fIsMC;
239     fSelectedPdg=copy.fSelectedPdg;
240     fP=copy.fP;
241     fPt=copy.fPt;
242     fEta=copy.fEta;
243     fPhi=copy.fPhi;
244     fTPCOuterPhi=copy.fTPCOuterPhi;
245     fL=copy.fL;
246     fMatchingMomCut=copy.fMatchingMomCut;
247     fTof=copy.fTof;
248     fHlist=copy.fHlist;
249     fHlistTimeZero=copy.fHlistTimeZero;
250     fHlistPID=copy.fHlistPID;
251     fHlistTRD=copy.fHlistTRD;
252     fHlistTrigger=copy.fHlistTrigger;
253   }
254   return *this;
255 }
256
257 //___________________________________________________________________________
258 AliAnalysisTaskTOFqaID::~AliAnalysisTaskTOFqaID() {
259   //
260   //destructor
261   //
262
263   Info("~AliAnalysisTaskTOFqaID","Calling Destructor");
264   if (fESDpid) delete fESDpid;
265   if (fTOFHeader) delete fTOFHeader;
266   if (fTOFT0v1) delete fTOFT0v1;
267   if (fVertex) delete fVertex;
268   if (fTrackFilter) delete fTrackFilter;
269   if (AliAnalysisManager::GetAnalysisManager()->IsProofMode()) return;  
270
271   if (fHlist) {
272     delete fHlist;
273     fHlist = 0;
274   }
275   if (fHlistTimeZero) {
276     delete fHlistTimeZero;
277     fHlistTimeZero = 0;
278   }
279   if (fHlistPID){
280     delete fHlistPID;
281     fHlistPID = 0;
282   }
283   if (fHlistTRD){
284     delete fHlistTRD;
285     fHlistTRD = 0;
286   }
287   if (fHlistTrigger){
288     delete fHlistTrigger;
289     fHlistTrigger = 0;
290   }
291 }
292
293 //________________________________________________________________________
294 void AliAnalysisTaskTOFqaID::UserCreateOutputObjects()
295 {
296   //Defines output objects and histograms
297   Info("CreateOutputObjects","CreateOutputObjects (TList) of task %s", GetName());
298   OpenFile(1);
299
300   if (!fHlist) fHlist = new TList();    
301   fHlist->SetOwner(kTRUE);
302   fHlist->SetName("base");
303
304   if (!fHlistTimeZero) fHlistTimeZero = new TList();    
305   fHlistTimeZero->SetOwner(kTRUE);
306   fHlistTimeZero->SetName("startTime");
307
308   if (!fHlistPID) fHlistPID = new TList();      
309   fHlistPID->SetOwner(kTRUE);
310   fHlistPID->SetName("pid");
311
312   if (!fHlistTRD) fHlistTRD = new TList();      
313   fHlistTRD->SetOwner(kTRUE);  
314   fHlistTRD->SetName("TRD");
315
316   if (!fHlistTrigger) fHlistTrigger = new TList();      
317   fHlistTrigger->SetOwner(kTRUE);  
318   fHlistTrigger->SetName("trigger");
319
320   if (fExpTimeRangeMax<fExpTimeRangeMin) {
321     SetExpTimeHistoRange(-25010.,25010.);
322   }
323   fnExpTimeBins = TMath::Nint((fExpTimeRangeMax - fExpTimeRangeMin)/fExpTimeBinWidth);//ps
324   fExpTimeRangeMax=fExpTimeRangeMin+fnExpTimeBins*fExpTimeBinWidth;//ps
325   
326   if (fExpTimeSmallRangeMax<fExpTimeSmallRangeMin) {
327     SetExpTimeHistoSmallRange(-5002.,5002.);
328   }
329   fnExpTimeSmallBins = TMath::Nint((fExpTimeSmallRangeMax - fExpTimeSmallRangeMin)/fExpTimeBinWidth);//ps
330   fExpTimeSmallRangeMax=fExpTimeSmallRangeMin+fnExpTimeSmallBins*fExpTimeBinWidth;//ps
331   
332   //add plots for start time QA
333   AddStartTimeHisto(fHlistTimeZero,"");
334   
335   //add plots for base TOF quantities
336   AddTofBaseHisto(fHlist,  1, "");
337   AddTofBaseHisto(fHlist, -1, "");
338   
339   //add plots for matching efficiency
340   AddMatchingEffHisto(fHlist,  1, "");
341   AddMatchingEffHisto(fHlist, -1, "");
342
343   //add plots for PID checks
344   AddPidHisto(fHlistPID,  1, ""); 
345   AddPidHisto(fHlistPID, -1, ""); 
346
347   //add trd plots
348   if (fEnableAdvancedCheck) {
349     AddTrdHisto();
350     AddTofTrgHisto("");
351   }
352   PostData(1, fHlist);
353   PostData(2, fHlistTimeZero);
354   PostData(3, fHlistPID);
355   PostData(4, fHlistTRD);
356   PostData(5, fHlistTrigger);
357 }
358 //________________________________________________________________________
359 void AliAnalysisTaskTOFqaID::UserExec(Option_t *) 
360
361   /* Main - executed for each event.
362     It extracts event information and track information after selecting 
363     primary tracks via standard cuts. */
364   
365   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
366   if (!esdH) {
367     AliError("ERROR: Could not get ESDInputHandler");
368     return;
369   } else {
370     fESD = (AliESDEvent*) esdH->GetEvent();
371   } 
372   
373   if (!fESD) {
374     AliError("ERROR: fESD not available");
375     return;
376   }
377
378   //access MC event handler for MC truth information
379   if (fIsMC) {
380     AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
381     if (!mcH) {
382       AliError("ERROR: Could not get MCeventHandler");
383       return;
384     } else {
385       fMCevent = (AliMCEvent *) mcH->MCEvent();
386       if (!fMCevent) {
387         AliError("Trying to retrieve an invalid MC event.");
388         return;
389       }
390     }
391   }
392   
393   // get run number
394   Int_t runNb = fESD->GetRunNumber();
395   if (runNb>0) fRunNumber = runNb;
396   //reset matched track counters
397   for (Int_t j=0;j<3;j++){fNTOFtracks[j]=0;}  
398   
399   //Get vertex info and apply vertex cut
400   if (!IsEventSelected(fESD)) return;
401   
402   //set response tof_t0 for all other checks
403   fESDpid->SetTOFResponse(fESD,AliESDpid::kTOF_T0);//(fill_t0, tof_t0, t0_t0, best_t0)
404   
405   Printf("Momentum cut for eta and phi distributions set: Pt>%3.2f", fMatchingMomCut);
406
407   // loop over ESD tracks
408   for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
409     AliESDtrack* track = fESD->GetTrack(iTracks);
410     if (!track) {
411       Printf("ERROR: Could not receive track %d", iTracks);
412       continue;
413     }
414     
415     //primary tracks selection: kTPCrefit and std cuts
416     if (fTrackFilter){
417       if(!fTrackFilter->IsSelected(track)) continue;
418     } else {
419       Printf("No track filter found, skipping the track loop");
420       return;
421     }
422     
423     //select specie if MC
424     if ( fIsMC && 
425          (!SelectMCspecies(fMCevent, track))) {
426       AliDebug(4, Form("MC tracks selection: Track=%i  label=%i  Not Accepted", iTracks, track->GetLabel()));
427       continue;
428     }
429     
430     //apply cut for eta acceptance
431     fEta=track->Eta();
432     if (TMath::Abs(fEta)>0.8) continue; 
433     
434     //get other track variables
435     fP = track->P();
436     fPt = track->Pt();
437     fPhi = track->Phi()*TMath::RadToDeg();
438     fTPCOuterPhi = GetPhiAtTPCouterRadius(track);
439     fL = track->GetIntegratedLength();
440     track->GetIntegratedTimes(fTrkExpTimes);
441     
442     Int_t charge = track->Charge();
443     
444     FillPrimaryTrkHisto(charge,"");
445     if (IsTPCTOFMatched(track)) {     
446       fTof=track->GetTOFsignal()*1E-3;//in ps
447       //increment track counters
448       fNTOFtracks[0]++;
449       if (charge>0) fNTOFtracks[1]++;
450       if (charge<0) fNTOFtracks[2]++;
451       //fill histos
452       FillTofBaseHisto(track, charge,"");
453       FillMatchedTrkHisto(charge,"");
454       FillPidHisto(track, charge, "");
455     }    
456     if (fEnableAdvancedCheck) FillTrdHisto(track, charge);
457   }//end loop on tracks
458    
459   //fill time zero histos  
460   FillStartTimeHisto("");  
461   ((TH1F*)fHlist->FindObject("hTOFmulti_pos"))->Fill(fNTOFtracks[1]);
462   ((TH1F*)fHlist->FindObject("hTOFmulti_neg"))->Fill(fNTOFtracks[2]);
463   
464   //fill TOF trg histos from infos in TOF header
465   fTOFHeader=(AliTOFHeader*)fESD->GetTOFHeader();
466   if (!fTOFHeader) {
467     AliWarning("Cannot get TOF header: no TOF trigger info available");
468   } else {
469     FillTofTrgHisto("");
470   }
471   
472   PostData(1, fHlist);
473   PostData(2, fHlistTimeZero);
474   PostData(3, fHlistPID);
475   PostData(4, fHlistTRD);
476   PostData(5, fHlistTrigger);
477
478 }      
479
480 //________________________________________________________________________
481 void AliAnalysisTaskTOFqaID::Terminate(Option_t *) 
482 {
483   //check on output validity
484   fHlist = dynamic_cast<TList*> (GetOutputData(1));
485   if (!fHlist) {
486     Printf("ERROR: lists not available");
487     return;   
488   } 
489   
490   // TH1D*hDummy = ((TH1D*)fHlist->FindObject("hTOFmatchedESDPt"));
491   // TH1D*hMatchingEff = (TH1D*) hDummy->Clone("hMatchingEff");
492   // hMatchingEff->SetTitle("Matching efficiency");
493   // hMatchingEff->Divide((TH1F*) fHlist->FindObject("hESDprimaryTrackPt"));
494   // TCanvas *c1 = new TCanvas("AliAnalysisTaskTOFqaID","Matching vs Pt",10,10,510,510);
495   // c1->cd(1)->SetLogy();
496   // hMatchingEff->DrawCopy("E");
497   // fHlist->AddLast(hMatchingEff);
498   ComputeMatchingEfficiency(fHlist, "pt");
499   ComputeMatchingEfficiency(fHlist, "eta");
500   ComputeMatchingEfficiency(fHlist, "phi");
501
502   PostData(1, fHlist);
503 }
504
505 //---------------------------------------------------------------
506 Int_t AliAnalysisTaskTOFqaID::GetStripIndex(const Int_t * const in)
507 {
508   /* return tof strip index between 0 and 91 */
509   
510   Int_t nStripA = AliTOFGeometry::NStripA();
511   Int_t nStripB = AliTOFGeometry::NStripB();
512   Int_t nStripC = AliTOFGeometry::NStripC();
513
514   Int_t iplate = in[1];
515   Int_t istrip = in[2];
516   
517   Int_t stripOffset = 0;
518   switch (iplate) {
519   case 0:
520     stripOffset = 0;
521       break;
522   case 1:
523     stripOffset = nStripC;
524     break;
525   case 2:
526     stripOffset = nStripC+nStripB;
527     break;
528   case 3:
529     stripOffset = nStripC+nStripB+nStripA;
530     break;
531   case 4:
532     stripOffset = nStripC+nStripB+nStripA+nStripB;
533     break;
534   default:
535     stripOffset=-1;
536     break;
537   };
538   
539   if (stripOffset<0 || stripOffset>92) return -1;
540   else 
541     return (stripOffset+istrip);
542 }
543
544 //-----------------------------------------------------------------
545 Double_t AliAnalysisTaskTOFqaID::GetPhiAtTPCouterRadius(AliESDtrack * track)
546 {
547   //get track phi at TPC outer radius
548   if (!track) return 1e10;
549   Double_t tpcoutcoord[3]={0.,0.,0.};
550   track->GetOuterXYZ(tpcoutcoord);
551   Double_t phiOuterTPC=TMath::ATan2(tpcoutcoord[1],tpcoutcoord[0])*TMath::RadToDeg();
552   if (phiOuterTPC<0) 
553     phiOuterTPC+= (2*TMath::Pi()*TMath::RadToDeg());
554   return phiOuterTPC;
555 }
556 //-----------------------------------------------------------------
557 Bool_t  AliAnalysisTaskTOFqaID::IsEventSelected(AliESDEvent * event)
558 {
559   //select event based on primary vertex
560   if (!event) {
561     AliError("Invalid ESD event");
562     return kFALSE;
563   }
564   fVertex = (AliESDVertex*) event->GetPrimaryVertexTracks(); 
565   if(fVertex->GetNContributors()<1) { 
566     // SPD vertex
567     fVertex = (AliESDVertex*) event->GetPrimaryVertexSPD(); 
568     if(fVertex->GetNContributors()<1) fVertex = 0x0;
569   }
570   if (!fVertex) return kFALSE; 
571   if (TMath::Abs(fVertex->GetZv())<10.0) return kTRUE;
572   else return kFALSE;
573 }
574
575 //-----------------------------------------------------------------
576 Bool_t  AliAnalysisTaskTOFqaID::IsTPCTOFMatched(AliESDtrack * track)
577 {
578   //defines TOF matching
579   if (!track){
580     AliWarning("Invalid track object");
581     return kFALSE;
582   }
583   
584   if ( (track->IsOn(AliESDtrack::kTOFout)) &&
585        (track->IsOn(AliESDtrack::kTIME)) &&
586        (track->IsOn(AliESDtrack::kTPCout))  )
587     return kTRUE;
588   else 
589     return kFALSE;
590 }
591 //-----------------------------------------------------------------
592 Bool_t  AliAnalysisTaskTOFqaID::IsInTRD(AliESDtrack * track)
593 {
594   //defines cut to select particles in/out TRD
595   if (!track){
596     AliWarning("Invalid track object");
597     return kFALSE;
598   }
599   
600   if ( track->IsOn(AliESDtrack::kTPCout) 
601        && track->IsOn(AliESDtrack::kTRDout) )
602     return kTRUE;
603   else 
604     return kFALSE; 
605 }
606 //-----------------------------------------------------------------
607 void AliAnalysisTaskTOFqaID::FillStartTimeMaskHisto(TString suffix)
608 {
609   /* set pid response to use best_T0 and for each
610      accepted track fills the histogram with the 
611      used start time 
612   */
613
614   //set response best_t0 
615   //fESDpid->SetTOFResponse(fESD,AliESDpid::kBest_T0);
616
617   for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
618     AliESDtrack* track = fESD->GetTrack(iTracks);
619     if (!track) {
620       Printf("ERROR: Could not receive track %d", iTracks);
621       continue;
622     }    
623     //primary tracks selection: kTPCrefit and std cuts
624     if (fTrackFilter){
625       if(!fTrackFilter->IsSelected(track)) continue;
626     }
627     else{
628       Printf("No track filter found, skipping the track loop");
629       break;
630     }
631     if (TMath::Abs(track->Eta())>0.8) continue; //cut for acceptance  
632     
633     Int_t StartTimeBit = fESDpid->GetTOFResponse().GetStartTimeMask(track->P());
634     ((TH2F*)fHlistTimeZero->FindObject(Form("hStartTimeMask%s",suffix.Data())))->Fill(track->P(),StartTimeBit);
635     
636     //matched tracks selection: kTOFout and kTIME
637     if ( (track->IsOn(AliESDtrack::kTOFout)) &&
638          (track->IsOn(AliESDtrack::kTIME)) &&
639          (track->IsOn(AliESDtrack::kTPCout))  ) {
640       ((TH2F*)fHlistTimeZero->FindObject(Form("hStartTimeMaskMatched%s",suffix.Data())))->Fill(track->P(),StartTimeBit);
641     }
642   }
643   return;
644 }
645
646 //----------------------------------------------------
647 Bool_t AliAnalysisTaskTOFqaID::ComputeTimeZeroByTOF1GeV()
648 {
649   /* compute T0-TOF for tracks within momentum range [0.95, 1.05] */
650   /* init T0-TOF */
651   fTOFT0v1->Init(fESD);
652   //AliTOFT0v1 *fTOFT0v1 = new AliTOFT0v1(fESDpid);
653   fTOFT0v1->DefineT0("all", 0.95, 1.05);
654   fMyTimeZeroTOF = -1000. * fTOFT0v1->GetResult(0);
655   fMyTimeZeroTOFsigma = 1000. * fTOFT0v1->GetResult(1);
656   fMyTimeZeroTOFtracks = fTOFT0v1->GetResult(3);
657   Bool_t hasTimeZeroTOF = kFALSE;
658   /* check T0-TOF sigma */
659   if (fMyTimeZeroTOFsigma < 250.)
660     hasTimeZeroTOF = kTRUE;  
661   return hasTimeZeroTOF;
662 }
663
664 //------------------------------------------------------
665 TString AliAnalysisTaskTOFqaID::GetSpeciesName(Int_t absPdgCode)
666 {
667   //returns name of selected specie 
668   TString name;
669   switch (absPdgCode){
670   case 11:
671     name = "electron";
672     break;
673   case 13:
674     name = "muon";
675     break;
676   case 211:
677     name = "pion";
678     break;
679   case 321:
680     name = "kaon";
681       break;
682   case 2212:
683     name = "proton";
684     break;
685   default:
686     name = "noPID";
687     break;
688   }
689   return name.Data();
690 }
691
692 //-----------------------------------------------
693 Bool_t AliAnalysisTaskTOFqaID::SelectMCspecies(AliMCEvent * ev, AliESDtrack * track)
694 {
695   //
696   //Retrieves particle true ID from MC and selects the desired species
697   //
698   if ((!ev) || (!track)) {
699     AliError("SelectMCspecies - Invalid object set as argument");
700     return kFALSE;
701   }
702   
703   if (fSelectedPdg==0) return kTRUE;   //if fSelectedPdg==0, no species selection is applied
704
705   Long_t label = track->GetLabel();
706   if (label<0) return kFALSE;
707   
708   // get number of particles
709   Long_t nMC = ev->GetNumberOfTracks();
710   // if label too large --> failed
711   if (label>= nMC) {
712     AliWarning(Form("Stack overflow: track label = %li -- stack maximum = %li", label, nMC));
713     return kFALSE;
714   } 
715   // retrieve particle
716   AliMCParticle *mcPart = (AliMCParticle *)ev->GetTrack(label);
717   if (!mcPart) {// if particle = NULL --> failed
718     AliWarning(Form("Stack discontinuity: label %li refers to a NULL object", label));
719     return kFALSE;
720   }
721
722   Int_t pdgCode = mcPart->PdgCode();
723   if (!(TMath::Abs(pdgCode)==fSelectedPdg))
724     return kFALSE;
725   else  //Printf("Found particle with selected pdg code = %i", pdgCode);
726   return  kTRUE;
727 }
728
729 //----------------------------------------------------------------------------------
730 Bool_t  AliAnalysisTaskTOFqaID::ComputeMatchingEfficiency(TList* list, TString variable)
731 {
732   //computes matching efficiency from previously filled histos
733   // to be called in terminate function  
734   if (!list) return kFALSE;
735
736   TString matchedName, primaryName, xAxisTitle;
737   if (variable.Contains("pt")) {
738     matchedName = "hTOFmatchedESDPt";
739     primaryName = "hESDprimaryTrackPt";
740     xAxisTitle = "p_{T} (GeV/c)";
741   }
742   if (variable.Contains("eta")) {
743     matchedName = "hTOFmatchedESDeta";
744     primaryName = "hTOFprimaryESDeta";
745     xAxisTitle = "#eta";
746   }
747   if (variable.Contains("phi")) {
748     matchedName = "hTOFmatchedESDphi";
749     primaryName = "hTOFprimaryESDphi";
750     xAxisTitle = "#phi_vtx (deg)";
751   }
752   
753   TH1D*hDummy = ((TH1D*)list->FindObject(matchedName.Data()));
754   if (!hDummy) return 0;
755
756   TH1D*hMatchingEff = (TH1D*) hDummy->Clone("hMatchingEff");
757   hMatchingEff->SetNameTitle(Form("hMatchingEff_%s", variable.Data()),Form("Matching efficiency vs %s", variable.Data()));
758   hMatchingEff->Divide((TH1F*) list->FindObject(primaryName.Data()));
759   hMatchingEff->GetXaxis()->SetTitle(xAxisTitle.Data());
760   hMatchingEff->GetYaxis()->SetRangeUser(0.0,1.0);
761   hMatchingEff->GetYaxis()->SetTitle("#epsilon_{match}");
762   list->AddLast(hMatchingEff);
763   return 1;
764 }
765 //----------------------------------------------------------------------------------
766 void AliAnalysisTaskTOFqaID::HistogramMakeUp(TH1* hist, Color_t color, Int_t markerStyle, TString drawOpt, TString newName, TString newTitle, TString xTitle, TString yTitle)
767 {
768   //set histogram style and axes style at once
769   if (!hist) return;  
770   if (!newName.IsNull()) hist->SetName(newName.Data());
771   if (!newTitle.IsNull()) hist->SetTitle(newTitle.Data());
772   if (!xTitle.IsNull()) hist->GetXaxis()->SetTitle(xTitle.Data());
773   if (!yTitle.IsNull()) hist->GetYaxis()->SetTitle(yTitle.Data());
774   hist->SetLineColor(color);
775   hist->SetMarkerColor(color);
776   hist->SetMarkerStyle(markerStyle);
777   hist->SetMarkerSize(0.7);
778   hist->SetDrawOption(drawOpt.Data());
779   hist->Sumw2();
780   return;
781 }
782
783 //----------------------------------------------------------------------------------
784 void AliAnalysisTaskTOFqaID::AddTofBaseHisto(TList *list, Int_t charge, TString suffix)
785 {
786   //Creates histograms for monitoring TOF signal, time alignement and matching-related quantities
787   if (!list){
788     AliError("Invalid list passed as argument.");
789     return;
790   }
791  
792   TString cLabel; 
793   if (charge<0) cLabel.Form("neg"); 
794   else 
795     if (charge>0) cLabel.Form("pos"); 
796     else cLabel.Form("all");
797   
798   TH1I* hTOFmulti = new TH1I(Form("hTOFmulti%s_%s",suffix.Data(), cLabel.Data()), Form("%s matched trk per event (|#eta|#leq0.8, p_{T}#geq0.3 GeV/c)", cLabel.Data()), 100, 0, 100);  
799   HistogramMakeUp(hTOFmulti, ((charge>0)? kRed : kBlue+2), 1, "E1", "","", "N","events");
800   list->AddLast(hTOFmulti);
801   
802   TH1F* hTOFtime = new TH1F(Form("hTime%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk TOF signal", cLabel.Data()), 250, 0., 610. ) ; 
803   HistogramMakeUp(hTOFtime,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "t (ns)","tracks");
804   list->AddLast(hTOFtime);
805   
806   TH1F* hTOFrawTime = new TH1F(Form("hRawTime%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk TOF raw signal", cLabel.Data()), 250, 0., 610. ) ; 
807   HistogramMakeUp(hTOFrawTime,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "t_{raw} (ns)","tracks");
808   list->AddLast(hTOFrawTime);
809
810   TH1F* hTOFtot = new TH1F(Form("hTot%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk ToT", cLabel.Data()), 50, 0., 50. ) ; 
811   HistogramMakeUp(hTOFtot,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "ToT (ns)","tracks");
812   list->AddLast(hTOFtot);
813     
814   TH1F* hMatchedL  = new TH1F(Form("hMatchedL%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk lenght", cLabel.Data()), 1200, -400., 800) ; 
815   HistogramMakeUp(hMatchedL,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "L (cm)","tracks");
816   list->AddLast(hMatchedL);
817      
818   TH2F* hMatchedDxVsPt = new TH2F(Form("hMatchedDxVsPt%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk dx vs.p_{T}", cLabel.Data()), 1000,0.,10.,200, -10., 10.) ; 
819   HistogramMakeUp(hMatchedDxVsPt,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","", "GeV/c (p_{T})","dx (cm)");
820   list->AddLast(hMatchedDxVsPt); 
821
822   TH2F* hMatchedDzVsStrip = new TH2F(Form("hMatchedDzVsStrip%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk dz vs. strip (#eta)", cLabel.Data()), 92,0.,92.,200, -10., 10.) ; 
823   HistogramMakeUp(hMatchedDzVsStrip,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","", "strip","dz (cm)");
824   list->AddLast(hMatchedDzVsStrip) ; 
825
826   TProfile *hMatchedDxVsCh = new TProfile(Form("hMatchedDxVsCh%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk dx vs. channel", cLabel.Data()), 157248., 0.,157248.);
827   HistogramMakeUp(hMatchedDxVsCh,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "channel","dx (cm)");
828   list->AddLast(hMatchedDxVsCh);
829    
830   TProfile *hMatchedDzVsCh = new TProfile(Form("hMatchedDzVsCh%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk dz vs. channel", cLabel.Data()), 157248., 0.,157248.);
831   HistogramMakeUp(hMatchedDzVsCh,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "channel","dz (cm)");
832   list->AddLast(hMatchedDzVsCh);    
833
834   return;
835 }
836
837 //----------------------------------------------------------------------------------
838 void    AliAnalysisTaskTOFqaID::AddMatchingEffHisto(TList *list, Int_t charge, TString suffix)
839 {
840   if (!list){
841     AliError("Invalid list passed as argument.");
842     return;
843   }
844   TString cLabel; 
845   if (charge<0) cLabel.Form("neg"); else if (charge>0) cLabel.Form("pos"); else cLabel.Form("all");
846   
847   TH1F* hMatchedP  = new TH1F(Form("hMatchedP%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p", cLabel.Data()), 1000,0.,10.) ;  
848   HistogramMakeUp(hMatchedP,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "L (cm)","tracks");
849   list->AddLast(hMatchedP) ; 
850
851   TH1F* hMatchedPt  = new TH1F(Form("hMatchedPt%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p_{T}", cLabel.Data()), 1000,0.,10.) ;  
852   HistogramMakeUp(hMatchedPt,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "p_{T} (GeV/c)","tracks");
853   list->AddLast(hMatchedPt) ; 
854
855   TH1F* hMatchedEta = new TH1F(Form("hMatchedEta%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk #eta", cLabel.Data()), 200, -1., 1.) ; 
856   HistogramMakeUp(hMatchedEta,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "#eta","tracks");
857   list->AddLast(hMatchedEta) ; 
858
859   TH1F* hMatchedPhi = new TH1F(Form("hMatchedPhi%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk #phi_{vtx}", cLabel.Data()), 72, 0., 360.) ; 
860   HistogramMakeUp(hMatchedPhi,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "#phi_{vtx} (deg)","tracks");
861   list->AddLast(hMatchedPhi) ; 
862
863   TH2F* hMatchedPtVsOutPhi  = new TH2F(Form("hMatchedPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p_{T} vs. #phi_{TPC out}", cLabel.Data()), 72, 0.0, 360.0, 1000,0.,10.) ;  
864   HistogramMakeUp(hMatchedPtVsOutPhi,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","", "#phi_{TPC out} (deg)","p_{T} (GeV/c)");
865   list->AddLast(hMatchedPtVsOutPhi) ;
866    
867   TH1F* hPrimaryP  = new TH1F(Form("hPrimaryP%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p", cLabel.Data()), 1000,0.,10.) ;  
868   HistogramMakeUp(hPrimaryP,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "L (cm)","tracks");
869   list->AddLast(hPrimaryP) ; 
870
871   TH1F* hPrimaryPt  = new TH1F(Form("hPrimaryPt%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p_{T}", cLabel.Data()), 1000,0.,10.) ;  
872   HistogramMakeUp(hPrimaryPt,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "p_{T} (GeV/c)","tracks");
873   list->AddLast(hPrimaryPt) ; 
874
875   TH1F* hPrimaryEta = new TH1F(Form("hPrimaryEta%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk #eta", cLabel.Data()), 200, -1., 1.) ; 
876   HistogramMakeUp(hPrimaryEta,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "#eta","tracks");
877   list->AddLast(hPrimaryEta) ; 
878
879   TH1F* hPrimaryPhi = new TH1F(Form("hPrimaryPhi%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk #phi_{vtx}", cLabel.Data()), 72, 0., 360.) ; 
880   HistogramMakeUp(hPrimaryPhi,((charge>0)? kRed+2 : kBlue+2), 1, "E1", "","", "#phi_{vtx} (deg)","tracks");
881   list->AddLast(hPrimaryPhi) ; 
882    
883   TH2F* hPrimaryPtVsOutPhi  = new TH2F(Form("hPrimaryPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk p_{T} vs. #phi_{TPC out}", cLabel.Data()), 72, 0.0, 360.0, 1000,0.,10.) ;  
884   HistogramMakeUp(hPrimaryPtVsOutPhi,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","", "#phi_{TPC out} (deg)","p_{T} (GeV/c)");
885   list->AddLast(hPrimaryPtVsOutPhi) ; 
886   return;
887 }
888   
889 //----------------------------------------------------------------------------------
890 void  AliAnalysisTaskTOFqaID::AddPidHisto(TList *list, Int_t charge, TString suffix)
891 {
892   //Creates histograms for monitoring TOF PID
893   if (!list){
894     AliError("Invalid list passed as argument.");
895     return;
896   }
897   TString cLabel; 
898   if (charge<0) cLabel.Form("neg"); else if (charge>0) cLabel.Form("pos"); else cLabel.Form("all");
899     
900   TH2F* hMatchedBetaVsP  = new TH2F(Form("hMatchedBetaVsP%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched trk #beta vs. p", cLabel.Data()), 1000, 0.0, 10.0, 150, 0., 1.5) ; 
901   HistogramMakeUp(hMatchedBetaVsP,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","", "p (GeV/c)","#beta");
902   list->AddLast(hMatchedBetaVsP);
903     
904   TH1F* hMatchedMass= new TH1F(Form("hMatchedMass%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched p.le M", cLabel.Data()), 500, 0., 5. );
905   HistogramMakeUp(hMatchedMass,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "M (GeV/c^{2})","entries");
906   list->AddLast(hMatchedMass);
907     
908   TH1F* hMatchedMass2= new TH1F(Form("hMatchedMass2%s_%s",suffix.Data(),cLabel.Data()), Form("%s matched p.le M^{2}", cLabel.Data()), 500, 0., 10. );
909   HistogramMakeUp(hMatchedMass2,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "M^{2} (GeV^{2}/c^{4})","entries");
910   list->AddLast(hMatchedMass2);
911
912   TH2F* hExpTimePiVsStrip = new TH2F(Form("hExpTimePiVsStrip%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{#pi,exp} vs strip",cLabel.Data()), 92, 0, 92,  fnExpTimeSmallBins, fExpTimeSmallRangeMin, fExpTimeSmallRangeMax) ; 
913   HistogramMakeUp(hExpTimePiVsStrip,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "strip (#eta)","t_{TOF}-t_{#pi,exp} [ps]");
914   list->AddLast(hExpTimePiVsStrip);
915   
916   TH2F* hExpTimePiT0Sub1GeV = new TH2F(Form("hExpTimePiT0Sub1GeV%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk (0.95#leq p_{T}#leq 1.05 GeV/c) t_{TOF}-t_{#pi,exp}-t_{0}^{TOF}",cLabel.Data()), 500, 0., 500., fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
917   HistogramMakeUp(hExpTimePiT0Sub1GeV,((charge>0)? kRed+2 : kBlue+2), 1, "colz", "","","n. tracks used for t_{0}^{TOF}","t_{TOF}-t_{#pi,exp}-t_{0}^{TOF}");    
918   list->AddLast(hExpTimePiT0Sub1GeV) ;
919   
920   TH1F* hExpTimePi = new TH1F(Form("hExpTimePi%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{#pi,exp}",cLabel.Data()), fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
921   HistogramMakeUp(hExpTimePi,((charge>0)? kRed+2 : kBlue+2), 1, "", "","", "t_{TOF}-t_{#pi,exp} [ps]","tracks");
922   list->AddLast(hExpTimePi);
923   
924   TH2F* hExpTimePiVsP = new TH2F(Form("hExpTimePiVsP%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{#pi,exp}",cLabel.Data()), 1000, 0.0, 10.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
925   HistogramMakeUp(hExpTimePiVsP,kRed+2, 1, "colz", "","", "p (GeV/c)","t_{TOF}-t_{#pi,exp} [ps]");
926   list->AddLast(hExpTimePiVsP);
927   
928   TH2F* hExpTimeKaVsP = new TH2F(Form("hExpTimeKaVsP%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{K,exp}",cLabel.Data()), 1000, 0.0, 10.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
929   HistogramMakeUp(hExpTimeKaVsP,kBlue+2, 1, "colz", "","", "p (GeV/c)","t_{TOF}-t_{K,exp} [ps]");
930   list->AddLast(hExpTimeKaVsP);
931   
932   TH2F* hExpTimeProVsP = new TH2F(Form("hExpTimeProVsP%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{p,exp}",cLabel.Data()),1000, 0.0, 10.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
933   HistogramMakeUp(hExpTimeProVsP,kGreen+2, 1, "colz", "","", "p (GeV/c)","t_{TOF}-t_{p,exp} [ps]");
934   list->AddLast(hExpTimeProVsP);
935   
936   TH2F* hTOFpidSigmaPi = new TH2F(Form("hTOFpidSigmaPi%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk n#sigma^{TOF}_{#pi} vs p_{T}",cLabel.Data()), 500,0.,5.,200, -10., 10. ) ; 
937   HistogramMakeUp(hTOFpidSigmaPi,kRed+2, 1, "colz", "","", "p (GeV/c)","n#sigma_{#pi,exp} [ps]");
938   list->AddLast(hTOFpidSigmaPi) ;
939   
940   TH2F* hTOFpidSigmaKa = new TH2F(Form("hTOFpidSigmaKa%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk n#sigma^{TOF}_{K} vs p_{T}",cLabel.Data()), 500, 0.,5.,200, -10., 10. ) ; 
941   HistogramMakeUp(hTOFpidSigmaKa,kBlue+2, 1, "colz", "","", "p (GeV/c)","n#sigma_{K,exp} [ps]");
942   list->AddLast(hTOFpidSigmaKa) ;
943     
944   TH2F* hTOFpidSigmaPro = new TH2F(Form("hTOFpidSigmaPro%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk TOF n#sigma^{TOF}_{p} vs p_{T}",cLabel.Data()), 500, 0.,5.,200, -10., 10. ) ; 
945   HistogramMakeUp(hTOFpidSigmaPro,kGreen+2, 1, "colz", "","","p (GeV/c)","n#sigma_{p,exp} [ps]");
946   list->AddLast(hTOFpidSigmaPro);
947     
948   TH2F* hExpTimePiT0SubVsP = new TH2F(Form("hExpTimePiT0SubVsP%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk t_{TOF}-t_{#pi,exp}-t_{0}^{TOF}",cLabel.Data()), 1000, 0.,10., fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
949   HistogramMakeUp(hExpTimePiT0SubVsP,kRed+2, 1, "colz", "","","p (GeV/c)","t_{TOF}-t_{#pi,exp}-t_{0}^{TOF}");
950   list->AddLast(hExpTimePiT0SubVsP) ;
951     
952   TH2F* hExpTimeKaT0SubVsP = new TH2F(Form("hExpTimeKaT0SubVsP%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk t_{TOF}-t_{K,exp}-t_{0}^{TOF}",cLabel.Data()), 1000, 0.,10., fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
953   HistogramMakeUp(hExpTimeKaT0SubVsP,kBlue+2, 1, "colz", "","","p (GeV/c)","t_{TOF}-t_{K,exp}-t_{0}^{TOF}");
954   list->AddLast(hExpTimeKaT0SubVsP) ;
955     
956   TH2F* hExpTimeProT0SubVsP = new TH2F(Form("hExpTimeProT0SubVsP%s_%s",suffix.Data(),cLabel.Data()), Form("%s trk t_{TOF}-t_{p,exp}-t_{0}^{TOF}",cLabel.Data()), 1000, 0.,10., fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
957   HistogramMakeUp(hExpTimeProT0SubVsP,kGreen+2, 1, "colz", "","","p (GeV/c)","t_{TOF}-t_{p,exp}-t_{0}^{TOF}");
958   list->AddLast(hExpTimeProT0SubVsP) ;
959   
960   TH2F* hExpTimePiVsOutPhi = new TH2F(Form("hExpTimePiVsOutPhi%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{#pi,exp} vs #phi_{TPC out}",cLabel.Data()), 72, 0.0, 360.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
961   HistogramMakeUp(hExpTimePiVsOutPhi,kRed+2, 1, "colz", "","", "#phi_{TPC out} (deg)","t_{TOF}-t_{#pi,exp} [ps]");
962   list->AddLast(hExpTimePiVsOutPhi);
963   
964   TH2F* hExpTimeKaVsOutPhi = new TH2F(Form("hExpTimeKaVsOutPhi%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{K,exp} vs #phi_{TPC out}",cLabel.Data()), 72, 0.0, 360.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
965   HistogramMakeUp(hExpTimeKaVsOutPhi,kBlue+2, 1, "colz", "","", "#phi_{TPC out} (deg)","t_{TOF}-t_{K,exp} [ps]");
966   list->AddLast(hExpTimeKaVsOutPhi);
967   
968   TH2F* hExpTimeProVsOutPhi = new TH2F(Form("hExpTimeProVsOutPhi%s_%s",suffix.Data(),cLabel.Data()),Form("%s matched trk t_{TOF}-t_{p,exp} vs #phi_{TPC out}",cLabel.Data()), 72, 0.0, 360.0, fnExpTimeBins, fExpTimeRangeMin, fExpTimeRangeMax) ; 
969   HistogramMakeUp(hExpTimeProVsOutPhi,kGreen+2, 1, "colz", "","", "#phi_{TPC out} (deg)","t_{TOF}-t_{p,exp} [ps]");
970   list->AddLast(hExpTimeProVsOutPhi);
971   
972   return;
973 }
974 //----------------------------------------------------------------------------------
975 void    AliAnalysisTaskTOFqaID::AddStartTimeHisto(TList *list, TString suffix)
976 {  
977   //Creates histograms for monitoring T0 signal and start-time related quantities
978   if (!list){
979     AliError("Invalid list passed as argument.");
980     return;
981   }
982   TH1D* hT0AC = new TH1D(Form("hT0AC%s",suffix.Data()), "Event timeZero from T0A&C; t_{0,AC} [ps]; events", 1000, -25000., 25000. ) ; 
983   HistogramMakeUp(hT0AC, kRed+2, 20, "", "","","","");    
984   list->AddLast(hT0AC);
985
986   TH1D* hT0A = new TH1D(Form("hT0A%s",suffix.Data()), "Event timeZero from T0A; t_{0,A} [ps]; events", 1000, -25000., 25000. ) ; 
987   HistogramMakeUp(hT0A, kBlue+2, 25, "", "","","","");    
988   list->AddLast(hT0A);
989     
990   TH1D* hT0C = new TH1D(Form("hT0C%s",suffix.Data()), "Event timeZero from T0C; t_{0,C} [ps]; events", 1000, -25000., 25000. ) ; 
991   HistogramMakeUp(hT0C, kGreen+2, 28, "", "","","","");    
992   list->AddLast(hT0C);
993     
994   TH1D* hT0DetRes = new TH1D(Form("hT0DetRes%s",suffix.Data()), "T0 detector (T0A-T0C)/2; (T0A-T0C)/2 [ps]; events", 200, -500.,500. ) ; 
995   HistogramMakeUp(hT0DetRes, kMagenta+1, 1, "", "","","","");    
996   list->AddLast(hT0DetRes) ; 
997
998   TH1F* hT0fill = new TH1F(Form("hT0fill%s",suffix.Data()), "Event timeZero of fill; t_{0,fill} [ps]; events", 1000, -25000., 25000. ) ; 
999   HistogramMakeUp(hT0fill, kOrange+1, 25, "", "","","","");    
1000   list->AddLast(hT0fill) ; 
1001     
1002   TH1F* hT0TOF = new TH1F(Form("hT0TOF%s",suffix.Data()), "Event timeZero estimated by TOF; t0 [ps]; events", 1000, -25000., 25000. ) ; 
1003   HistogramMakeUp(hT0TOF, kTeal-5, 21, "", "","","","");    
1004   list->AddLast(hT0TOF) ; 
1005     
1006   TH1F* hT0T0 = new TH1F(Form("hT0T0%s",suffix.Data()), "Best timeZero between AC, A, C; t_{0} [ps]; events", 1000, -25000.,25000. ) ; 
1007   HistogramMakeUp(hT0T0, kAzure+7, 26, "", "","","","");    
1008   list->AddLast(hT0T0) ; 
1009     
1010   TH1F* hT0best = new TH1F(Form("hT0best%s",suffix.Data()), "Event timeZero estimated as T0best; t0 [ps]; events", 1000, -25000.,25000. ) ; 
1011   HistogramMakeUp(hT0best, kBlack, 20, "", "","","","");    
1012   list->AddLast(hT0best) ; 
1013
1014   TH1F* hT0fillRes = new TH1F(Form("hT0fillRes%s",suffix.Data()), "Resolution of fillT0; #sigma_{fillT0} [ps];events", 250, 0.,250. ) ; 
1015   HistogramMakeUp(hT0fillRes, kOrange+1, 25, "", "","","","");    
1016   list->AddLast(hT0fillRes) ; 
1017     
1018   TH1F* hT0TOFRes = new TH1F(Form("hT0TOFRes%s",suffix.Data()), "Resolution of timeZero from TOF; #sigma_{TOFT0} [ps];events", 250, 0.,250. ) ; 
1019   HistogramMakeUp(hT0TOFRes, kTeal-5, 21, "", "","","","");    
1020   list->AddLast(hT0TOFRes) ; 
1021
1022   TH1F* hT0T0Res = new TH1F(Form("hT0T0Res%s",suffix.Data()), "Resolution of timeZero from T0;#sigma_{T0T0}  [ps];events", 250, -0., 250. ) ; 
1023   HistogramMakeUp(hT0T0Res, kAzure+7, 26, "", "","","","");    
1024   list->AddLast(hT0T0Res) ; 
1025   
1026   TH1F* hT0bestRes = new TH1F(Form("hT0bestRes%s",suffix.Data()), "Resolution of bestT0; #sigma_{bestT0} [ps];events", 250, 0.,250. ) ; 
1027   HistogramMakeUp(hT0bestRes, kBlack, 20, "", "","","","");    
1028   list->AddLast(hT0bestRes) ; 
1029
1030   TH2F* hT0TOFvsNtrk = new TH2F(Form("hT0TOFvsNtrk%s",suffix.Data()), "Event timeZero estimated by TOF vs. number of tracks in event;TOF-matching tracks; t0 [ps]", 100, 0., 100.,1000,-25000.,25000. ) ; 
1031   HistogramMakeUp(hT0TOFvsNtrk, kTeal-5, 1, "colz", "","","","");    
1032   list->AddLast(hT0TOFvsNtrk) ;
1033   
1034   TH2F* hEventT0MeanVsVtx = new TH2F(Form("hEventT0MeanVsVtx%s",suffix.Data()), "T0 detector: mean vs vertex ; (t0_{A}-t0_{C})/2 [ns]; (t0_{A}+t0_{C})/2 [ns]; events", 500, -25., 25., 500, -25., 25. ) ; 
1035   HistogramMakeUp(hEventT0MeanVsVtx, kBlue+2, 1, "colz", "","","","");    
1036   list->AddLast(hEventT0MeanVsVtx) ;
1037     
1038   TH2F* hEventV0MeanVsVtx = new TH2F(Form("hEventV0MeanVsVtx%s",suffix.Data()), "V0 detector: mean vs vertex ; (V0_{A}-V0_{C})/2 [ns]; (V0_{A}+V0_{C})/2 [ns]; events", 500, -50., 50., 500, -50., 50. ) ; 
1039   HistogramMakeUp(hEventV0MeanVsVtx, kBlack, 1, "colz", "","","","");    
1040   list->AddLast(hEventV0MeanVsVtx) ;
1041   
1042   const Double_t startTimeMomBins[13]={ 0.0, 0.3, 0.5, 0.6, 0.7, 0.8, 0.9, 1., 1.2, 1.5, 2., 3., 10.};
1043     
1044   TH2F* hStartTimeMaskMatched = new TH2F(Form("hStartTimeMaskMatched%s",suffix.Data()),"Start Time Mask vs p bin for matched tracks; p(GeV/c);", 12, startTimeMomBins, 8,0.,8.);
1045   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(1,"fill_t0");
1046   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(2,"tof_t0");
1047   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(3,"T0AC");
1048   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(4,"T0AC & tof_t0");
1049   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(5,"T0A");
1050   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(6,"T0A & tof_t0");
1051   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(7,"T0C");
1052   hStartTimeMaskMatched->GetYaxis()->SetBinLabel(8,"T0C & tof_t0");
1053   HistogramMakeUp(hStartTimeMaskMatched, kRed+2, 1, "colz", "","","","");    
1054   list->AddLast(hStartTimeMaskMatched);
1055     
1056   TH2F* hStartTimeMask = new TH2F(Form("hStartTimeMask%s",suffix.Data()),"Start Time Mask vs p bin for primary tracks; p(GeV/c);", 12, startTimeMomBins, 8,0.,8.);
1057   hStartTimeMask->GetYaxis()->SetBinLabel(1,"fill_t0");
1058   hStartTimeMask->GetYaxis()->SetBinLabel(2,"tof_t0");
1059   hStartTimeMask->GetYaxis()->SetBinLabel(3,"T0AC");
1060   hStartTimeMask->GetYaxis()->SetBinLabel(4,"T0AC & tof_t0");
1061   hStartTimeMask->GetYaxis()->SetBinLabel(5,"T0A");
1062   hStartTimeMask->GetYaxis()->SetBinLabel(6,"T0A & tof_t0");
1063   hStartTimeMask->GetYaxis()->SetBinLabel(7,"T0C");
1064   hStartTimeMask->GetYaxis()->SetBinLabel(8,"T0C & tof_t0");
1065   HistogramMakeUp(hStartTimeMask, kRed+2, 1, "colz", "","","","");    
1066   list->AddLast(hStartTimeMask);
1067
1068   return;
1069 }
1070 //----------------------------------------------------------------------------------
1071 void AliAnalysisTaskTOFqaID::AddTrdHisto()
1072 {
1073   //Creates histograms for monitoring TOF base quantities wrt TRD/no TRD selection
1074   if (!fHlistTRD){
1075     AliError("Invalid TRD list");
1076     return;
1077   }
1078
1079   AddMatchingEffHisto(fHlistTRD, 1, "_noTrd");
1080   AddMatchingEffHisto(fHlistTRD, -1, "_noTrd");
1081   AddMatchingEffHisto(fHlistTRD, 1, "_Trd");
1082   AddMatchingEffHisto(fHlistTRD, -1, "_Trd");
1083
1084   AddPidHisto(fHlistTRD, 1, "_noTrd");
1085   AddPidHisto(fHlistTRD, -1, "_noTrd");
1086   AddPidHisto(fHlistTRD, 1, "_Trd");
1087   AddPidHisto(fHlistTRD, -1, "_Trd");
1088   return;
1089 }
1090  
1091
1092 //----------------------------------------------------------------------------------
1093 void AliAnalysisTaskTOFqaID::AddTofTrgHisto(TString suffix)
1094 {
1095   //defines histo with trigger info
1096   if (!fHlistTrigger){
1097     AliError("Invalid TOF trigger list");
1098     return;
1099   }
1100   
1101   TH1I* hFiredMaxipad = new TH1I(Form("hFiredMaxipad%s",suffix.Data()), Form("Fired maxipad per event"), 1584, 0, 1584);  
1102   HistogramMakeUp(hFiredMaxipad, kBlue+2, 1, "E1", "","", "N_{maxipad}","events");
1103   fHlistTrigger->AddLast(hFiredMaxipad);
1104   
1105   TH1I* hFiredReadoutPad = new TH1I(Form("hFiredReadoutPad%s",suffix.Data()), Form("Fired readout pad per event"), 153000, 0, 153000);  
1106   HistogramMakeUp(hFiredReadoutPad, kRed+2, 1, "E1", "","", "N_{pad}","events");
1107   fHlistTrigger->AddLast(hFiredReadoutPad);
1108   
1109   TH1I* hFiredReadoutTrgPad = new TH1I(Form("hFiredReadoutTrgPad%s",suffix.Data()), Form("Fired readout pad in trg window"), 153000, 0, 153000);  
1110   HistogramMakeUp(hFiredReadoutTrgPad, kBlack, 1, "E1", "","", "N_{pad} in trg window","events");
1111   fHlistTrigger->AddLast(hFiredReadoutTrgPad);
1112   
1113   TH2I* hFiredMaxipadVsTrgPad = new TH2I(Form("hFiredMaxipadVsTrgPad%s",suffix.Data()), Form("Fired maxipad vs pads in trg window per event"), 100, 0, 100, 100, 0, 100);  
1114   HistogramMakeUp(hFiredMaxipadVsTrgPad, kBlue+2, 1, "colz", "","", "N_{pad} in trg window","N_{maxipad}");
1115   fHlistTrigger->AddLast(hFiredMaxipadVsTrgPad);
1116   
1117   TH2I* hTrgMap = new TH2I(Form("hTrgMap%s",suffix.Data()), Form("Map of fired maxipads"), 72, 0, 72, 22, 0, 22);  
1118   HistogramMakeUp(hTrgMap, kBlue+2, 1, "colz", "","", "LTM","maxipad");
1119   fHlistTrigger->AddLast(hTrgMap);
1120   
1121   return;
1122 }
1123
1124 //----------------------------------------------------------------------------------
1125 void AliAnalysisTaskTOFqaID::FillTofBaseHisto(AliESDtrack * track, Int_t charge, TString suffix)
1126 {
1127   //fill histo with TOF base quantities
1128   if (!track) return;
1129   
1130   //  Double_t tofTime=track->GetTOFsignal();//in ps
1131   Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
1132   Double_t tofToT=track->GetTOFsignalToT(); //in ps
1133   Int_t channel=track->GetTOFCalChannel(); 
1134   Int_t volId[5]; //(sector, plate,strip,padZ,padX)
1135   AliTOFGeometry::GetVolumeIndices(channel,volId);
1136   TString cLabel; if (charge<0) cLabel.Form("neg"); else if (charge>0) cLabel.Form("pos"); else cLabel.Form("all");
1137
1138   ((TH1F*)fHlist->FindObject(Form("hTime%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTof); //ns
1139   ((TH1F*)fHlist->FindObject(Form("hRawTime%s_%s",suffix.Data(),cLabel.Data())))->Fill(tofTimeRaw*1E-3); //ns
1140   ((TH1F*)fHlist->FindObject(Form("hTot%s_%s",suffix.Data(),cLabel.Data())))->Fill(tofToT);
1141   ((TH1F*)fHlist->FindObject(Form("hMatchedL%s_%s", suffix.Data(), cLabel.Data())))->Fill(fL);      
1142   ((TH2F*)fHlist->FindObject(Form("hMatchedDxVsPt%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt,track->GetTOFsignalDx());
1143   ((TH2F*)fHlist->FindObject(Form("hMatchedDzVsStrip%s_%s",suffix.Data(),cLabel.Data())))->Fill((Int_t)GetStripIndex(volId),track->GetTOFsignalDz());
1144   ((TProfile*)fHlist->FindObject(Form("hMatchedDxVsCh%s_%s",suffix.Data(),cLabel.Data())))->Fill(channel,track->GetTOFsignalDx());
1145   ((TProfile*)fHlist->FindObject(Form("hMatchedDzVsCh%s_%s",suffix.Data(),cLabel.Data())))->Fill(channel,track->GetTOFsignalDz());
1146   
1147   return;
1148 }
1149 //----------------------------------------------------------------------------------
1150 void AliAnalysisTaskTOFqaID::FillPrimaryTrkHisto(Int_t charge, TString suffix)
1151 {
1152   // fill histos with primary tracks info
1153   // => denominator for matching efficiency
1154   TString cLabel; 
1155   
1156   if (charge<0) cLabel.Form("neg"); else if (charge>0) cLabel.Form("pos"); else cLabel.Form("all");
1157   if (suffix.Contains("Trd")) {
1158     ((TH1F*)fHlistTRD->FindObject(Form("hPrimaryP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP); 
1159     ((TH1F*)fHlistTRD->FindObject(Form("hPrimaryPt%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt); 
1160     if (fPt>=fMatchingMomCut) {
1161       ((TH1F*)fHlistTRD->FindObject(Form("hPrimaryEta%s_%s",suffix.Data(),cLabel.Data())))->Fill(fEta);
1162       ((TH1F*)fHlistTRD->FindObject(Form("hPrimaryPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPhi);
1163       ((TH2F*)fHlistTRD->FindObject(Form("hPrimaryPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,fPt);
1164     }
1165   } else {
1166     ((TH1F*)fHlist->FindObject(Form("hPrimaryP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP); 
1167     ((TH1F*)fHlist->FindObject(Form("hPrimaryPt%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt); 
1168     if (fPt>=fMatchingMomCut) {
1169       ((TH1F*)fHlist->FindObject(Form("hPrimaryEta%s_%s",suffix.Data(),cLabel.Data())))->Fill(fEta);
1170       ((TH1F*)fHlist->FindObject(Form("hPrimaryPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPhi);
1171       ((TH2F*)fHlist->FindObject(Form("hPrimaryPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,fPt); 
1172     }
1173   }
1174   return;
1175 }
1176 //----------------------------------------------------------------------------------
1177 void AliAnalysisTaskTOFqaID::FillMatchedTrkHisto(Int_t charge, TString suffix)
1178 {
1179   //get matched tracks variables (matching cut to be applied externally)
1180   //=> numerator for matching efficiency
1181   TString cLabel; 
1182
1183   if (charge<0) cLabel.Form("neg"); 
1184   else 
1185     if (charge>0) cLabel.Form("pos"); 
1186     else cLabel.Form("all");
1187
1188   if (suffix.Contains("Trd")) { 
1189     ((TH1F*)fHlistTRD->FindObject(Form("hMatchedP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP); 
1190     ((TH1F*)fHlistTRD->FindObject(Form("hMatchedPt%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt); 
1191     if (fPt>=fMatchingMomCut) {
1192       ((TH1F*)fHlistTRD->FindObject(Form("hMatchedEta%s_%s",suffix.Data(),cLabel.Data())))->Fill(fEta);
1193       ((TH1F*)fHlistTRD->FindObject(Form("hMatchedPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPhi);
1194       ((TH2F*)fHlistTRD->FindObject(Form("hMatchedPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,fPt);
1195     }
1196   } else {
1197     ((TH1F*)fHlist->FindObject(Form("hMatchedP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP); 
1198     ((TH1F*)fHlist->FindObject(Form("hMatchedPt%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt); 
1199     if (fPt>=fMatchingMomCut) {
1200       ((TH1F*)fHlist->FindObject(Form("hMatchedEta%s_%s",suffix.Data(),cLabel.Data())))->Fill(fEta);
1201       ((TH1F*)fHlist->FindObject(Form("hMatchedPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPhi);
1202       ((TH2F*)fHlist->FindObject(Form("hMatchedPtVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,fPt);
1203     }
1204   }
1205   return;
1206 }
1207
1208 //----------------------------------------------------------------------------------
1209 void AliAnalysisTaskTOFqaID::FillPidHisto(AliESDtrack * track, Int_t charge, TString suffix)
1210 {
1211   //basic PID performance check
1212   if (fTof<=0) {
1213     printf("WARNING: track with negative TOF time found! Skipping this track for PID checks\n");
1214     return;
1215   }
1216   if (fL<=0){
1217     printf("WARNING: track with negative length found!Skipping this track for PID checks\n");
1218     return;
1219   }
1220   if (!track) return;
1221   
1222   TString cLabel; 
1223   if (charge<0) cLabel.Form("neg"); 
1224   else 
1225     if (charge>0) cLabel.Form("pos");
1226     else cLabel.Form("all");
1227   
1228   //calculate beta
1229   Double_t c=TMath::C()*1.E-9;// m/ns
1230   Double_t mass=0.; //GeV
1231   Double_t length=fL*0.01; // in meters
1232   Double_t tof=fTof*c;
1233   Double_t beta=length/tof;
1234   Double_t fact= (tof/length)*(tof/length) -1.;
1235   Double_t fP2 = fP*fP;
1236   
1237   if(fact<=0) {
1238     mass = -fP*TMath::Sqrt(-fact);
1239   } else { 
1240     mass = fP*TMath::Sqrt(fact); 
1241   }
1242   
1243   if (suffix.Contains("Trd")) {
1244     ((TH2F*) fHlistTRD->FindObject(Form("hMatchedBetaVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,beta);
1245     ((TH1F*) fHlistTRD->FindObject(Form("hMatchedMass%s_%s",suffix.Data(),cLabel.Data())))->Fill(mass);
1246     ((TH1F*) fHlistTRD->FindObject(Form("hMatchedMass2%s_%s",suffix.Data(),cLabel.Data())))->Fill(mass*mass);
1247   } else {
1248     ((TH2F*) fHlist->FindObject(Form("hMatchedBetaVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,beta);
1249     ((TH1F*) fHlist->FindObject(Form("hMatchedMass%s_%s",suffix.Data(),cLabel.Data())))->Fill(mass);
1250     ((TH1F*) fHlist->FindObject(Form("hMatchedMass2%s_%s",suffix.Data(),cLabel.Data())))->Fill(mass*mass);
1251   }
1252   
1253   //PID sigmas
1254   Bool_t isValidBeta[AliPID::kSPECIES]={0,0,0,0,0};
1255   for (Int_t specie = 0; specie < AliPID::kSPECIES; specie++){
1256     fSigmaSpecie[specie] = fESDpid->GetTOFResponse().GetExpectedSigma(fP, fTrkExpTimes[specie], AliPID::ParticleMass(specie));
1257     beta=1/TMath::Sqrt(1+AliPID::ParticleMass(specie)*AliPID::ParticleMass(specie)/(fP2));
1258     if (beta>0) {
1259       fThExpTimes[specie]=length*1.E3/(beta*c);//ps
1260       isValidBeta[specie]=kTRUE;
1261     } else {
1262       fThExpTimes[specie]=1E-10;
1263       isValidBeta[specie]=kFALSE;
1264     }
1265   }
1266   Float_t timeZeroTOF = (Float_t) fESDpid->GetTOFResponse().GetStartTime(fPt);
1267   Double_t tofps=fTof*1E3;//ps for t-texp
1268   Int_t channel=track->GetTOFCalChannel(); 
1269   Int_t volId[5]; //(sector, plate,strip,padZ,padX)
1270   AliTOFGeometry::GetVolumeIndices(channel,volId);
1271   
1272 if (suffix.Contains("Trd")) {
1273   if (isValidBeta[AliPID::kPion]){
1274     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimePiVsStrip%s_%s",suffix.Data(),cLabel.Data())))->Fill((Int_t)GetStripIndex(volId),tofps-fTrkExpTimes[AliPID::kPion]);//ps
1275     ((TH1F*)fHlistTRD->FindObject(Form("hExpTimePi%s_%s",suffix.Data(),cLabel.Data())))->Fill(tofps-fTrkExpTimes[AliPID::kPion]);//ps
1276     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimePiVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kPion]);
1277     ((TH2F*)fHlistTRD->FindObject(Form("hTOFpidSigmaPi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kPion])/fSigmaSpecie[AliPID::kPion]);
1278     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimePiT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kPion]-timeZeroTOF);  
1279     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimePiVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kPion]-timeZeroTOF);
1280     
1281     if (ComputeTimeZeroByTOF1GeV()){
1282       if ((fPt>0.95)&&(fPt<1.05)){
1283         ((TH2F*)fHlistTRD->FindObject(Form("hExpTimePiT0Sub1GeV%s_%s",suffix.Data(),cLabel.Data())))->Fill(fMyTimeZeroTOFtracks,tofps-fMyTimeZeroTOF-fTrkExpTimes[AliPID::kPion]);
1284       }
1285     } 
1286   }//end pion
1287   
1288   if (isValidBeta[AliPID::kKaon]){
1289     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeKaVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kKaon]);
1290     ((TH2F*)fHlistTRD->FindObject(Form("hTOFpidSigmaKa%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kKaon])/fSigmaSpecie[AliPID::kKaon]);
1291     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeKaT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kKaon]-timeZeroTOF);  
1292     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeKaVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kKaon]-timeZeroTOF);    
1293   }//end kaon  
1294   
1295   if (isValidBeta[AliPID::kProton]){
1296     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeProVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kProton]);
1297     ((TH2F*)fHlistTRD->FindObject(Form("hTOFpidSigmaPro%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kProton])/fSigmaSpecie[AliPID::kProton]);
1298     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeProT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kProton]-timeZeroTOF);  
1299     ((TH2F*)fHlistTRD->FindObject(Form("hExpTimeProVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kProton]-timeZeroTOF);
1300   }//end proton  
1301
1302  } else { 
1303
1304   if (isValidBeta[AliPID::kPion]){
1305     ((TH2F*)fHlist->FindObject(Form("hExpTimePiVsStrip%s_%s",suffix.Data(),cLabel.Data())))->Fill((Int_t)GetStripIndex(volId),tofps-fTrkExpTimes[AliPID::kPion]);//ps
1306     ((TH1F*)fHlist->FindObject(Form("hExpTimePi%s_%s",suffix.Data(),cLabel.Data())))->Fill(tofps-fTrkExpTimes[AliPID::kPion]);//ps
1307     ((TH2F*)fHlist->FindObject(Form("hExpTimePiVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kPion]);
1308     ((TH2F*)fHlist->FindObject(Form("hTOFpidSigmaPi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kPion])/fSigmaSpecie[AliPID::kPion]);
1309     ((TH2F*)fHlist->FindObject(Form("hExpTimePiT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kPion]-timeZeroTOF);  
1310     ((TH2F*)fHlist->FindObject(Form("hExpTimePiVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kPion]-timeZeroTOF);
1311     
1312     if (ComputeTimeZeroByTOF1GeV()){
1313       if ((fPt>0.95)&&(fPt<1.05)){
1314         ((TH2F*)fHlist->FindObject(Form("hExpTimePiT0Sub1GeV%s_%s",suffix.Data(),cLabel.Data())))->Fill(fMyTimeZeroTOFtracks,tofps-fMyTimeZeroTOF-fTrkExpTimes[AliPID::kPion]);
1315       }
1316     } 
1317   }//end pion
1318   
1319   if (isValidBeta[AliPID::kKaon]){
1320     ((TH2F*)fHlist->FindObject(Form("hExpTimeKaVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kKaon]);
1321     ((TH2F*)fHlist->FindObject(Form("hTOFpidSigmaKa%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kKaon])/fSigmaSpecie[AliPID::kKaon]);
1322     ((TH2F*)fHlist->FindObject(Form("hExpTimeKaT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kKaon]-timeZeroTOF);  
1323     ((TH2F*)fHlist->FindObject(Form("hExpTimeKaVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kKaon]-timeZeroTOF);    
1324   }//end kaon  
1325   
1326   if (isValidBeta[AliPID::kProton]){
1327     ((TH2F*)fHlist->FindObject(Form("hExpTimeProVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP, tofps-fTrkExpTimes[AliPID::kProton]);
1328     ((TH2F*)fHlist->FindObject(Form("hTOFpidSigmaPro%s_%s",suffix.Data(),cLabel.Data())))->Fill(fPt, (tofps-fTrkExpTimes[AliPID::kProton])/fSigmaSpecie[AliPID::kProton]);
1329     ((TH2F*)fHlist->FindObject(Form("hExpTimeProT0SubVsP%s_%s",suffix.Data(),cLabel.Data())))->Fill(fP,tofps-fTrkExpTimes[AliPID::kProton]-timeZeroTOF);  
1330     ((TH2F*)fHlist->FindObject(Form("hExpTimeProVsOutPhi%s_%s",suffix.Data(),cLabel.Data())))->Fill(fTPCOuterPhi,tofps-fTrkExpTimes[AliPID::kProton]-timeZeroTOF);
1331   }//end proton  
1332  }
1333   
1334   // if (fEnableAdvancedCheck && (fPt<1.)) {
1335   //   Double_t pos[3]={0.,0.,0.};
1336   //   track->GetXYZAt(378.,5.,pos);
1337   //   if ((pos[0]==0.)&&(pos[1]==0.)&&(pos[2]==0.))continue;
1338   
1339   //   Double_t phiTOF=TMath::ATan2(pos[1],pos[0])*TMath::RadToDeg();
1340   //   if (phiTOF<0) phiTOF+= (2*TMath::Pi()*TMath::RadToDeg());
1341   
1342   //   //fill t-texp vs phi@TOF
1343   //    if ((phiOuterTPC<=75) || ((phiOuterTPC>=125)&&(phiOuterTPC<=235)) || (phiOuterTPC>=305) ) { //TRD sectors 2012
1344   //    if ( ((phiOuterTPC>=85)&&(phiOuterTPC<=115)) || ((phiOuterTPC>=245)&&(phiOuterTPC<=295)) ) {//no TRD sectors 2012
1345   // }
1346   return;      
1347 }
1348 //----------------------------------------------------------------------------------
1349 void AliAnalysisTaskTOFqaID::FillStartTimeHisto(TString suffix)
1350 {
1351   //fill start time histo
1352   if (!fESD) {
1353     AliError("Invalid event object");
1354     return;
1355   }
1356   // info from V0 detector QA 
1357   AliESDVZERO * vzero = fESD->GetVZEROData();
1358   Float_t V0Atime = vzero->GetV0ATime();
1359   Float_t V0Ctime = vzero->GetV0CTime(); 
1360   ((TH2F*)fHlistTimeZero->FindObject(Form("hEventV0MeanVsVtx%s",suffix.Data())))->Fill((V0Atime-V0Ctime)*0.5,(V0Atime+V0Ctime)*0.5);
1361   
1362   // info from T0 detector QA 
1363   for (Int_t j=0;j<3;j++){
1364     fT0[j]= (Float_t) fESD->GetT0TOF(j);//ps
1365     if (fT0[j]>90000.) fT0[j]=99999.;//fix old default values to the new one
1366   }
1367   
1368   Float_t t0cut = 90000.; 
1369   //Float_t t0cut =3 * t0spread; //use this cut to check t0 used in tof response
1370   // if(t0cut < 500) t0cut = 500;
1371   
1372   if(TMath::Abs(fT0[1]) < t0cut && TMath::Abs(fT0[2]) < t0cut ) {
1373     //&& TMath::Abs(fT0[2]-fT0[1]) < 500)  //add this condition to check t0 used in tof response
1374     ((TH1F*)fHlistTimeZero->FindObject(Form("hT0DetRes%s",suffix.Data())))->Fill((fT0[2]-fT0[1])*0.5);
1375     ((TH1F*)fHlistTimeZero->FindObject(Form("hT0AC%s",suffix.Data())))->Fill(fT0[0]);  
1376     ((TH2F*)fHlistTimeZero->FindObject(Form("hEventT0MeanVsVtx%s",suffix.Data())))->Fill((fT0[2]-fT0[1])*0.5,(fT0[2]+fT0[1])*0.5);
1377   } 
1378   if(TMath::Abs(fT0[1]) < t0cut){
1379     ((TH1F*)fHlistTimeZero->FindObject(Form("hT0A%s",suffix.Data())))->Fill(fT0[1]);   
1380   }
1381   if(TMath::Abs(fT0[2]) < t0cut){
1382     ((TH1F*)fHlistTimeZero->FindObject(Form("hT0C%s",suffix.Data())))->Fill(fT0[2]);
1383   }
1384   
1385   //  event timeZero QA via AliESDpid::SetTOFResponse() 
1386   Double_t timeZero[4]={99999.,99999.,99999.,99999.};
1387   Double_t timeZeroRes[4]={99999.,99999.,99999.,99999.}; 
1388   
1389   TString timeZeroHisto[4]={"hT0fill","hT0TOF","hT0T0","hT0best"};
1390   TString timeZeroHistoRes[4]={"hT0fillRes","hT0TOFRes","hT0T0Res","hT0bestRes"};
1391   for (Int_t j=0;j<4;j++){
1392     timeZeroHisto[j].Append(suffix.Data());
1393     timeZeroHistoRes[j].Append(suffix.Data());
1394     fESDpid->SetTOFResponse(fESD, (AliESDpid::EStartTimeType_t) j);//(fill_t0, tof_t0, t0_t0, best_t0)
1395     timeZero[j]=fESDpid->GetTOFResponse().GetStartTime(10.); //timeZero for bin pT>10GeV/c
1396     timeZeroRes[j]=fESDpid->GetTOFResponse().GetStartTimeRes(10.); //timeZero for bin pT>10GeV/c
1397     ((TH1D*)(fHlistTimeZero->FindObject(timeZeroHisto[j].Data())))->Fill(timeZero[j]);
1398     ((TH1D*)(fHlistTimeZero->FindObject(timeZeroHistoRes[j].Data())))->Fill(timeZeroRes[j]);
1399   }
1400   ((TH2F*)fHlistTimeZero->FindObject("hT0TOFvsNtrk"))->Fill(fNTOFtracks[0],timeZero[AliESDpid::kTOF_T0]);
1401
1402   //response set to best_t0 by previous loop
1403   FillStartTimeMaskHisto(suffix.Data());
1404
1405   return;
1406 }
1407 //----------------------------------------------------------------------------------
1408 void AliAnalysisTaskTOFqaID::FillTrdHisto(AliESDtrack * track, Int_t charge)
1409 {
1410   //fill histograms for TRD/noTRD
1411   if (!track){
1412     AliError("Invalid track object");
1413     return;
1414   }
1415   
1416   if (IsInTRD(track)){
1417     FillPrimaryTrkHisto(charge,"_Trd");
1418     if (IsTPCTOFMatched(track)) {      
1419       FillMatchedTrkHisto(charge,"_Trd");
1420       FillPidHisto(track,charge, "_Trd");
1421     }
1422   } else {
1423     FillPrimaryTrkHisto(charge,"_noTrd");    
1424     if (IsTPCTOFMatched(track)) {      
1425       FillMatchedTrkHisto(charge,"_noTrd");
1426       FillPidHisto(track, charge, "_noTrd");
1427     }
1428   }
1429   return;
1430 }
1431 //----------------------------------------------------------------------------------
1432 void AliAnalysisTaskTOFqaID::FillTofTrgHisto(TString suffix)
1433 {
1434   //fills histo with trigger info
1435   if (!fHlistTrigger){
1436     AliError("Invalid TOF trigger list");
1437     return;
1438   }
1439   if (!fTOFHeader) {
1440     AliWarning("Invalid AliTOFHeader object - cannot fill trg histo");
1441     return;    
1442   }
1443   
1444   Int_t nPad = fTOFHeader->GetNumberOfTOFclusters(); // tutti i pad di readout che hanno sparato
1445   Int_t nTrgPad = fTOFHeader->GetNumberOfTOFtrgPads();//i pad di readout che hanno sparato nella finesta di trigger
1446   Int_t nMaxiPad = fTOFHeader->GetNumberOfTOFmaxipad(); // numero di maxipad che hanno sparato
1447   
1448   // aggiornare un histo con le mappe di maxipad che hanno sparato
1449   AliTOFTriggerMask *trgMask = fTOFHeader->GetTriggerMask();
1450   for(Int_t j=0;j<72;j++){
1451     for(Int_t i=22;i>=0;i--){
1452       if(trgMask->IsON(j,i))
1453         ((TH1I*)fHlistTrigger->FindObject(Form("hTrgMap%s",suffix.Data())))->Fill(j+1,i+1); 
1454     }   
1455   }
1456   ((TH1I*)fHlistTrigger->FindObject(Form("hFiredMaxipad%s",suffix.Data())))->Fill(nMaxiPad);
1457   ((TH1I*)fHlistTrigger->FindObject(Form("hFiredReadoutPad%s",suffix.Data())))->Fill(nPad);
1458   ((TH1I*)fHlistTrigger->FindObject(Form("hFiredReadoutTrgPad%s",suffix.Data())))->Fill(nTrgPad);
1459   ((TH2I*)fHlistTrigger->FindObject(Form("hFiredMaxipadVsTrgPad%s",suffix.Data())))->Fill(nTrgPad,nMaxiPad);
1460   
1461   return;
1462 }
1463
1464 #endif