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