]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerSim.cxx
added verbosity to QA histograms (Yves)
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMakerSim.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 <TTree.h>
30 #include <TMath.h>
31
32 #include "AliLog.h"
33
34 #include "AliTOFdigit.h"
35 #include "AliTOFSDigit.h"
36 #include "AliTOFhitT0.h"
37 #include "AliTOFQADataMakerSim.h"
38 #include "AliQAChecker.h"
39 #include "AliTOFGeometry.h"
40
41
42 ClassImp(AliTOFQADataMakerSim)
43            
44 //____________________________________________________________________________ 
45   AliTOFQADataMakerSim::AliTOFQADataMakerSim() : 
46   AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker")
47 {
48   //
49   // ctor
50   //
51 }
52
53 //____________________________________________________________________________ 
54 AliTOFQADataMakerSim::AliTOFQADataMakerSim(const AliTOFQADataMakerSim& qadm) :
55   AliQADataMakerSim()
56 {
57   //
58   //copy ctor 
59   //
60   SetName((const char*)qadm.GetName()) ; 
61   SetTitle((const char*)qadm.GetTitle()); 
62 }
63
64 //__________________________________________________________________
65 AliTOFQADataMakerSim& AliTOFQADataMakerSim::operator = (const AliTOFQADataMakerSim& qadm )
66 {
67   //
68   // assignment operator.
69   //
70   this->~AliTOFQADataMakerSim();
71   new(this) AliTOFQADataMakerSim(qadm);
72   return *this;
73 }
74  
75 //____________________________________________________________________________ 
76 void AliTOFQADataMakerSim::InitHits()
77 {
78   //
79   // create Hits histograms in Hits subdir
80   //
81
82   const Bool_t expert   = kTRUE ; 
83   const Bool_t image    = kTRUE ; 
84   
85   TH1F * h0 = new TH1F("hTOFHits",    "Number of TOF Hits;TOF hit number [10 power];Counts ",301, -1.02, 5.) ; 
86   h0->Sumw2() ;
87   Add2HitsList(h0, 0, !expert, image) ;
88
89   TH1F * h1  = new TH1F("hTOFHitsTime", "Hits Time Spectrum in TOF (ns);Simulated TOF time [ns];Counts", 2000, 0., 200) ; 
90   h1->Sumw2() ;
91   Add2HitsList(h1, 1, !expert, image) ;
92
93   TH1F * h2  = new TH1F("hTOFHitsLength", "Length Spectrum in TOF (cm);Track length from primary vertex till hit TOF pad [cm];Counts", 500, 0., 500) ; 
94   h2->Sumw2() ;
95   Add2HitsList(h2, 2, !expert, image) ;
96
97   TH2F * h3  = new TH2F("hTOFHitsClusMap","Hits vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ; 
98   h3->Sumw2() ;
99   h3->GetYaxis()->SetTitleOffset(1.15);
100   Add2HitsList(h3, 3, !expert, image) ;
101 }
102
103 //____________________________________________________________________________ 
104 void AliTOFQADataMakerSim::InitDigits()
105 {
106   //
107   // create Digits histograms in Digits subdir
108   //
109
110   const Bool_t expert   = kTRUE ; 
111   const Bool_t image    = kTRUE ; 
112   
113   TH1F * h0 = new TH1F("hTOFDigits",    "Number of TOF Digit;TOF digit number [10 power];Counts ",301, -1.02, 5.) ;
114   h0->Sumw2() ;
115   Add2DigitsList(h0, 0, !expert, image) ;
116
117   TH1F * h1  = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 2000, 0., 200) ; 
118   h1->Sumw2() ;
119   Add2DigitsList(h1, 1, !expert, image) ;
120
121   TH1F * h2  = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 500, 0., 50) ; 
122   h2->Sumw2() ;
123   Add2DigitsList(h2, 2, !expert, image) ;
124
125   TH2F * h3  = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ; 
126   h3->Sumw2() ;
127   h3->GetYaxis()->SetTitleOffset(1.15);
128   Add2DigitsList(h3, 3, !expert, image) ;
129
130 }
131
132 //____________________________________________________________________________ 
133 void AliTOFQADataMakerSim::InitSDigits()
134 {
135   //
136   // create SDigits histograms in SDigits subdir
137   //
138
139   const Bool_t expert   = kTRUE ; 
140   const Bool_t image    = kTRUE ; 
141   
142   TH1F * h0 = new TH1F("hTOFSDigits",    "Number of TOF SDigits;TOF sdigit number [10 power];Counts ",301, -1.02, 5.) ;
143   h0->Sumw2() ;
144   Add2SDigitsList(h0, 0, !expert, image) ;
145
146   TH1F * h1  = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns);SDigitized TOF time [ns];Counts", 2000, 0., 200) ; 
147   h1->Sumw2() ;
148   Add2SDigitsList(h1, 1, !expert, image) ;
149
150   TH2F * h2  = new TH2F("hTOFSDigitsClusMap","SDigits vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ; 
151   h2->Sumw2() ;
152   h2->GetYaxis()->SetTitleOffset(1.15);
153   Add2SDigitsList(h2, 2, !expert, image) ;
154
155 }
156
157 //____________________________________________________________________________
158 void AliTOFQADataMakerSim::MakeHits(TClonesArray * hits)
159 {
160   //
161   //make QA data from Hits
162   //
163
164   Int_t in[5];
165   Int_t out[5];
166
167   Int_t nentries=hits->GetEntriesFast();
168   if(nentries<=0) {
169     GetHitsData(0)->Fill(-1.) ; 
170   } else{
171     GetHitsData(0)->Fill(TMath::Log10(nentries)) ; 
172   }
173   TIter next(hits) ; 
174   AliTOFhitT0 * hit ; 
175   while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
176
177     GetHitsData(1)->Fill( hit->GetTof()*1.E9) ;//in ns
178     GetHitsData(2)->Fill( hit->GetLen()) ;//in cm
179   
180     in[0] = hit->GetSector();
181     in[1] = hit->GetPlate();
182     in[2]= hit->GetStrip();
183     in[3]= hit->GetPadx();
184     in[4]= hit->GetPadz();
185     GetMapIndeces(in,out);
186     GetHitsData(3)->Fill( out[0],out[1]) ;//hit map
187   }
188
189 }
190
191
192 //____________________________________________________________________________
193 void AliTOFQADataMakerSim::MakeHits(TTree * hitTree)
194 {
195   //
196   // make QA data from Hit Tree
197   //
198   if(!hitTree){
199     AliError("can't get the tree with TOF hits !");
200     return;
201   }     
202
203   TBranch * branch = hitTree->GetBranch("TOF") ;
204
205   if (!branch ) {
206     AliError("TOF branch in Hit Tree not found") ; 
207     return;
208   }
209
210   static TClonesArray statichits("AliTOFhitT0", 1000);
211   statichits.Clear();
212   TClonesArray *hits = &statichits;
213   static TClonesArray staticdummy("AliTOFhitT0", 1000);
214   staticdummy.Clear();
215   TClonesArray *dummy = &staticdummy;
216   branch->SetAddress(&dummy);
217   Int_t index = 0 ;  
218   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
219     branch->GetEntry(ientry) ; 
220     for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
221       AliTOFhitT0 * hit = dynamic_cast<AliTOFhitT0 *> (dummy->At(ihit)) ; 
222       new((*hits)[index]) AliTOFhitT0(*hit) ; 
223       index++ ;
224     } 
225   }     
226
227   MakeHits(hits) ; 
228
229 }
230
231 //____________________________________________________________________________
232 void AliTOFQADataMakerSim::MakeDigits(TClonesArray * digits)
233 {
234   //
235   // makes data from Digits
236   //
237   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
238   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
239   Int_t in[5];
240   Int_t out[5];
241
242   Int_t nentries=digits->GetEntriesFast();
243   if(nentries<=0){
244     GetDigitsData(0)->Fill(-1.) ; 
245   }else{
246     GetDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
247   } 
248
249   TIter next(digits) ; 
250   AliTOFdigit * digit ; 
251   while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
252     
253     GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
254     GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
255
256     in[0] = digit->GetSector();
257     in[1] = digit->GetPlate();
258     in[2] = digit->GetStrip();
259     in[3] = digit->GetPadx();
260     in[4]= digit->GetPadz();
261     GetMapIndeces(in,out);
262     GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
263   }
264
265 }
266
267
268 //____________________________________________________________________________
269 void AliTOFQADataMakerSim::MakeDigits(TTree * digitTree)
270 {
271   //
272   // makes data from Digit Tree
273   //
274   TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ; 
275   
276   TBranch * branch = digitTree->GetBranch("TOF") ;
277   if ( ! branch ) {
278     AliError("TOF branch in Digit Tree not found") ; 
279     return;
280   }
281   branch->SetAddress(&digits) ;
282   branch->GetEntry(0) ; 
283   MakeDigits(digits) ; 
284 }
285
286 //____________________________________________________________________________
287 void AliTOFQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
288 {
289   //
290   // makes data from SDigits
291   //
292
293   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
294   Int_t in[5];
295   Int_t out[5];
296
297   Int_t nentries=sdigits->GetEntriesFast();
298   if(nentries<=0){
299     GetSDigitsData(0)->Fill(-1.) ; 
300   }else{
301     GetSDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
302   } 
303
304   TIter next(sdigits) ; 
305   AliTOFSDigit * sdigit ; 
306   while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
307     
308     for(Int_t i=0;i<sdigit->GetNDigits();i++){
309       GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
310     }
311
312     in[0] = sdigit->GetSector();
313     in[1] = sdigit->GetPlate();
314     in[2] = sdigit->GetStrip();
315     in[3] = sdigit->GetPadx();
316     in[4]= sdigit->GetPadz();
317     GetMapIndeces(in,out);
318     GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
319   }
320 }
321
322 //____________________________________________________________________________
323 void AliTOFQADataMakerSim::MakeSDigits(TTree * sdigitTree)
324 {
325   //
326   // makes data from SDigit Tree
327   //
328   TClonesArray * sdigits = new TClonesArray("AliTOFSDigit", 1000) ; 
329   
330   TBranch * branch = sdigitTree->GetBranch("TOF") ;
331   if ( ! branch ) {
332     AliError("TOF branch in SDigit Tree not found") ; 
333     return;
334   }
335   branch->SetAddress(&sdigits) ;
336   branch->GetEntry(0) ; 
337   MakeSDigits(sdigits) ; 
338 }
339
340 //____________________________________________________________________________ 
341 void AliTOFQADataMakerSim::StartOfDetectorCycle()
342 {
343   //
344   //Detector specific actions at start of cycle
345   //to be implemented  
346 }
347
348 //____________________________________________________________________________ 
349 void AliTOFQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
350 {
351   //Detector specific actions at end of cycle
352   // do the QA checking
353
354   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
355 }
356 //____________________________________________________________________________
357 void AliTOFQADataMakerSim::GetMapIndeces(Int_t* in , Int_t* out)
358 {
359   //
360   //return appropriate indeces for the theta-phi map
361   //
362
363   Int_t npadX = AliTOFGeometry::NpadX();
364   Int_t npadZ = AliTOFGeometry::NpadZ();
365   Int_t nStripA = AliTOFGeometry::NStripA();
366   Int_t nStripB = AliTOFGeometry::NStripB();
367   Int_t nStripC = AliTOFGeometry::NStripC();
368
369   Int_t isector = in[0];
370   Int_t iplate = in[1];
371   Int_t istrip = in[2];
372   Int_t ipadX = in[3]; 
373   Int_t ipadZ = in[4]; 
374   
375   Int_t stripOffset = 0;
376   switch (iplate) {
377   case 0:
378     stripOffset = 0;
379       break;
380   case 1:
381     stripOffset = nStripC;
382     break;
383   case 2:
384     stripOffset = nStripC+nStripB;
385     break;
386   case 3:
387     stripOffset = nStripC+nStripB+nStripA;
388     break;
389   case 4:
390     stripOffset = nStripC+nStripB+nStripA+nStripB;
391     break;
392   default:
393     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
394     break;
395   };
396   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
397   Int_t phiindex=npadX*isector+ipadX+1;
398   out[0]=zindex;  
399   out[1]=phiindex;  
400   
401 }