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