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