]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMaker.cxx
Improved memory leaks in AliTOFtracker* classes: additional solution for bug #66136
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMaker.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 ///////////////////////////////////////////////////////////////////////
17 //                                                                   //
18 //  Produces the data needed to calculate the TOF quality assurance. //
19 //  QA objects are 1 & 2 Dimensional histograms.                     //
20 //  author: S.Arcelli                                                //
21 //                                                                   //
22 ///////////////////////////////////////////////////////////////////////
23
24 #include <TClonesArray.h>
25 #include <TFile.h> 
26 #include <TH1I.h> 
27 #include <TH1F.h> 
28 #include <TH2F.h> 
29 #include "AliESDEvent.h"
30 #include "AliESDtrack.h"
31 #include "AliTOFcluster.h"
32 #include "AliTOFdigit.h"
33 #include "AliTOFSDigit.h"
34 #include "AliTOFhitT0.h"
35 #include "AliTOFQADataMaker.h"
36 #include "AliQAChecker.h"
37 #include "AliRawReader.h"
38 #include "AliTOFRawStream.h"
39 #include "AliTOFrawData.h"
40 #include "AliLog.h"
41
42 ClassImp(AliTOFQADataMaker)
43            
44 //____________________________________________________________________________ 
45   AliTOFQADataMaker::AliTOFQADataMaker() : 
46   AliQADataMaker(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker")
47 {
48   //
49   // ctor
50   //
51 }
52
53 //____________________________________________________________________________ 
54 AliTOFQADataMaker::AliTOFQADataMaker(const AliTOFQADataMaker& qadm) :
55   AliQADataMaker()
56 {
57   //
58   //copy ctor 
59   //
60   SetName((const char*)qadm.GetName()) ; 
61   SetTitle((const char*)qadm.GetTitle()); 
62 }
63
64 //__________________________________________________________________
65 AliTOFQADataMaker& AliTOFQADataMaker::operator = (const AliTOFQADataMaker& qadm )
66 {
67   //
68   // assignment operator.
69   //
70   this->~AliTOFQADataMaker();
71   new(this) AliTOFQADataMaker(qadm);
72   return *this;
73 }
74  
75 //____________________________________________________________________________ 
76 void AliTOFQADataMaker::InitHits()
77 {
78   //
79   // create Hits histograms in Hits subdir
80   //
81
82   Bool_t expert = kFALSE;
83
84   TH1F * h0 = new TH1F("hTOFHits",    "Number of TOF Hits ",301, -1.02, 5.) ;
85   h0->Sumw2() ;
86   h0->GetXaxis()->SetTitle("TOF hit number [10 power]");
87   Add2HitsList(h0, 0, expert) ;
88
89   TH1F * h1  = new TH1F("hTOFHitsTime", "Hits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
90   h1->Sumw2() ;
91   h1->GetXaxis()->SetTitle("Simulated TOF time [ns]");
92   Add2HitsList(h1, 1, expert) ;
93
94   TH1F * h2  = new TH1F("hTOFHitsLength", "Length Spectrum in TOF (cm)", 500, 0., 500) ;
95   h2->Sumw2() ;
96   h2->GetXaxis()->SetTitle("Track length from primary vertex till hit TOF pad [cm]");
97   Add2HitsList(h2, 2, expert) ;
98
99   TH2F * h3  = new TH2F("hTOFHitsClusMap","Hits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
100   h3->Sumw2() ;
101   h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
102   h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
103   h3->GetYaxis()->SetTitleOffset(1.15);
104   Add2HitsList(h3, 3, expert) ;
105 }
106
107 //____________________________________________________________________________ 
108 void AliTOFQADataMaker::InitDigits()
109 {
110   //
111   // create Digits histograms in Digits subdir
112   //
113
114   Bool_t expert = kFALSE;
115
116   TH1F * h0 = new TH1F("hTOFDigits",    "Number of TOF Digits ",301, -1.02, 5.) ;
117   h0->Sumw2() ;
118   h0->GetXaxis()->SetTitle("TOF digit number [10 power]");
119   Add2DigitsList(h0, 0, expert) ;
120
121   TH1F * h1  = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
122   h1->Sumw2() ;
123   h1->GetXaxis()->SetTitle("Digitized TOF time [ns]");
124   Add2DigitsList(h1, 1, expert) ;
125
126   TH1F * h2  = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns)", 500, 0., 50) ;
127   h2->Sumw2() ;
128   h2->GetXaxis()->SetTitle("Digitized TOT [ns]");
129   h2->GetYaxis()->SetTitle("Digitized TOF time [ns]");
130   Add2DigitsList(h2, 2, expert) ;
131
132   TH2F * h3  = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
133   h3->Sumw2() ;
134   h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
135   h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
136   h3->GetYaxis()->SetTitleOffset(1.15);
137   Add2DigitsList(h3, 3, expert) ;
138
139 }
140
141 //____________________________________________________________________________ 
142 void AliTOFQADataMaker::InitSDigits()
143 {
144   //
145   // create SDigits histograms in SDigits subdir
146   //
147
148   Bool_t expert = kFALSE;
149
150   TH1F * h0 = new TH1F("hTOFSDigits",    "Number of TOF SDigits ",301, -1.02, 5.) ;
151   h0->Sumw2() ;
152   h0->GetXaxis()->SetTitle("TOF sdigit number [10 power]");
153   Add2SDigitsList(h0, 0, expert) ;
154
155   TH1F * h1  = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
156   h1->Sumw2() ;
157   h1->GetXaxis()->SetTitle("SDigitized TOF time [ns]");
158   Add2SDigitsList(h1, 1, expert) ;
159
160   TH2F * h2  = new TH2F("hTOFSDigitsClusMap","SDigits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
161   h2->Sumw2() ;
162   h2->GetXaxis()->SetTitle("2*strip+padz (eta)");
163   h2->GetYaxis()->SetTitle("48*sector+padx (phi)");
164   h2->GetYaxis()->SetTitleOffset(1.15);
165   Add2SDigitsList(h2, 2, expert) ;
166
167 }
168
169 //____________________________________________________________________________ 
170 void AliTOFQADataMaker::InitRaws()
171 {
172   //
173   // create Raws histograms in Raws subdir
174   //
175
176   Bool_t expert = kFALSE;
177
178   TH1F * h0 = new TH1F("hTOFRaws",    "Number of TOF Raws ",301, -1.02, 5.) ;
179   h0->Sumw2() ;
180   h0->GetXaxis()->SetTitle("TOF raw number [10 power]");
181   Add2RawsList(h0, 0, expert) ;
182
183   TH1F * h1  = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
184   h1->Sumw2() ;
185   h1->GetXaxis()->SetTitle("Measured TOF time [ns]");
186   Add2RawsList(h1, 1, expert) ;
187
188   TH1F * h2  = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns)", 500, 0., 50) ;
189   h2->Sumw2() ;
190   h2->GetXaxis()->SetTitle("Measured TOT [ns]");
191   h2->GetYaxis()->SetTitle("Measured TOF time [ns]");
192   Add2RawsList(h2, 2, expert) ;
193
194   TH2F * h3  = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
195   h3->Sumw2() ;
196   h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
197   h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
198   h3->GetYaxis()->SetTitleOffset(1.15);
199   Add2RawsList(h3, 3, expert) ;
200
201 }
202
203 //____________________________________________________________________________ 
204 void AliTOFQADataMaker::InitRecPoints()
205 {
206   //
207   // create RecPoints histograms in RecPoints subdir
208   //
209
210   Bool_t expert = kFALSE;
211
212   TH1F * h0 = new TH1F("hTOFRecPoints",    "Number of TOF RecPoints ",301, -1.02, 5.) ;
213   h0->Sumw2() ;
214   h0->GetXaxis()->SetTitle("TOF recPoint number [10 power]");
215   Add2RecPointsList(h0, 0, expert) ;
216
217   TH1F * h1  = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
218   h1->Sumw2() ;
219   h1->GetXaxis()->SetTitle("Calibrated TOF time [ns]");
220   Add2RecPointsList(h1, 1, expert) ;
221
222   TH1F * h2  = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
223   h2->Sumw2() ;
224   h2->GetXaxis()->SetTitle("Measured TOT [ns]");
225   h2->GetYaxis()->SetTitle("Measured TOF time [ns]");
226   Add2RecPointsList(h2, 2, expert) ;
227
228   TH1F * h3  = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns)", 500, 0., 50) ;
229   h3->Sumw2() ;
230   h3->GetXaxis()->SetTitle("Measured TOT [ns]");
231   Add2RecPointsList(h3, 3, expert) ;
232
233   TH2F * h4  = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta",183, -0.5, 182.5,865,-0.5,864.5) ;
234   h4->Sumw2() ;
235   h4->GetXaxis()->SetTitle("2*strip+padz (eta)");
236   h4->GetYaxis()->SetTitle("48*sector+padx (phi)");
237   h4->GetYaxis()->SetTitleOffset(1.15);
238   Add2RecPointsList(h4, 4, expert) ;
239
240 }
241
242 //____________________________________________________________________________ 
243 void AliTOFQADataMaker::InitESDs()
244 {
245   //
246   //create ESDs histograms in ESDs subdir
247   //
248
249   Bool_t expert = kFALSE;
250
251   TH1F * h0 = new TH1F("hTOFESDs",    "Number of matched TOF tracks over ESDs",       250, -1., 4.) ;
252   h0->Sumw2() ; 
253   h0->GetXaxis()->SetTitle("Number of TOF matched ESD tracks [10 power]");
254   Add2ESDsList(h0, 0, expert) ;
255
256   TH1F * h1  = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
257   h1->Sumw2() ;
258   h1->GetXaxis()->SetTitle("Calibrated TOF time [ns]");
259   Add2ESDsList(h1, 1, expert) ;
260
261   TH1F * h2  = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
262   h2->Sumw2() ;
263   h2->GetXaxis()->SetTitle("Measured TOF time [ns]");
264   Add2ESDsList(h2, 2, expert) ;
265
266   TH1F * h3  = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 500, 0., 50) ;
267   h3->Sumw2() ;
268   h3->GetXaxis()->SetTitle("Measured TOT [ns]");
269   Add2ESDsList(h3, 3, expert) ;
270
271   TH1F * h4 = new TH1F("hTOFESDsPID",    "Fraction of matched TOF tracks with good PID flag", 101, 0., 101.) ;
272   h4->Sumw2() ;
273   h4->GetXaxis()->SetTitle("Fraction of TOF matched ESD tracks with good flag [%]");
274   Add2ESDsList(h4, 4, expert) ;
275 }
276
277
278 //____________________________________________________________________________
279 void AliTOFQADataMaker::MakeHits(TClonesArray * hits)
280 {
281   //
282   //make QA data from Hits
283   //
284
285   Int_t in[5];
286   Int_t out[5];
287
288   Int_t nentries=hits->GetEntriesFast();
289   if(nentries<=0) {
290     GetHitsData(0)->Fill(-1.) ; 
291   } else{
292     GetHitsData(0)->Fill(TMath::Log10(nentries)) ; 
293   }
294   TIter next(hits) ; 
295   AliTOFhitT0 * hit ; 
296   while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
297
298     GetHitsData(1)->Fill( hit->GetTof()*1.E9) ;//in ns
299     GetHitsData(2)->Fill( hit->GetLen()) ;//in cm
300   
301     in[0] = hit->GetSector();
302     in[1] = hit->GetPlate();
303     in[2]= hit->GetStrip();
304     in[3]= hit->GetPadx();
305     in[4]= hit->GetPadz();
306     GetMapIndeces(in,out);
307     GetHitsData(3)->Fill( out[0],out[1]) ;//hit map
308   }
309
310 }
311
312
313 //____________________________________________________________________________
314 void AliTOFQADataMaker::MakeHits(TTree * hitTree)
315 {
316   //
317   // make QA data from Hit Tree
318   //
319   if(!hitTree){
320     AliError("can't get the tree with TOF hits !");
321     return;
322   }     
323
324   TBranch * branch = hitTree->GetBranch("TOF") ;
325
326   if (!branch ) {
327     AliError("TOF branch in Hit Tree not found") ; 
328     return;
329   }
330
331   TClonesArray * hits = new TClonesArray("AliTOFhitT0", 1000);
332   TClonesArray * dummy = new TClonesArray("AliTOFhitT0", 1000);
333   branch->SetAddress(&dummy);
334   Int_t index = 0 ;  
335   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
336     branch->GetEntry(ientry) ; 
337     for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
338       AliTOFhitT0 * hit = dynamic_cast<AliTOFhitT0 *> (dummy->At(ihit)) ; 
339       new((*hits)[index]) AliTOFhitT0(*hit) ; 
340       index++ ;
341     } 
342   }     
343
344   dummy->Delete();
345   delete dummy;
346   MakeHits(hits) ; 
347
348 }
349
350 //____________________________________________________________________________
351 void AliTOFQADataMaker::MakeDigits(TClonesArray * digits)
352 {
353   //
354   // makes data from Digits
355   //
356   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
357   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
358   Int_t in[5];
359   Int_t out[5];
360
361   Int_t nentries=digits->GetEntriesFast();
362   if(nentries<=0){
363     GetDigitsData(0)->Fill(-1.) ; 
364   }else{
365     GetDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
366   } 
367
368   TIter next(digits) ; 
369   AliTOFdigit * digit ; 
370   while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
371     
372     GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
373     GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
374
375     in[0] = digit->GetSector();
376     in[1] = digit->GetPlate();
377     in[2] = digit->GetStrip();
378     in[3] = digit->GetPadx();
379     in[4]= digit->GetPadz();
380     GetMapIndeces(in,out);
381     GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
382   }
383
384 }
385
386
387 //____________________________________________________________________________
388 void AliTOFQADataMaker::MakeDigits(TTree * digitTree)
389 {
390   //
391   // makes data from Digit Tree
392   //
393   TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ; 
394   
395   TBranch * branch = digitTree->GetBranch("TOF") ;
396   if ( ! branch ) {
397     AliError("TOF branch in Digit Tree not found") ; 
398     return;
399   }
400   branch->SetAddress(&digits) ;
401   branch->GetEntry(0) ; 
402   MakeDigits(digits) ; 
403 }
404
405 //____________________________________________________________________________
406 void AliTOFQADataMaker::MakeSDigits(TClonesArray * sdigits)
407 {
408   //
409   // makes data from SDigits
410   //
411
412   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
413   Int_t in[5];
414   Int_t out[5];
415
416   Int_t nentries=sdigits->GetEntriesFast();
417   if(nentries<=0){
418     GetSDigitsData(0)->Fill(-1.) ; 
419   }else{
420     GetSDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
421   } 
422
423   TIter next(sdigits) ; 
424   AliTOFSDigit * sdigit ; 
425   while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
426     
427     for(Int_t i=0;i<sdigit->GetNDigits();i++){
428       GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
429     }
430
431     in[0] = sdigit->GetSector();
432     in[1] = sdigit->GetPlate();
433     in[2] = sdigit->GetStrip();
434     in[3] = sdigit->GetPadx();
435     in[4]= sdigit->GetPadz();
436     GetMapIndeces(in,out);
437     GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
438   }
439 }
440
441 //____________________________________________________________________________
442 void AliTOFQADataMaker::MakeSDigits(TTree * sdigitTree)
443 {
444   //
445   // makes data from SDigit Tree
446   //
447   TClonesArray * sdigits = new TClonesArray("AliTOFSDigit", 1000) ; 
448   
449   TBranch * branch = sdigitTree->GetBranch("TOF") ;
450   if ( ! branch ) {
451     AliError("TOF branch in SDigit Tree not found") ; 
452     return;
453   }
454   branch->SetAddress(&sdigits) ;
455   branch->GetEntry(0) ; 
456   MakeSDigits(sdigits) ; 
457 }
458
459 //____________________________________________________________________________
460 void AliTOFQADataMaker::MakeRaws(AliRawReader* rawReader)
461 {
462   //
463   // makes data from Raws
464   //
465
466   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
467   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
468
469
470   Int_t ntof = 0 ; 
471   Int_t in[5];
472   Int_t out[5];
473
474   TClonesArray * clonesRawData;
475   AliTOFRawStream tofInput(rawReader);
476   for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
477     rawReader->Reset();
478     tofInput.LoadRawData(iDDL);
479     clonesRawData = (TClonesArray*)tofInput.GetRawData();
480     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
481       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
482       //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
483       if (!tofRawDatum->GetTOF()) continue;
484       ntof++;
485       GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
486       GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
487
488       tofInput.EquipmentId2VolumeId(iDDL, 
489                                     tofRawDatum->GetTRM(), 
490                                     tofRawDatum->GetTRMchain(),
491                                     tofRawDatum->GetTDC(), 
492                                     tofRawDatum->GetTDCchannel(), 
493                                     in);
494     
495       GetMapIndeces(in,out);
496       GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
497       
498     } // while loop
499     
500     clonesRawData->Clear();
501     
502   } // DDL Loop
503   
504   Int_t nentries=ntof;
505   if(nentries<=0){
506     GetRawsData(0)->Fill(-1.) ; 
507   }else{
508     GetRawsData(0)->Fill(TMath::Log10(nentries)) ; 
509   }
510 }
511
512 //____________________________________________________________________________
513 void AliTOFQADataMaker::MakeRecPoints(TTree * clustersTree)
514 {
515   //
516   // Make data from Clusters
517   //
518
519   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
520   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
521
522   Int_t in[5];
523   Int_t out[5];
524
525   TBranch *branch=clustersTree->GetBranch("TOF");
526   if (!branch) { 
527     AliError("can't get the branch with the TOF clusters !");
528     return;
529   }
530
531   TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
532   branch->SetAddress(&clusters);
533
534   // Import the tree
535   clustersTree->GetEvent(0);  
536   
537   Int_t nentries=clusters->GetEntriesFast();
538   if(nentries<=0){
539     GetRecPointsData(0)->Fill(-1.) ; 
540   }else{
541     GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ; 
542   } 
543  
544   TIter next(clusters) ; 
545   AliTOFcluster * c ; 
546   while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
547     GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
548     GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
549     GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
550     
551     in[0] = c->GetDetInd(0);
552     in[1] = c->GetDetInd(1);
553     in[2] = c->GetDetInd(2);
554     in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
555     in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
556     
557     GetMapIndeces(in,out);
558     GetRecPointsData(4)->Fill(out[0],out[1]);
559     
560   }
561 }
562
563 //____________________________________________________________________________
564 void AliTOFQADataMaker::MakeESDs(AliESDEvent * esd)
565 {
566   //
567   // make QA data from ESDs
568   //  
569   Int_t ntrk = esd->GetNumberOfTracks() ; 
570   Int_t ntof=0;
571   Int_t ntofpid=0;
572   while (ntrk--) {
573     AliESDtrack *track=esd->GetTrack(ntrk);
574     Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
575     Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
576     Double_t tofToT=track->GetTOFsignalToT(); //in ns
577     if(!(tofTime>0))continue;
578     ntof++;
579     GetESDsData(1)->Fill(tofTime);
580     GetESDsData(2)->Fill(tofTimeRaw); 
581     GetESDsData(3)->Fill(tofToT);
582     //check how many tracks where ESD PID is ok 
583     UInt_t status=track->GetStatus();
584     if (((status&AliESDtrack::kESDpid)==0) || 
585         ((status&AliESDtrack::kTOFpid)==0)) continue;
586     ntofpid++;
587   }
588   
589   Int_t nentries=ntof;
590   if(nentries<=0){
591     GetESDsData(0)->Fill(-1.) ;
592   }else{
593     GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
594   }
595
596   Float_t ratio = 0.;
597   if(ntof>0) ratio=(Float_t)ntofpid/(Float_t)ntof*100.;
598   GetESDsData(4)->Fill(ratio) ;
599
600 }
601
602 //____________________________________________________________________________ 
603 void AliTOFQADataMaker::StartOfDetectorCycle()
604 {
605   //
606   //Detector specific actions at start of cycle
607   //to be implemented  
608 }
609
610 //____________________________________________________________________________ 
611 void AliTOFQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX task, TObjArray * list)
612 {
613   //Detector specific actions at end of cycle
614   // do the QA checking
615
616   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
617 }
618 //____________________________________________________________________________
619 void AliTOFQADataMaker::GetMapIndeces(Int_t* in , Int_t* out)
620 {
621   //
622   //return appropriate indeces for the theta-phi map
623   //
624
625   Int_t npadX = AliTOFGeometry::NpadX();
626   Int_t npadZ = AliTOFGeometry::NpadZ();
627   Int_t nStripA = AliTOFGeometry::NStripA();
628   Int_t nStripB = AliTOFGeometry::NStripB();
629   Int_t nStripC = AliTOFGeometry::NStripC();
630
631   Int_t isector = in[0];
632   Int_t iplate = in[1];
633   Int_t istrip = in[2];
634   Int_t ipadX = in[3]; 
635   Int_t ipadZ = in[4]; 
636   
637   Int_t stripOffset = 0;
638   switch (iplate) {
639   case 0:
640     stripOffset = 0;
641       break;
642   case 1:
643     stripOffset = nStripC;
644     break;
645   case 2:
646     stripOffset = nStripC+nStripB;
647     break;
648   case 3:
649     stripOffset = nStripC+nStripB+nStripA;
650     break;
651   case 4:
652     stripOffset = nStripC+nStripB+nStripA+nStripB;
653     break;
654   default:
655     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
656     break;
657   };
658   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
659   Int_t phiindex=npadX*isector+ipadX+1;
660   out[0]=zindex;  
661   out[1]=phiindex;  
662   
663 }