]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliTOFQATask.cxx
Updates on Lambda_c decays (S. Masciocchi)
[u/mrichter/AliRoot.git] / ESDCheck / AliTOFQATask.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 //_________________________________________________________________________
16 // An analysis task to check the TOF in simulated data
17 //
18 //*-- Silvia Arcelli
19 //-Distributions of the matching performance
20 //-TOF Time info and (TOF - expected time) plots
21 //-Summary Plots on TOF Pid
22 //////////////////////////////////////////////////////////////////////////////
23 #include <TMath.h>
24 #include <TChain.h>
25 #include <TFile.h> 
26 #include <TObject.h> 
27 #include <TCanvas.h>
28 #include <TStyle.h>
29 #include <TGaxis.h>
30 #include <TLatex.h>
31 #include <TPad.h>
32 #include <TSpectrum.h>
33 #include <TF1.h>
34 #include <TROOT.h>
35
36 #include "AliTOFQATask.h" 
37 #include "AliESD.h" 
38 #include "AliESDtrack.h" 
39 #include "AliLog.h"
40
41 //______________________________________________________________________________
42 AliTOFQATask::AliTOFQATask(const char *name) : 
43   AliAnalysisTask(name,""),  
44   fChain(0),
45   fESD(0), 
46   fOutputContainer(0), 
47   fhTOFMatch(0),
48   fhESDeffPhi(0),
49   fhESDeffTheta(0),
50   fhESDeffMom(0),
51   fhTOFeffPhi(0),
52   fhTOFeffTheta(0),
53   fhTOFeffMom(0),
54   fhTOFeffPhiMT(0),
55   fhTOFeffThetaMT(0),
56   fhTOFeffMomMT(0),
57   fhTOFsector(0),
58   fhTOFsectorMT(0),
59   fhTOFTime(0),
60   fhTOFDeltaTime(0),
61   fhTOFDeltaTimeMT(0),
62   fhTOFIDSpecies(0),
63   fhTOFMassVsMom(0),
64   fhTOFMass(0),
65   fmatchFracMin(0.15),
66   fmatchEffMin(0.28),
67   ftimePeakMax(0.05),
68   fmassPeakMax(0.05)
69 {
70   // Constructor.
71   // Input slot #0 works with an Ntuple
72   DefineInput(0, TChain::Class());
73   // Output slot #0 writes into a TH1 container
74   DefineOutput(0,  TObjArray::Class()) ; 
75 }
76
77 //______________________________________________________________________________
78 AliTOFQATask::AliTOFQATask(const AliTOFQATask &qatask) : 
79   AliAnalysisTask("AliTOFQATask",""),  
80   fChain(0),
81   fESD(0), 
82   fOutputContainer(0), 
83   fhTOFMatch(0),
84   fhESDeffPhi(0),
85   fhESDeffTheta(0),
86   fhESDeffMom(0),
87   fhTOFeffPhi(0),
88   fhTOFeffTheta(0),
89   fhTOFeffMom(0),
90   fhTOFeffPhiMT(0),
91   fhTOFeffThetaMT(0),
92   fhTOFeffMomMT(0),
93   fhTOFsector(0),
94   fhTOFsectorMT(0),
95   fhTOFTime(0),
96   fhTOFDeltaTime(0),
97   fhTOFDeltaTimeMT(0),
98   fhTOFIDSpecies(0),
99   fhTOFMassVsMom(0),
100   fhTOFMass(0),
101   fmatchFracMin(0.15),
102   fmatchEffMin(0.28),
103   ftimePeakMax(0.05),
104   fmassPeakMax(0.05)
105 {
106   // Copy Constructor.
107   fChain=qatask.fChain;
108   fESD=qatask.fESD; 
109   fOutputContainer=qatask.fOutputContainer; 
110   fhTOFMatch=qatask.fhTOFMatch;
111   fhESDeffPhi=qatask.fhESDeffPhi;
112   fhESDeffTheta=qatask.fhESDeffTheta;
113   fhESDeffMom=qatask.fhESDeffMom;
114   fhTOFeffPhi=qatask.fhTOFeffPhi;
115   fhTOFeffTheta=qatask.fhTOFeffTheta;
116   fhTOFeffMom=qatask.fhTOFeffMom;
117   fhTOFeffPhiMT=qatask.fhTOFeffPhiMT;
118   fhTOFeffThetaMT=qatask.fhTOFeffThetaMT;
119   fhTOFeffMomMT=qatask.fhTOFeffMomMT;
120   fhTOFsector=qatask.fhTOFsector;
121   fhTOFsectorMT=qatask.fhTOFsectorMT;
122   fhTOFTime=qatask.fhTOFTime;
123   fhTOFDeltaTime=qatask.fhTOFDeltaTime;
124   fhTOFDeltaTimeMT=qatask.fhTOFDeltaTimeMT;
125   fhTOFIDSpecies=qatask.fhTOFIDSpecies;
126   fhTOFMassVsMom=qatask.fhTOFMassVsMom;
127   fhTOFMass=qatask.fhTOFMass;
128   fmatchFracMin=qatask.fmatchFracMin;
129   fmatchEffMin=qatask.fmatchEffMin;
130   ftimePeakMax=qatask.ftimePeakMax;
131   fmassPeakMax=qatask.fmassPeakMax;
132 }
133 //______________________________________________________________________________
134 AliTOFQATask:: ~AliTOFQATask() 
135 {
136   delete fOutputContainer;
137   delete fhTOFMatch;
138   delete fhESDeffPhi;
139   delete fhESDeffTheta;
140   delete fhESDeffMom;
141   delete fhTOFeffPhi;
142   delete fhTOFeffTheta;
143   delete fhTOFeffMom;
144   delete fhTOFeffPhiMT;
145   delete fhTOFeffThetaMT;
146   delete fhTOFeffMomMT;
147   delete fhTOFsector;
148   delete fhTOFsectorMT;
149   delete fhTOFTime;
150   delete fhTOFDeltaTime;
151   delete fhTOFDeltaTimeMT;
152   delete fhTOFIDSpecies;
153   delete fhTOFMassVsMom;
154   delete fhTOFMass;
155   }
156 //______________________________________________________________________________
157 AliTOFQATask& AliTOFQATask::operator=(const AliTOFQATask &qatask)  
158
159    //assignment operator
160   this->fChain=qatask.fChain;
161   this->fESD=qatask.fESD; 
162   this->fOutputContainer=qatask.fOutputContainer; 
163   this->fhTOFMatch=qatask.fhTOFMatch;
164   this->fhESDeffPhi=qatask.fhESDeffPhi;
165   this->fhESDeffTheta=qatask.fhESDeffTheta;
166   this->fhESDeffMom=qatask.fhESDeffMom;
167   this->fhTOFeffPhi=qatask.fhTOFeffPhi;
168   this->fhTOFeffTheta=qatask.fhTOFeffTheta;
169   this->fhTOFeffMom=qatask.fhTOFeffMom;
170   this->fhTOFeffPhiMT=qatask.fhTOFeffPhiMT;
171   this->fhTOFeffThetaMT=qatask.fhTOFeffThetaMT;
172   this->fhTOFeffMomMT=qatask.fhTOFeffMomMT;
173   this->fhTOFsector=qatask.fhTOFsector;
174   this->fhTOFsectorMT=qatask.fhTOFsectorMT;
175   this->fhTOFTime=qatask.fhTOFTime;
176   this->fhTOFDeltaTime=qatask.fhTOFDeltaTime;
177   this->fhTOFDeltaTimeMT=qatask.fhTOFDeltaTimeMT;
178   this->fhTOFIDSpecies=qatask.fhTOFIDSpecies;
179   this->fhTOFMassVsMom=qatask.fhTOFMassVsMom;
180   this->fhTOFMass=qatask.fhTOFMass;
181   this->fmatchFracMin=qatask.fmatchFracMin;
182   this->fmatchEffMin=qatask.fmatchEffMin;
183   this->ftimePeakMax=qatask.ftimePeakMax;
184   this->fmassPeakMax=qatask.fmassPeakMax;
185   return *this;
186 }
187 //______________________________________________________________________________
188 void AliTOFQATask::ConnectInputData(const Option_t*)
189 {
190   // Initialisation of branch container and histograms 
191     
192   AliInfo(Form("*** Initialization of %s", GetName())) ; 
193   
194   // Get input data
195   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
196   if (!fChain) {
197     AliError(Form("Input 0 for %s not found\n", GetName()));
198     return ;
199   }
200   
201   // One should first check if the branch address was taken by some other task
202   char ** address = (char **)GetBranchAddress(0, "ESD") ;
203   if (address) 
204     fESD = (AliESD *)(*address) ; 
205   else {
206     fESD = new AliESD() ; 
207     SetBranchAddress(0, "ESD", &fESD) ;  
208   }
209 }
210
211 //______________________________________________________________________________
212 void AliTOFQATask::Exec(Option_t *) 
213 {
214
215 //******* The loop over events --------------------------------------------------
216
217   Int_t nselESD=0;
218   Int_t nmatchTOF=0;
219   Int_t npidok=0;
220   Int_t pisel=0,kasel=0,prsel=0,elsel=0,musel=0;
221   //Set equal a-priori weights (just charged hadrions)
222   const Int_t nCalinSec=8736;
223   Double_t c[5]={0, 0, 1, 1, 1};
224   // Processing of one event
225   Long64_t entry = fChain->GetReadEntry() ;  
226   if (!fESD) {
227     AliError("fESD is not connected to the input!") ; 
228     return ; 
229   }
230   
231   if ( !((entry-1)%100) ) 
232     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
233   
234   // ************************  TOF *************************************
235
236
237   Int_t ntrk = fESD->GetNumberOfTracks() ;
238
239   while ( ntrk-- ) {
240     
241     AliESDtrack * t = fESD->GetTrack(ntrk) ;
242     if ( (t->GetStatus() & AliESDtrack::kTIME)==0 )continue;    
243
244     nselESD++;
245
246     Double_t mom   = t->GetP() ; 
247     Double_t phi   = TMath::ATan2(t->GetX(),t->GetY()) ; 
248     Double_t theta = TMath::ACos(t->GetZ()/
249     TMath::Sqrt(t->GetX()*t->GetX()+t->GetY()*t->GetY()+t->GetZ()*t->GetZ())); 
250     phi*=180/TMath::Pi();
251     theta*=180/TMath::Pi();
252
253     fhESDeffPhi->Fill(phi);
254     fhESDeffTheta->Fill(theta);
255     fhESDeffMom->Fill(mom);
256
257
258     if(t->GetTOFsignal()<0)continue;
259
260     nmatchTOF++;
261
262     Double_t time=t->GetTOFsignal();//TOF time in ps    
263     Int_t detid=t->GetTOFCalChannel();//which pad was hit
264     Int_t sector = detid/nCalinSec;
265     fhTOFTime->Fill(time*1.E-3);
266     fhTOFeffPhi->Fill(phi);
267     fhTOFeffTheta->Fill(theta);
268     fhTOFeffMom->Fill(mom);
269     fhTOFsector->Fill(sector);
270     //Well matched
271
272     Int_t label=TMath::Abs(t->GetLabel());
273     Int_t clab[3]; t->GetTOFLabel(clab);    
274     if(label==clab[0] || label==clab[1] || label==clab[2]) {
275       fhTOFeffPhiMT->Fill(phi);
276       fhTOFeffThetaMT->Fill(theta);
277       fhTOFeffMomMT->Fill(mom);
278       fhTOFsectorMT->Fill(sector);
279     }           
280
281     //Look at TOF PID 
282
283     UInt_t status=AliESDtrack::kESDpid;status|=AliESDtrack::kTOFpid; 
284     if (!((t->GetStatus()&status) == status))continue;
285     npidok++;
286     Double_t times[10];
287     t->GetIntegratedTimes(times);//ps
288     Double_t l   =t->GetIntegratedLength()/100.; // (m)
289     Double_t mass= -1.;
290     Double_t invBetaGamma= (0.299*time*1.E-3/l)*(0.299*time*1.E-3/l) -1.;
291     if(invBetaGamma<0){mass = -mom*TMath::Sqrt(-invBetaGamma);}
292     else{mass = mom*TMath::Sqrt(invBetaGamma);}
293
294     //The Mass/ vs Momentum Plot:
295     fhTOFMassVsMom->Fill(mass,mom);
296     fhTOFMass->Fill(mass);
297     
298     //PID weights 
299     Double_t r[10]; t->GetTOFpid(r);
300     Double_t rcc=0.;
301     Int_t i;
302     for (i=0; i<AliPID::kSPECIES; i++) rcc+=(c[i]*r[i]);
303     if (rcc==0.) continue;
304     Double_t w[10];
305     for (i=0; i<AliPID::kSPECIES; i++) w[i]=c[i]*r[i]/rcc;
306     
307     fhTOFDeltaTime->Fill((time-times[2])*1.E-3);
308
309     if(label==clab[0] || label==clab[1] || label==clab[2]) {
310       fhTOFDeltaTimeMT->Fill((time-times[2])*1.E-3);
311     }           
312     
313     if (w[4]>w[0] && w[4]>w[1] && w[4]>w[2] && w[4]>w[3]){
314       prsel++;
315       fhTOFIDSpecies->Fill(4);
316     }
317     if (w[3]>w[0] && w[3]>w[1] && w[3]>w[2] && w[3]>w[4]){
318       kasel++;
319       fhTOFIDSpecies->Fill(3);
320     }
321     if (w[2]>w[0] && w[2]>w[1] && w[2]>w[3] && w[2]>w[4]){
322       pisel++;
323       fhTOFIDSpecies->Fill(2);
324     }
325     if (w[1]>w[0] && w[1]>w[2] && w[1]>w[3] && w[1]>w[4]){
326       musel++;
327       fhTOFIDSpecies->Fill(1);
328     }
329     if (w[0]>w[1] && w[0]>w[2] && w[0]>w[3] && w[0]>w[4]){
330       elsel++;   
331       fhTOFIDSpecies->Fill(0);
332     }
333   }
334
335   Float_t fracM= -1;
336   if(nselESD>10)fracM=((Float_t) nmatchTOF)/((Float_t) nselESD);
337   fhTOFMatch->Fill(fracM);
338
339   PostData(0, fOutputContainer);  
340
341 }
342 //______________________________________________________________________________
343 void AliTOFQATask::CreateOutputObjects()
344 {  
345   // Construct histograms:
346   
347   OpenFile(0) ; 
348
349   fhTOFMatch= 
350     new TH1F("hTOFMatch","Fraction of Matched TOF tracks",101,-0.005,1.005);
351   fhESDeffPhi= 
352     new TH1F("hESDeffPhi","ESD tracks Phi(vtx)",    180, -180., 180.) ;
353   fhESDeffTheta= 
354     new TH1F("hESDeffTheta","ESD tracks Theta (vtx)",90, 45., 135.) ;
355   fhESDeffMom= 
356     new TH1F("hESDeffMom","ESD tracks Momentum (vtx)",40, 0., 6.) ;
357   fhTOFeffPhi = 
358     new TH1F("hTOFeffPhi","TOF, Matched vs Phi(vtx)", 180,-180, 180.);
359   fhTOFeffPhiMT= 
360     new TH1F("hTOFeffPhiMT","TOF, Well Matched vs Phi(vtx)",180,-180,180.);
361   fhTOFeffTheta= 
362     new TH1F("hTOFeffTheta","TOF, Matched vs Theta(vtx)",90,45.,135.);  
363   fhTOFeffThetaMT= 
364     new TH1F("hTOFeffThetaMT","TOF, Well Matched vs Theta(vtx)",90,45.,135.);
365   fhTOFeffMom = 
366     new TH1F("hTOFeffMom","TOF, Matched vs Momentum ", 40, 0.,6.);
367   fhTOFeffMomMT = 
368     new TH1F("hTOFeffMomMT","TOF, Well Matched vs Momentum ", 40, 0.,6.); 
369   fhTOFsector = 
370     new TH1F("hTOFsector","TOF, Matched vs Sector ", 18,0.,18.); 
371   fhTOFsectorMT = 
372     new TH1F("hTOFsectorMT","TOF, Well Matched vs Sector", 18, 0.,18.); 
373
374   fhESDeffMom->Sumw2(); fhTOFeffMom->Sumw2();  fhTOFeffMomMT->Sumw2();
375   fhESDeffTheta->Sumw2();  fhTOFeffTheta->Sumw2();  fhTOFeffThetaMT->Sumw2();
376   fhESDeffPhi->Sumw2();  fhTOFeffPhi->Sumw2();  fhTOFeffPhiMT->Sumw2();
377   fhTOFsector->Sumw2(); fhTOFsectorMT->Sumw2();
378
379   fhTOFTime = 
380     new TH1F("hTOFTime","TOF, t(TOF)in ns ",1000,0,100.); 
381   fhTOFDeltaTime = 
382     new TH1F("hTOFDeltaTime","TOF,t(TOF)-t(exp,pion), ns ",1000,-4.4,20.);   
383   fhTOFDeltaTimeMT = 
384     new TH1F("hTOFDeltaTimeMT","TOF, t(TOF)-t(exp,pion) for Well Matched, ns ",1000,-4.4,20.); 
385   fhTOFIDSpecies = 
386     new TH1F("hTOFIDSpecies","TOF, ID Sample Composition ",5,-0.5,4.5);
387   fhTOFMassVsMom = 
388     new TH2F("hTOFMassVsMom","TOF, Mass Vs Momentum ",280,-0.2,1.2,600,0.,6.);
389   fhTOFMass = 
390     new TH1F("hTOFMass","TOF, reconstructed mass ",240,0.,1.2);
391
392   
393   // create the output container
394   
395   fOutputContainer = new TObjArray(18) ; 
396   fOutputContainer->SetName(GetName()) ; 
397
398   fOutputContainer->AddAt(fhTOFMatch,             0) ; 
399   fOutputContainer->AddAt(fhESDeffPhi,            1) ; 
400   fOutputContainer->AddAt(fhESDeffTheta,          2) ; 
401   fOutputContainer->AddAt(fhESDeffMom,            3) ; 
402   fOutputContainer->AddAt(fhTOFeffPhi,            4) ; 
403   fOutputContainer->AddAt(fhTOFeffPhiMT,          5) ; 
404   fOutputContainer->AddAt(fhTOFeffTheta,          6) ; 
405   fOutputContainer->AddAt(fhTOFeffThetaMT,        7) ; 
406   fOutputContainer->AddAt(fhTOFeffMom,            8) ; 
407   fOutputContainer->AddAt(fhTOFeffMomMT,          9) ; 
408   fOutputContainer->AddAt(fhTOFsector,           10) ; 
409   fOutputContainer->AddAt(fhTOFsectorMT,         11) ; 
410   fOutputContainer->AddAt(fhTOFTime,             12) ; 
411   fOutputContainer->AddAt(fhTOFDeltaTime ,       13) ; 
412   fOutputContainer->AddAt(fhTOFDeltaTimeMT,      14) ; 
413   fOutputContainer->AddAt(fhTOFIDSpecies,        15) ; 
414   fOutputContainer->AddAt(fhTOFMassVsMom,        16) ; 
415   fOutputContainer->AddAt(fhTOFMass,             17) ; 
416   
417 }
418 //______________________________________________________________________________
419 void AliTOFQATask::GetEfficiency() 
420 {
421   // calculates the efficiency
422   
423   
424   fhTOFeffPhiMT->Divide(fhTOFeffPhiMT,fhTOFeffPhi,1,1,"B");
425   fhTOFeffThetaMT->Divide(fhTOFeffThetaMT,fhTOFeffTheta,1,1,"B");
426   fhTOFeffMomMT->Divide(fhTOFeffMomMT,fhTOFeffMom,1,1,"B");
427
428   fhTOFeffPhi->Divide(fhTOFeffPhi,fhESDeffPhi,1,1,"B");
429   fhTOFeffTheta->Divide(fhTOFeffTheta,fhESDeffTheta,1,1,"B");
430   fhTOFeffMom->Divide(fhTOFeffMom,fhESDeffMom,1,1,"B");
431   
432 }
433
434 //______________________________________________________________________________
435 Bool_t AliTOFQATask::DrawHistos() 
436 {
437
438
439   Int_t TOFsectors[18]={0,1,1,0,0,0,1,1,0,1,1,1,1,0,0,1,1,1};//TOF sectors which are supposed to be present
440   const char* part[3]={"pions","kaons","protons"};// pi,ka,pr
441   const Float_t masses[3]={0.1396,0.494,0.938};// particle masses pi,ka,pr
442   // Makes a few plots
443
444   gROOT->SetStyle("Plain");
445   gStyle->SetPalette(1);
446   gStyle->SetOptStat(111110);
447   gStyle->SetPalette(1);
448   gStyle->SetCanvasColor(0);
449   gStyle->SetFrameFillColor(0);
450   //
451   TGaxis::SetMaxDigits(3);  
452   gStyle->SetLabelFont(52, "XYZ");
453   gStyle->SetTitleFont(62, "XYZ");
454   gStyle->SetPadRightMargin(0.02);
455
456   
457   TCanvas * cTOFeff = new TCanvas("cTOFeff", "TOF ESD General", 400, 30, 550, 630) ;
458
459   cTOFeff->Divide(1,2) ;
460   cTOFeff->cd(1);
461   fhTOFMatch->GetXaxis()->SetTitle("Fraction of matched ESD tracks/event");
462   fhTOFMatch->GetYaxis()->SetTitle("Events");
463   fhTOFMatch->SetFillColor(4);
464   fhTOFMatch->Draw();  
465   cTOFeff->cd(2);
466   TPad *b = (TPad*)gPad;
467   b->Divide(2,1);
468   b->cd(1) ;
469   fhTOFsector->GetXaxis()->SetTitle("Sector index of matched TOF Cluster");
470   fhTOFsector->GetYaxis()->SetTitle("Entries");
471   fhTOFsector->SetFillColor(4);
472   fhTOFsector->SetMinimum(0.);
473   fhTOFsector->GetYaxis()->SetNdivisions(205);
474   fhTOFsector->GetYaxis()->SetTitleOffset(1.2);
475   fhTOFsector->Draw("histo");
476   b->cd(2) ;
477   fhTOFeffMom->SetMaximum(1.2);
478   fhTOFeffMom->GetXaxis()->SetTitle("Track momentum (GeV/c)");
479   fhTOFeffMom->GetYaxis()->SetNdivisions(205);
480   fhTOFeffMom->GetYaxis()->SetTitleOffset(1.2);
481   fhTOFeffMom->GetYaxis()->SetTitle("Fraction of matched ESD tracks");
482   fhTOFeffMom->SetFillColor(4);
483   fhTOFeffMom->Draw();
484   fhTOFeffMom->Fit("pol0","Q","",0.5,3.);
485   TF1 *fitMom = fhTOFeffMom->GetFunction("pol0");
486   fhTOFeffMom->Draw("histo,same");
487   
488   cTOFeff->Print("TOF_eff.eps");
489
490
491
492   TCanvas * cTOFtime = new TCanvas("cTOFtime", "TOF measured Times ", 400, 30, 550, 630) ;  
493   cTOFtime->Divide(1,2) ;
494   cTOFtime->cd(1);
495   if ( fhTOFTime->GetMaximum() > 0 ) 
496     cTOFtime->GetPad(1)->SetLogy(1);
497   fhTOFTime->GetXaxis()->SetTitle("TOF time (ns)");
498   fhTOFTime->GetYaxis()->SetTitle("Entries");
499   fhTOFTime->SetFillColor(4);
500   fhTOFTime->Draw();
501   cTOFtime->cd(2);
502   if ( fhTOFDeltaTime->GetMaximum() > 0 ) 
503     cTOFtime->GetPad(2)->SetLogy(1);
504   fhTOFDeltaTime->GetXaxis()->SetTitle("t^{TOF}-t^{exp}_{#pi} (ns)");
505   fhTOFDeltaTime->GetYaxis()->SetTitle("Entries");
506   fhTOFDeltaTime->SetFillColor(4);
507   Int_t ntimepeak=1;
508   TSpectrum *timeDiff = new TSpectrum(ntimepeak);
509   Int_t ntime = timeDiff->Search(fhTOFDeltaTime,ntimepeak,"new",0.1);
510   fhTOFDeltaTime->Draw();
511
512   cTOFtime->Print("TOF_time.eps");
513
514
515   TCanvas * cTOFpid = new TCanvas("cTOFpid", "TOF PID ", 400, 30, 550, 630) ;  
516   cTOFpid->Divide(1,3) ;
517   cTOFpid->cd(1);
518   if ( fhTOFMass->GetMaximum() > 0 ) 
519     cTOFpid->SetLogy(1);
520   fhTOFMass->GetXaxis()->SetTitle("Reconstructed Mass (GeV/c^{2})");
521   fhTOFMass->GetYaxis()->SetTitle("Entries");
522   fhTOFMass->SetFillColor(4);
523   Int_t npmass=1;
524   if(fhTOFMass->GetEntries()>1000)npmass=3;
525   TSpectrum *mass = new TSpectrum(npmass);
526   Int_t nmass = mass->Search(fhTOFMass,npmass,"new",0.02);
527   fhTOFMass->Draw();
528   if ( fhTOFMassVsMom->GetMaximum() > 0 ) 
529     cTOFpid->SetLogy(0);
530   cTOFpid->cd(2);
531   fhTOFMassVsMom->GetYaxis()->SetRange(0,400);
532   fhTOFMassVsMom->GetXaxis()->SetTitle("Reconstructed Mass (GeV/c^{2})");
533   fhTOFMassVsMom->GetYaxis()->SetTitle("Track Momentum (GeV/c)");
534   fhTOFMassVsMom->GetXaxis()->SetTitleSize(0.05);
535   fhTOFMassVsMom->GetYaxis()->SetTitleSize(0.05);
536   fhTOFMassVsMom->SetMarkerStyle(20);
537   fhTOFMassVsMom->SetMarkerSize(0.05);
538   fhTOFMassVsMom->SetMarkerColor(2);
539   fhTOFMassVsMom->Draw();
540   cTOFpid->cd(3);
541
542   TLatex *   tex = new TLatex(1., 1.25, "Bayesian PID: a-priori concentrations: [0,0,1,1,1]");
543   tex->SetTextColor(1);
544   tex->SetTextSize(0.045);
545   tex->SetLineWidth(2);
546
547   Float_t norm=1./fhTOFIDSpecies->GetEntries();
548   fhTOFIDSpecies->Scale(norm);
549   fhTOFIDSpecies->SetMaximum(1.2);
550   fhTOFIDSpecies->GetXaxis()->SetTitle("Particle Type");
551   fhTOFIDSpecies->GetYaxis()->SetTitle("ID Fractions");
552   fhTOFIDSpecies->GetXaxis()->SetTitleSize(0.05);
553   fhTOFIDSpecies->GetYaxis()->SetTitleSize(0.05);
554   fhTOFIDSpecies->SetFillColor(4);
555   fhTOFIDSpecies->Draw();
556   tex->Draw();
557  
558   char ch[10];
559
560   Float_t pifrac=fhTOFIDSpecies->GetBinContent(3);
561   Float_t kafrac=fhTOFIDSpecies->GetBinContent(4);
562   Float_t prfrac=fhTOFIDSpecies->GetBinContent(5);
563
564   sprintf(ch,"[2]:pion fraction   = %5.3f",pifrac);    
565   TLatex *   texpi = new TLatex(-0.3, 0.9, ch);
566   texpi->SetTextColor(1);
567   texpi->SetTextSize(0.05);
568   texpi->SetLineWidth(2);
569   texpi->Draw();
570   sprintf(ch,"[3]:kaon fraction   = %5.3f",kafrac);    
571   TLatex *   texka = new TLatex(-0.3, 0.8, ch);
572   texka->SetTextColor(1);
573   texka->SetTextSize(0.05);
574   texka->SetLineWidth(2);
575   texka->Draw();
576   sprintf(ch,"[4]:proton fraction = %5.3f",prfrac);    
577   TLatex *   texpr = new TLatex(-0.3, 0.7, ch);
578   texpr->SetTextColor(1);
579   texpr->SetTextSize(0.05);
580   texpr->SetLineWidth(2);
581   texpr->Draw();
582
583
584   cTOFpid->Print("TOF_pid.eps");
585
586   char line[1024] ; 
587   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
588   gROOT->ProcessLine(line);
589   sprintf(line, ".!rm -fR *.eps"); 
590   gROOT->ProcessLine(line);
591   AliInfo(Form("*** TOF QA plots saved in %s.tar.gz...", GetName())) ;
592
593   Bool_t hardProblem=kFALSE;
594   Bool_t softProblem=kFALSE;
595   //------------------------------Matching Efficiency
596
597   //Overall Fraction:
598   Float_t matchFrac= fhTOFMatch->GetMean();
599   if(matchFrac<fmatchFracMin){
600     AliWarning(Form("*** Overall Fraction of matched tracks too low! Fraction = %f", matchFrac)) ;
601     softProblem=kTRUE;
602   }else{
603     AliInfo(Form("*** Fraction of matched tracks  = %f", matchFrac)) ; 
604   } 
605
606   if(fhTOFeffMom->GetEntries()<1.){
607     AliWarning(Form("*** No tracks matching with TOF! Fraction is = %f", matchFrac)) ; 
608     hardProblem=kTRUE;
609     return hardProblem;
610   }
611
612   
613   //The efficiency as a function of momentum:
614   Float_t eff=  fitMom->GetParameter(0);
615   Float_t deff=  fitMom->GetParError(0);
616   if(eff+3*deff<fmatchEffMin){
617     AliWarning(Form("*** Fraction of matched tracks vs Momentum is too low! Fraction= %f", eff)) ;
618     softProblem=kTRUE;
619   }else{
620     AliInfo(Form("*** Fraction of matched tracks for p>0.5 GeV is = %f", eff)) ;
621   } 
622
623   //Matched tracks vs TOF Sector:
624   for(Int_t isec=1;isec<=18;isec++){
625     if(fhTOFsector->GetBinContent(isec)<1 && TOFsectors[isec-1]>0){
626       AliWarning(Form("*** Missing Entries in sector %i", isec)); 
627       softProblem=kTRUE;
628     }
629     if(fhTOFsector->GetBinContent(isec)>0 && TOFsectors[isec-1]==0){
630       AliWarning(Form("*** Unexpected Entries in sector %i", isec)); 
631       softProblem=kTRUE;
632     }
633   }
634
635   //-----------------------------Pid Quality
636
637   // Look at the time - expected time: 
638   if(ntime==0){
639     AliWarning("*** No peak was found in time difference spectrum!");
640   }
641   else if (ntime>0){
642     Float_t *timePos = timeDiff->GetPositionX();
643     if(TMath::Abs(timePos[0])>3*ftimePeakMax){
644       AliWarning(Form("*** Main Peak position in tTOF-TEXP spectrum is sitting very far from where expected! Tpeak = %f ps",timePos[0]*1.E3));  
645       hardProblem=kTRUE;
646     }
647     else if(TMath::Abs(timePos[0])>ftimePeakMax){
648       AliWarning(Form("*** Main Peak position in tTOF-TEXP spectrum is sitting far from where expected! Tpeak = %f ps",timePos[0]*1.E3));  
649       softProblem=kTRUE;
650     }else{
651       AliInfo(Form("*** Main Peak position in tTOF-TEXP found at = %f ps",timePos[0]*1.E3));  
652     }
653   }
654
655   // Look at the Mass Spectrum: 
656   if(nmass==0){
657     AliWarning("*** No peak was found in Mass difference spectrum!");
658     softProblem=kTRUE;
659   }
660   else if(nmass>0){
661     Int_t massind[3];
662     Float_t *massPos = mass->GetPositionX();
663     //check the found peaks (sorted in ascending mass order)
664     TMath::Sort(nmass,massPos,massind,kFALSE);
665     for(Int_t imass=0;imass<TMath::Min(nmass,3);imass++){   
666       AliInfo(Form("*** the Mass peak for %s found at  = %f GeV/c^2",part[imass],massPos[massind[imass]]));
667       if(TMath::Abs( massPos[massind[imass]]-masses[imass])> fmassPeakMax){
668         AliWarning(Form("*** the Mass peak position for %s is not in the right place, found at  = %f GeV/c^2",part[imass],massPos[massind[imass]]));
669         softProblem=kTRUE;
670       }
671     }
672     // harder check on the pion peak (more statistically significant)
673     if(TMath::Abs( massPos[massind[0]]-masses[0])> 3*fmassPeakMax){
674       hardProblem=kTRUE;
675     }
676   }
677   // Look at the ID Species: 
678
679   if(fhTOFIDSpecies->GetEntries()>1000){
680     if(pifrac<0.8 || (kafrac<0.01 || kafrac>0.2) || (prfrac<0.01 || prfrac>0.2)){
681       AliWarning(Form("*** Unexpected Id fractions: pions = %f, kaons = %f, protons %f", pifrac,kafrac,prfrac));
682       softProblem=kTRUE;
683     }
684   }
685
686   delete mass;
687   delete timeDiff;
688   return hardProblem || softProblem ; 
689 }
690
691 //______________________________________________________________________________
692 void AliTOFQATask::Terminate(Option_t *)
693 {
694   // Processing when the event loop is ended, some plots+checks
695   
696   AliInfo(Form("*** %s Report:", GetName())) ; 
697   GetEfficiency();
698   Bool_t problem = DrawHistos() ; 
699
700   char * report ; 
701   if(problem)
702     report="Problems found, please check!!!";  
703   else 
704     report="OK";
705
706   AliInfo(Form("*** %s Summary Report: %s\n",GetName(), report)) ; 
707   
708 }