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