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