]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliTOFQATask.cxx
added copyright note
[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 <TChain.h>
24 #include <TFile.h> 
25 #include <TObject.h> 
26 #include <TCanvas.h>
27 #include <TStyle.h>
28 #include <TGaxis.h>
29 #include <TLatex.h>
30 #include <TPad.h>
31
32 #include "AliTOFQATask.h" 
33 #include "AliESD.h" 
34 #include "AliESDtrack.h" 
35 #include "AliLog.h"
36 #include "Riostream.h"
37
38 //______________________________________________________________________________
39 AliTOFQATask::AliTOFQATask(const char *name) : 
40   AliAnalysisTask(name,""),  
41   fChain(0),
42   fESD(0), 
43   fOutputContainer(0), 
44   fhTOFMatch(0),
45   fhESDeffPhi(0),
46   fhESDeffTheta(0),
47   fhESDeffMom(0),
48   fhTOFeffPhi(0),
49   fhTOFeffTheta(0),
50   fhTOFeffMom(0),
51   fhTOFeffPhiMT(0),
52   fhTOFeffThetaMT(0),
53   fhTOFeffMomMT(0),
54   fhTOFsector(0),
55   fhTOFsectorMT(0),
56   fhTOFTime(0),
57   fhTOFDeltaTime(0),
58   fhTOFDeltaTimeMT(0),
59   fhTOFIDSpecies(0),
60   fhTOFMassVsMom(0),
61   fhTOFMass(0)
62 {
63   // Constructor.
64   // Input slot #0 works with an Ntuple
65   DefineInput(0, TChain::Class());
66   // Output slot #0 writes into a TH1 container
67   DefineOutput(0,  TObjArray::Class()) ; 
68 }
69
70 //______________________________________________________________________________
71 AliTOFQATask::AliTOFQATask(const AliTOFQATask &qatask) : 
72   AliAnalysisTask("AliTOFQATask",""),  
73   fChain(0),
74   fESD(0), 
75   fOutputContainer(0), 
76   fhTOFMatch(0),
77   fhESDeffPhi(0),
78   fhESDeffTheta(0),
79   fhESDeffMom(0),
80   fhTOFeffPhi(0),
81   fhTOFeffTheta(0),
82   fhTOFeffMom(0),
83   fhTOFeffPhiMT(0),
84   fhTOFeffThetaMT(0),
85   fhTOFeffMomMT(0),
86   fhTOFsector(0),
87   fhTOFsectorMT(0),
88   fhTOFTime(0),
89   fhTOFDeltaTime(0),
90   fhTOFDeltaTimeMT(0),
91   fhTOFIDSpecies(0),
92   fhTOFMassVsMom(0),
93   fhTOFMass(0)
94 {
95   // Copy Constructor.
96   fChain=qatask.fChain;
97   fESD=qatask.fESD; 
98   fOutputContainer=qatask.fOutputContainer; 
99   fhTOFMatch=qatask.fhTOFMatch;
100   fhESDeffPhi=qatask.fhESDeffPhi;
101   fhESDeffTheta=qatask.fhESDeffTheta;
102   fhESDeffMom=qatask.fhESDeffMom;
103   fhTOFeffPhi=qatask.fhTOFeffPhi;
104   fhTOFeffTheta=qatask.fhTOFeffTheta;
105   fhTOFeffMom=qatask.fhTOFeffMom;
106   fhTOFeffPhiMT=qatask.fhTOFeffPhiMT;
107   fhTOFeffThetaMT=qatask.fhTOFeffThetaMT;
108   fhTOFeffMomMT=qatask.fhTOFeffMomMT;
109   fhTOFsector=qatask.fhTOFsector;
110   fhTOFsectorMT=qatask.fhTOFsectorMT;
111   fhTOFTime=qatask.fhTOFTime;
112   fhTOFDeltaTime=qatask.fhTOFDeltaTime;
113   fhTOFDeltaTimeMT=qatask.fhTOFDeltaTimeMT;
114   fhTOFIDSpecies=qatask.fhTOFIDSpecies;
115   fhTOFMassVsMom=qatask.fhTOFMassVsMom;
116   fhTOFMass=qatask.fhTOFMass;
117 }
118 //______________________________________________________________________________
119 AliTOFQATask:: ~AliTOFQATask() 
120 {
121   delete fOutputContainer;
122   delete fhTOFMatch;
123   delete fhESDeffPhi;
124   delete fhESDeffTheta;
125   delete fhESDeffMom;
126   delete fhTOFeffPhi;
127   delete fhTOFeffTheta;
128   delete fhTOFeffMom;
129   delete fhTOFeffPhiMT;
130   delete fhTOFeffThetaMT;
131   delete fhTOFeffMomMT;
132   delete fhTOFsector;
133   delete fhTOFsectorMT;
134   delete fhTOFTime;
135   delete fhTOFDeltaTime;
136   delete fhTOFDeltaTimeMT;
137   delete fhTOFIDSpecies;
138   delete fhTOFMassVsMom;
139   delete fhTOFMass;
140   }
141 //______________________________________________________________________________
142 AliTOFQATask& AliTOFQATask::operator=(const AliTOFQATask &qatask)  
143
144    //assignment operator
145   this->fChain=qatask.fChain;
146   this->fESD=qatask.fESD; 
147   this->fOutputContainer=qatask.fOutputContainer; 
148   this->fhTOFMatch=qatask.fhTOFMatch;
149   this->fhESDeffPhi=qatask.fhESDeffPhi;
150   this->fhESDeffTheta=qatask.fhESDeffTheta;
151   this->fhESDeffMom=qatask.fhESDeffMom;
152   this->fhTOFeffPhi=qatask.fhTOFeffPhi;
153   this->fhTOFeffTheta=qatask.fhTOFeffTheta;
154   this->fhTOFeffMom=qatask.fhTOFeffMom;
155   this->fhTOFeffPhiMT=qatask.fhTOFeffPhiMT;
156   this->fhTOFeffThetaMT=qatask.fhTOFeffThetaMT;
157   this->fhTOFeffMomMT=qatask.fhTOFeffMomMT;
158   this->fhTOFsector=qatask.fhTOFsector;
159   this->fhTOFsectorMT=qatask.fhTOFsectorMT;
160   this->fhTOFTime=qatask.fhTOFTime;
161   this->fhTOFDeltaTime=qatask.fhTOFDeltaTime;
162   this->fhTOFDeltaTimeMT=qatask.fhTOFDeltaTimeMT;
163   this->fhTOFIDSpecies=qatask.fhTOFIDSpecies;
164   this->fhTOFMassVsMom=qatask.fhTOFMassVsMom;
165   this->fhTOFMass=qatask.fhTOFMass;
166   return *this;
167 }
168 //______________________________________________________________________________
169 void AliTOFQATask::ConnectInputData(const Option_t*)
170 {
171   // Initialisation of branch container and histograms 
172     
173   AliInfo(Form("*** Initialization of %s", GetName())) ; 
174   
175   // Get input data
176   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
177   if (!fChain) {
178     AliError(Form("Input 0 for %s not found\n", GetName()));
179     return ;
180   }
181   
182   // One should first check if the branch address was taken by some other task
183   char ** address = (char **)GetBranchAddress(0, "ESD");
184   if (address) {
185     fESD = (AliESD*)(*address);
186   } else {
187     fESD = new AliESD();
188     SetBranchAddress(0, "ESD", &fESD);
189   }
190 }
191
192 //________________________________________________________________________
193 void AliTOFQATask::CreateOutputObjects()
194 {
195 // Create histograms
196   BookHistos();
197 }
198
199 //______________________________________________________________________________
200 void AliTOFQATask::Exec(Option_t *) 
201 {
202
203 //******* The loop over events --------------------------------------------------
204
205   Int_t nselESD=0;
206   Int_t nmatchTOF=0;
207   Int_t npidok=0;
208   Int_t pisel=0,kasel=0,prsel=0,elsel=0,musel=0;
209   //Set equal a-priori weights (just charged hadrions)
210   Int_t nCalinSec=8736;
211   Double_t c[5]={0, 0, 1, 1, 1};
212   // Processing of one event
213   Long64_t entry = fChain->GetReadEntry() ;  
214   if (!fESD) {
215     AliError("fESD is not connected to the input!") ; 
216     return ; 
217   }
218   
219   if ( !((entry-1)%100) ) 
220     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
221   
222   // ************************  TOF *************************************
223
224
225   Int_t ntrk = fESD->GetNumberOfTracks() ;
226
227   while ( ntrk-- ) {
228     
229     AliESDtrack * t = fESD->GetTrack(ntrk) ;
230     if ( (t->GetStatus() & AliESDtrack::kTIME)==0 )continue;    
231
232     nselESD++;
233
234     Double_t mom   = t->GetP() ; 
235     Double_t phi   = TMath::ATan2(t->GetX(),t->GetY()) ; 
236     Double_t theta = TMath::ACos(t->GetZ()/
237     TMath::Sqrt(t->GetX()*t->GetX()+t->GetY()*t->GetY()+t->GetZ()*t->GetZ())); 
238     phi*=180/TMath::Pi();
239     theta*=180/TMath::Pi();
240
241     fhESDeffPhi->Fill(phi);
242     fhESDeffTheta->Fill(theta);
243     fhESDeffMom->Fill(mom);
244
245
246     if(t->GetTOFsignal()<0)continue;
247
248     nmatchTOF++;
249
250     Double_t time=t->GetTOFsignal();//TOF time in ps    
251     Int_t detid=t->GetTOFCalChannel();//which pad was hit
252     Int_t sector = detid/nCalinSec;
253     fhTOFTime->Fill(time*1.E-3);
254     fhTOFeffPhi->Fill(phi);
255     fhTOFeffTheta->Fill(theta);
256     fhTOFeffMom->Fill(mom);
257     fhTOFsector->Fill(sector);
258     //Well matched
259
260     Int_t label=TMath::Abs(t->GetLabel());
261     Int_t clab[3]; t->GetTOFLabel(clab);    
262     if(label==clab[0] || label==clab[1] || label==clab[2]) {
263       fhTOFeffPhiMT->Fill(phi);
264       fhTOFeffThetaMT->Fill(theta);
265       fhTOFeffMomMT->Fill(mom);
266       fhTOFsectorMT->Fill(sector);
267     }           
268
269     //Look at TOF PID 
270
271     UInt_t status=AliESDtrack::kESDpid;status|=AliESDtrack::kTOFpid; 
272     if (!((t->GetStatus()&status) == status))continue;
273     npidok++;
274     Double_t times[10];
275     t->GetIntegratedTimes(times);//ps
276     Double_t l   =t->GetIntegratedLength()/100.; // (m)
277     Double_t mass= -1.;
278     Double_t invBetaGamma= (0.299*time*1.E-3/l)*(0.299*time*1.E-3/l) -1.;
279     if(invBetaGamma<0){mass = -mom*TMath::Sqrt(-invBetaGamma);}
280     else{mass = mom*TMath::Sqrt(invBetaGamma);}
281
282     //The Mass/ vs Momentum Plot:
283     fhTOFMassVsMom->Fill(mass,mom);
284     fhTOFMass->Fill(mass);
285     
286     //PID weights 
287     Double_t r[10]; t->GetTOFpid(r);
288     Double_t rcc=0.;
289     Int_t i;
290     for (i=0; i<AliPID::kSPECIES; i++) rcc+=(c[i]*r[i]);
291     if (rcc==0.) continue;
292     Double_t w[10];
293     for (i=0; i<AliPID::kSPECIES; i++) w[i]=c[i]*r[i]/rcc;
294     
295     fhTOFDeltaTime->Fill((time-times[2])*1.E-3);
296
297     if(label==clab[0] || label==clab[1] || label==clab[2]) {
298       fhTOFDeltaTimeMT->Fill((time-times[2])*1.E-3);
299     }           
300     
301     if (w[4]>w[0] && w[4]>w[1] && w[4]>w[2] && w[4]>w[3]){
302       prsel++;
303       fhTOFIDSpecies->Fill(4);
304     }
305     if (w[3]>w[0] && w[3]>w[1] && w[3]>w[2] && w[3]>w[4]){
306       kasel++;
307       fhTOFIDSpecies->Fill(3);
308     }
309     if (w[2]>w[0] && w[2]>w[1] && w[2]>w[3] && w[2]>w[4]){
310       pisel++;
311       fhTOFIDSpecies->Fill(2);
312     }
313     if (w[1]>w[0] && w[1]>w[2] && w[1]>w[3] && w[1]>w[4]){
314       musel++;
315       fhTOFIDSpecies->Fill(1);
316     }
317     if (w[0]>w[1] && w[0]>w[2] && w[0]>w[3] && w[0]>w[4]){
318       elsel++;   
319       fhTOFIDSpecies->Fill(0);
320     }
321   }
322
323   Float_t fracM= -1;
324   if(nselESD>10)fracM=((Float_t) nmatchTOF)/((Float_t) nselESD);
325   fhTOFMatch->Fill(fracM);
326
327   PostData(0, fOutputContainer);  
328   AliInfo("Finishing event processing...") ; 
329
330 }
331 //______________________________________________________________________________
332 void AliTOFQATask::BookHistos()
333 {  
334   // Construct histograms:
335   
336   fhTOFMatch= 
337     new TH1F("hTOFMatch","Fraction of Matched TOF tracks",101,-0.005,1.005);
338   fhESDeffPhi= 
339     new TH1F("hESDeffPhi","ESD tracks Phi(vtx)",    180, -180., 180.) ;
340   fhESDeffTheta= 
341     new TH1F("hESDeffTheta","ESD tracks Theta (vtx)",90, 45., 135.) ;
342   fhESDeffMom= 
343     new TH1F("hESDeffMom","ESD tracks Momentum (vtx)",40, 0., 6.) ;
344   fhTOFeffPhi = 
345     new TH1F("hTOFeffPhi","TOF, Matched vs Phi(vtx)", 180,-180, 180.);
346   fhTOFeffPhiMT= 
347     new TH1F("hTOFeffPhiMT","TOF, Well Matched vs Phi(vtx)",180,-180,180.);
348   fhTOFeffTheta= 
349     new TH1F("hTOFeffTheta","TOF, Matched vs Theta(vtx)",90,45.,135.);  
350   fhTOFeffThetaMT= 
351     new TH1F("hTOFeffThetaMT","TOF, Well Matched vs Theta(vtx)",90,45.,135.);
352   fhTOFeffMom = 
353     new TH1F("hTOFeffMom","TOF, Matched vs Momentum ", 40, 0.,6.);
354   fhTOFeffMomMT = 
355     new TH1F("hTOFeffMomMT","TOF, Well Matched vs Momentum ", 40, 0.,6.); 
356   fhTOFsector = 
357     new TH1F("hTOFsector","TOF, Matched vs Sector ", 18,0.,18.); 
358   fhTOFsectorMT = 
359     new TH1F("hTOFsectorMT","TOF, Well Matched vs Sector", 18, 0.,18.); 
360
361   fhESDeffMom->Sumw2(); fhTOFeffMom->Sumw2();  fhTOFeffMomMT->Sumw2();
362   fhESDeffTheta->Sumw2();  fhTOFeffTheta->Sumw2();  fhTOFeffThetaMT->Sumw2();
363   fhESDeffPhi->Sumw2();  fhTOFeffPhi->Sumw2();  fhTOFeffPhiMT->Sumw2();
364   fhTOFsector->Sumw2(); fhTOFsectorMT->Sumw2();
365
366   fhTOFTime = 
367     new TH1F("hTOFTime","TOF, t(TOF)in ns ",1000,0,100.); 
368   fhTOFDeltaTime = 
369     new TH1F("hTOFDeltaTime","TOF,t(TOF)-t(exp,pion), ns ",1000,-4.4,20.);   
370   fhTOFDeltaTimeMT = 
371     new TH1F("hTOFDeltaTimeMT","TOF, t(TOF)-t(exp,pion) for Well Matched, ns ",1000,-4.4,20.); 
372   fhTOFIDSpecies = 
373     new TH1F("hTOFIDSpecies","TOF, ID Sample Composition ",5,-0.5,4.5);
374   fhTOFMassVsMom = 
375     new TH2F("hTOFMassVsMom","TOF, Mass Vs Momentum ",280,-0.2,1.2,600,0.,6.);
376   fhTOFMass = 
377     new TH1F("hTOFMass","TOF, reconstructed mass ",280,-0.2,1.2);
378
379   
380   // create the output container
381   
382   fOutputContainer = new TObjArray(18) ; 
383   fOutputContainer->SetName(GetName()) ; 
384
385   fOutputContainer->AddAt(fhTOFMatch,             0) ; 
386   fOutputContainer->AddAt(fhESDeffPhi,            1) ; 
387   fOutputContainer->AddAt(fhESDeffTheta,          2) ; 
388   fOutputContainer->AddAt(fhESDeffMom,            3) ; 
389   fOutputContainer->AddAt(fhTOFeffPhi,            4) ; 
390   fOutputContainer->AddAt(fhTOFeffPhiMT,          5) ; 
391   fOutputContainer->AddAt(fhTOFeffTheta,          6) ; 
392   fOutputContainer->AddAt(fhTOFeffThetaMT,        7) ; 
393   fOutputContainer->AddAt(fhTOFeffMom,            8) ; 
394   fOutputContainer->AddAt(fhTOFeffMomMT,          9) ; 
395   fOutputContainer->AddAt(fhTOFsector,           10) ; 
396   fOutputContainer->AddAt(fhTOFsectorMT,         11) ; 
397   fOutputContainer->AddAt(fhTOFTime,             12) ; 
398   fOutputContainer->AddAt(fhTOFDeltaTime ,       13) ; 
399   fOutputContainer->AddAt(fhTOFDeltaTimeMT,      14) ; 
400   fOutputContainer->AddAt(fhTOFIDSpecies,        15) ; 
401   fOutputContainer->AddAt(fhTOFMassVsMom,        16) ; 
402   fOutputContainer->AddAt(fhTOFMass,             17) ; 
403   
404 }
405 //______________________________________________________________________________
406 void AliTOFQATask::GetEfficiency() 
407 {
408   // calculates the efficiency
409   
410   
411   fhTOFeffPhiMT->Divide(fhTOFeffPhiMT,fhTOFeffPhi,1,1,"B");
412   fhTOFeffThetaMT->Divide(fhTOFeffThetaMT,fhTOFeffTheta,1,1,"B");
413   fhTOFeffMomMT->Divide(fhTOFeffMomMT,fhTOFeffMom,1,1,"B");
414
415   fhTOFeffPhi->Divide(fhTOFeffPhi,fhESDeffPhi,1,1,"B");
416   fhTOFeffTheta->Divide(fhTOFeffTheta,fhESDeffTheta,1,1,"B");
417   fhTOFeffMom->Divide(fhTOFeffMom,fhESDeffMom,1,1,"B");
418   
419 }
420
421 //______________________________________________________________________________
422 void AliTOFQATask::DrawHistos() 
423 {
424   // Makes a few plots
425
426   AliInfo("Plotting....") ; 
427
428   gROOT->SetStyle("Plain");
429   gStyle->SetPalette(1);
430   gStyle->SetOptStat(110010);
431   gStyle->SetPalette(1);
432   gStyle->SetCanvasColor(0);
433   gStyle->SetFrameFillColor(0);
434   //
435   TGaxis::SetMaxDigits(3);  
436   gStyle->SetLabelFont(52, "XYZ");
437   gStyle->SetTitleFont(62, "XYZ");
438   gStyle->SetPadRightMargin(0.02);
439
440   
441   TCanvas * cTOFeff = new TCanvas("cTOFeff", "TOF ESD General", 400, 30, 550, 630) ;
442
443   cTOFeff->Divide(1,2) ;
444   cTOFeff->cd(1);
445   fhTOFMatch->GetXaxis()->SetTitle("Fraction of matched ESD tracks/event");
446   fhTOFMatch->GetYaxis()->SetTitle("Events");
447   fhTOFMatch->SetFillColor(4);
448   fhTOFMatch->Draw();  
449   cTOFeff->cd(2);
450   TPad *b = (TPad*)gPad;
451   b->Divide(2,1);
452   b->cd(1) ;
453   fhTOFsector->GetXaxis()->SetTitle("Sector index of matched TOF Cluster");
454   fhTOFsector->GetYaxis()->SetTitle("Entries");
455   fhTOFsector->SetFillColor(4);
456   fhTOFsector->SetMinimum(0.);
457   fhTOFsector->GetYaxis()->SetNdivisions(205);
458   fhTOFsector->GetYaxis()->SetTitleOffset(1.2);
459   fhTOFsector->Draw("histo");
460   b->cd(2) ;
461   fhTOFeffMom->SetMaximum(1.2);
462   fhTOFeffMom->GetXaxis()->SetTitle("Track momentum (GeV/c)");
463   fhTOFeffMom->GetYaxis()->SetNdivisions(205);
464   fhTOFeffMom->GetYaxis()->SetTitleOffset(1.2);
465   fhTOFeffMom->GetYaxis()->SetTitle("Fraction of matched ESD tracks");
466   fhTOFeffMom->SetFillColor(4);
467   fhTOFeffMom->Draw();
468   fhTOFeffMom->Draw("histo,same");
469   
470   cTOFeff->Print("TOF_eff.gif");
471   cTOFeff->Print("TOF_eff.eps");
472
473
474
475   TCanvas * cTOFtime = new TCanvas("cTOFtime", "TOF measured Times ", 400, 30, 550, 630) ;  
476   cTOFtime->Divide(1,2) ;
477   cTOFtime->cd(1);
478   cTOFtime->GetPad(1)->SetLogy(1);
479   fhTOFTime->GetXaxis()->SetTitle("TOF time (ns)");
480   fhTOFTime->GetYaxis()->SetTitle("Entries");
481   fhTOFTime->SetFillColor(4);
482   fhTOFTime->Draw();
483   cTOFtime->cd(2);
484   cTOFtime->GetPad(2)->SetLogy(1);
485   fhTOFDeltaTime->GetXaxis()->SetTitle("t^{TOF}-t^{exp}_{#pi} (ns)");
486   fhTOFDeltaTime->GetYaxis()->SetTitle("Entries");
487   fhTOFDeltaTime->SetFillColor(4);
488   fhTOFDeltaTime->Draw();
489
490   cTOFtime->Print("TOF_time.gif");
491   cTOFtime->Print("TOF_time.eps");
492
493
494   TCanvas * cTOFpid = new TCanvas("cTOFpid", "TOF PID ", 400, 30, 550, 630) ;  
495   cTOFpid->Divide(1,3) ;
496   cTOFpid->cd(1);
497   cTOFpid->SetLogy(1);
498   fhTOFMass->GetXaxis()->SetTitle("Reconstructed Mass (GeV/c^{2})");
499   fhTOFMass->GetYaxis()->SetTitle("Entries");
500   fhTOFMass->SetFillColor(4);
501   fhTOFMass->Draw();
502   cTOFpid->SetLogy(0);
503   cTOFpid->cd(2);
504   fhTOFMassVsMom->GetYaxis()->SetRange(0,400);
505   fhTOFMassVsMom->GetXaxis()->SetTitle("Reconstructed Mass (GeV/c^{2})");
506   fhTOFMassVsMom->GetYaxis()->SetTitle("Track Momentum (GeV/c)");
507   fhTOFMassVsMom->GetXaxis()->SetTitleSize(0.05);
508   fhTOFMassVsMom->GetYaxis()->SetTitleSize(0.05);
509   fhTOFMassVsMom->SetMarkerStyle(20);
510   fhTOFMassVsMom->SetMarkerSize(0.05);
511   fhTOFMassVsMom->SetMarkerColor(2);
512   fhTOFMassVsMom->Draw();
513   cTOFpid->cd(3);
514
515   TLatex *   tex = new TLatex(1., 1.25, "Bayesian PID: a-priori concentrations: [0,0,1,1,1]");
516   tex->SetTextColor(1);
517   tex->SetTextSize(0.045);
518   tex->SetLineWidth(2);
519
520   Float_t norm=1./fhTOFIDSpecies->GetEntries();
521   fhTOFIDSpecies->Scale(norm);
522   fhTOFIDSpecies->SetMaximum(1.2);
523   fhTOFIDSpecies->GetXaxis()->SetTitle("Particle Type");
524   fhTOFIDSpecies->GetYaxis()->SetTitle("ID Fractions");
525   fhTOFIDSpecies->GetXaxis()->SetTitleSize(0.05);
526   fhTOFIDSpecies->GetYaxis()->SetTitleSize(0.05);
527   fhTOFIDSpecies->SetFillColor(4);
528   fhTOFIDSpecies->Draw();
529   tex->Draw();
530  
531   char ch[10];
532
533   Float_t pifrac=fhTOFIDSpecies->GetBinContent(3);
534   Float_t kafrac=fhTOFIDSpecies->GetBinContent(4);
535   Float_t prfrac=fhTOFIDSpecies->GetBinContent(5);
536
537   sprintf(ch,"pion fraction   = %5.3f",pifrac);    
538   TLatex *   texpi = new TLatex(-0.3, 0.9, ch);
539   texpi->SetTextColor(1);
540   texpi->SetTextSize(0.05);
541   texpi->SetLineWidth(2);
542   texpi->Draw();
543   sprintf(ch,"kaon fraction   = %5.3f",kafrac);    
544   TLatex *   texka = new TLatex(-0.3, 0.8, ch);
545   texka->SetTextColor(1);
546   texka->SetTextSize(0.05);
547   texka->SetLineWidth(2);
548   texka->Draw();
549   sprintf(ch,"proton fraction = %5.3f",prfrac);    
550   TLatex *   texpr = new TLatex(-0.3, 0.7, ch);
551   texpr->SetTextColor(1);
552   texpr->SetTextSize(0.05);
553   texpr->SetLineWidth(2);
554   texpr->Draw();
555
556
557   cTOFpid->Print("TOF_pid.gif");
558   cTOFpid->Print("TOF_pid.eps");
559
560   // draw all 
561   //The General
562   
563
564   //The efficiency  
565 }
566
567 //______________________________________________________________________________
568 void AliTOFQATask::Terminate(Option_t *)
569 {
570   // Processing when the event loop is ended
571   
572   // some plots
573
574   AliInfo("TOF QA Task: End of events loop");
575   fOutputContainer = (TObjArray*)GetOutputData(0);
576   fhTOFMatch       = (TH1F*)fOutputContainer->At(0);
577   fhESDeffPhi      = (TH1F*)fOutputContainer->At(1);
578   fhESDeffTheta    = (TH1F*)fOutputContainer->At(2);
579   fhESDeffMom      = (TH1F*)fOutputContainer->At(3);
580   fhTOFeffPhi      = (TH1F*)fOutputContainer->At(4);
581   fhTOFeffPhiMT    = (TH1F*)fOutputContainer->At(5);
582   fhTOFeffTheta    = (TH1F*)fOutputContainer->At(6);
583   fhTOFeffThetaMT  = (TH1F*)fOutputContainer->At(7);
584   fhTOFeffMom      = (TH1F*)fOutputContainer->At(8);
585   fhTOFeffMomMT    = (TH1F*)fOutputContainer->At(9);
586   fhTOFsector      = (TH1F*)fOutputContainer->At(10);
587   fhTOFsectorMT    = (TH1F*)fOutputContainer->At(11);
588   fhTOFTime        = (TH1F*)fOutputContainer->At(12);
589   fhTOFDeltaTime   = (TH1F*)fOutputContainer->At(13);
590   fhTOFDeltaTimeMT = (TH1F*)fOutputContainer->At(14);
591   fhTOFIDSpecies   = (TH1F*)fOutputContainer->At(15);
592   fhTOFMassVsMom   = (TH1F*)fOutputContainer->At(16);
593   fhTOFMass        = (TH1F*)fOutputContainer->At(17);
594   GetEfficiency();
595 //  PostData(0, fOutputContainer);
596   DrawHistos() ; 
597 }