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