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