]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerSim.cxx
Coding rules
[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()
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= fHitsArray->GetEntriesFast();
168   if(nentries<=0) {
169     GetHitsData(0)->Fill(-1.) ; 
170   } else{
171     GetHitsData(0)->Fill(TMath::Log10(nentries)) ; 
172   }
173   TIter next(fHitsArray) ; 
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 void AliTOFQADataMakerSim::MakeHits(TTree * hitTree)
193 {
194   //
195   // make QA data from Hit Tree
196   //
197   if(!hitTree){
198     AliError("can't get the tree with TOF hits !");
199     return; 
200   }     
201
202   TBranch * branch = hitTree->GetBranch("TOF") ;
203
204   if (!branch ) {
205     AliError("TOF branch in Hit Tree not found") ; 
206     return;
207   }
208
209   if (fHitsArray) 
210     fHitsArray->Clear() ; 
211   else 
212     fHitsArray = new TClonesArray("AliTOFhitT0", 1000) ;
213   
214   branch->SetAddress(&fHitsArray);
215   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
216     branch->GetEntry(ientry) ; 
217     MakeHits() ; 
218     fHitsArray->Clear() ; 
219   }     
220 }
221
222 //____________________________________________________________________________
223 void AliTOFQADataMakerSim::MakeDigits()
224 {
225   //
226   // makes data from Digits
227   //
228   
229   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
230   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
231   Int_t in[5];
232   Int_t out[5];
233
234   Int_t nentries=fDigitsArray->GetEntriesFast();
235   if(nentries<=0){
236     GetDigitsData(0)->Fill(-1.) ; 
237   }else{
238     GetDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
239   } 
240
241   TIter next(fDigitsArray) ; 
242   AliTOFdigit * digit ; 
243   while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
244     
245     GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
246     GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
247
248     in[0] = digit->GetSector();
249     in[1] = digit->GetPlate();
250     in[2] = digit->GetStrip();
251     in[3] = digit->GetPadx();
252     in[4]= digit->GetPadz();
253     GetMapIndeces(in,out);
254     GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
255   }
256
257 }
258
259
260 //____________________________________________________________________________
261 void AliTOFQADataMakerSim::MakeDigits(TTree * digitTree)
262 {
263   //
264   // makes data from Digit Tree
265   //
266   if (fDigitsArray) 
267     fDigitsArray->Clear() ; 
268   else 
269     fDigitsArray = new TClonesArray("AliTOFdigit", 1000) ; 
270   
271   TBranch * branch = digitTree->GetBranch("TOF") ;
272   if ( ! branch ) {
273     AliError("TOF branch in Digit Tree not found") ; 
274     return;
275   }
276   branch->SetAddress(&fDigitsArray) ;
277   branch->GetEntry(0) ; 
278   MakeDigits() ; 
279 }
280
281 //____________________________________________________________________________
282 void AliTOFQADataMakerSim::MakeSDigits()
283 {
284   //
285   // makes data from SDigits
286   //
287   
288   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
289   Int_t in[5];
290   Int_t out[5];
291
292   Int_t nentries=fSDigitsArray->GetEntriesFast();
293   if(nentries<=0){
294     GetSDigitsData(0)->Fill(-1.) ; 
295   }else{
296     GetSDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
297   } 
298
299   TIter next(fSDigitsArray) ; 
300   AliTOFSDigit * sdigit ; 
301   while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
302     
303     for(Int_t i=0;i<sdigit->GetNDigits();i++){
304       GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
305     }
306
307     in[0] = sdigit->GetSector();
308     in[1] = sdigit->GetPlate();
309     in[2] = sdigit->GetStrip();
310     in[3] = sdigit->GetPadx();
311     in[4]= sdigit->GetPadz();
312     GetMapIndeces(in,out);
313     GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
314   }
315 }
316
317 //____________________________________________________________________________
318 void AliTOFQADataMakerSim::MakeSDigits(TTree * sdigitTree)
319 {
320   //
321   // makes data from SDigit Tree
322   //
323   if (fSDigitsArray) 
324     fSDigitsArray->Clear() ; 
325   else 
326     fSDigitsArray = new TClonesArray("AliTOFSDigit", 1000) ; 
327   
328   TBranch * branch = sdigitTree->GetBranch("TOF") ;
329   if ( ! branch ) {
330     AliError("TOF branch in SDigit Tree not found") ; 
331     return;
332   }
333   branch->SetAddress(&fSDigitsArray) ;
334   branch->GetEntry(0) ; 
335   MakeSDigits() ; 
336 }
337
338 //____________________________________________________________________________ 
339 void AliTOFQADataMakerSim::StartOfDetectorCycle()
340 {
341   //
342   //Detector specific actions at start of cycle
343   //to be implemented  
344 }
345
346 //____________________________________________________________________________ 
347 void AliTOFQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
348 {
349   //Detector specific actions at end of cycle
350   // do the QA checking
351
352   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
353 }
354 //____________________________________________________________________________
355 void AliTOFQADataMakerSim::GetMapIndeces(Int_t* in , Int_t* out)
356 {
357   //
358   //return appropriate indeces for the theta-phi map
359   //
360
361   Int_t npadX = AliTOFGeometry::NpadX();
362   Int_t npadZ = AliTOFGeometry::NpadZ();
363   Int_t nStripA = AliTOFGeometry::NStripA();
364   Int_t nStripB = AliTOFGeometry::NStripB();
365   Int_t nStripC = AliTOFGeometry::NStripC();
366
367   Int_t isector = in[0];
368   Int_t iplate = in[1];
369   Int_t istrip = in[2];
370   Int_t ipadX = in[3]; 
371   Int_t ipadZ = in[4]; 
372   
373   Int_t stripOffset = 0;
374   switch (iplate) {
375   case 0:
376     stripOffset = 0;
377       break;
378   case 1:
379     stripOffset = nStripC;
380     break;
381   case 2:
382     stripOffset = nStripC+nStripB;
383     break;
384   case 3:
385     stripOffset = nStripC+nStripB+nStripA;
386     break;
387   case 4:
388     stripOffset = nStripC+nStripB+nStripA+nStripB;
389     break;
390   default:
391     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
392     break;
393   };
394   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
395   Int_t phiindex=npadX*isector+ipadX+1;
396   out[0]=zindex;  
397   out[1]=phiindex;  
398   
399 }