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