]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerSim.cxx
Fixes for bug #49914: Compilation breaks in trunk, and bug #48629: Trunk cannot read...
[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   Bool_t expert = kFALSE;
83
84   TH1F * h0 = new TH1F("hTOFHits",    "Number of TOF Hits ",301, -1.02, 5.) ; 
85   h0->Sumw2() ;
86   Add2HitsList(h0, 0, expert) ;
87
88   TH1F * h1  = new TH1F("hTOFHitsTime", "Hits Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
89   h1->Sumw2() ;
90   Add2HitsList(h1, 1, expert) ;
91
92   TH1F * h2  = new TH1F("hTOFHitsLength", "Length Spectrum in TOF (cm)", 500, 0., 500) ; 
93   h2->Sumw2() ;
94   Add2HitsList(h2, 2, expert) ;
95
96   TH2F * h3  = new TH2F("hTOFHitsClusMap","Hits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ; 
97   h3->Sumw2() ;
98   Add2HitsList(h3, 3, expert) ;
99 }
100
101 //____________________________________________________________________________ 
102 void AliTOFQADataMakerSim::InitDigits()
103 {
104   //
105   // create Digits histograms in Digits subdir
106   //
107
108   Bool_t expert = kFALSE;
109
110   TH1F * h0 = new TH1F("hTOFDigits",    "Number of TOF Digits ",301, -1.02, 5.) ;   h0->Sumw2() ;
111   Add2DigitsList(h0, 0, expert) ;
112
113   TH1F * h1  = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
114   h1->Sumw2() ;
115   Add2DigitsList(h1, 1, expert) ;
116
117   TH1F * h2  = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns)", 500, 0., 50) ; 
118   h2->Sumw2() ;
119   Add2DigitsList(h2, 2, expert) ;
120
121   TH2F * h3  = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ; 
122   h3->Sumw2() ;
123   Add2DigitsList(h3, 3, expert) ;
124
125 }
126
127 //____________________________________________________________________________ 
128 void AliTOFQADataMakerSim::InitSDigits()
129 {
130   //
131   // create SDigits histograms in SDigits subdir
132   //
133
134   Bool_t expert = kFALSE;
135
136   TH1F * h0 = new TH1F("hTOFSDigits",    "Number of TOF SDigits ",301, -1.02, 5.) ;   h0->Sumw2() ;
137   Add2SDigitsList(h0, 0, expert) ;
138
139   TH1F * h1  = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
140   h1->Sumw2() ;
141   Add2SDigitsList(h1, 1, expert) ;
142
143   TH2F * h2  = new TH2F("hTOFSDigitsClusMap","SDigits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ; 
144   h2->Sumw2() ;
145   Add2SDigitsList(h2, 2, expert) ;
146
147 }
148
149 //____________________________________________________________________________
150 void AliTOFQADataMakerSim::MakeHits(TClonesArray * hits)
151 {
152   //
153   //make QA data from Hits
154   //
155
156   Int_t in[5];
157   Int_t out[5];
158
159   Int_t nentries=hits->GetEntriesFast();
160   if(nentries<=0) {
161     GetHitsData(0)->Fill(-1.) ; 
162   } else{
163     GetHitsData(0)->Fill(TMath::Log10(nentries)) ; 
164   }
165   TIter next(hits) ; 
166   AliTOFhitT0 * hit ; 
167   while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
168
169     GetHitsData(1)->Fill( hit->GetTof()*1.E9) ;//in ns
170     GetHitsData(2)->Fill( hit->GetLen()) ;//in cm
171   
172     in[0] = hit->GetSector();
173     in[1] = hit->GetPlate();
174     in[2]= hit->GetStrip();
175     in[3]= hit->GetPadx();
176     in[4]= hit->GetPadz();
177     GetMapIndeces(in,out);
178     GetHitsData(3)->Fill( out[0],out[1]) ;//hit map
179   }
180
181 }
182
183
184 //____________________________________________________________________________
185 void AliTOFQADataMakerSim::MakeHits(TTree * hitTree)
186 {
187   //
188   // make QA data from Hit Tree
189   //
190   if(!hitTree){
191     AliError("can't get the tree with TOF hits !");
192     return;
193   }     
194
195   TBranch * branch = hitTree->GetBranch("TOF") ;
196
197   if (!branch ) {
198     AliError("TOF branch in Hit Tree not found") ; 
199     return;
200   }
201
202   static TClonesArray statichits("AliTOFhitT0", 1000);
203   statichits.Clear();
204   TClonesArray *hits = &statichits;
205   static TClonesArray staticdummy("AliTOFhitT0", 1000);
206   staticdummy.Clear();
207   TClonesArray *dummy = &staticdummy;
208   branch->SetAddress(&dummy);
209   Int_t index = 0 ;  
210   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
211     branch->GetEntry(ientry) ; 
212     for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
213       AliTOFhitT0 * hit = dynamic_cast<AliTOFhitT0 *> (dummy->At(ihit)) ; 
214       new((*hits)[index]) AliTOFhitT0(*hit) ; 
215       index++ ;
216     } 
217   }     
218
219   MakeHits(hits) ; 
220
221 }
222
223 //____________________________________________________________________________
224 void AliTOFQADataMakerSim::MakeDigits(TClonesArray * digits)
225 {
226   //
227   // makes data from Digits
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=digits->GetEntriesFast();
235   if(nentries<=0){
236     GetDigitsData(0)->Fill(-1.) ; 
237   }else{
238     GetDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
239   } 
240
241   TIter next(digits) ; 
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   TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ; 
267   
268   TBranch * branch = digitTree->GetBranch("TOF") ;
269   if ( ! branch ) {
270     AliError("TOF branch in Digit Tree not found") ; 
271     return;
272   }
273   branch->SetAddress(&digits) ;
274   branch->GetEntry(0) ; 
275   MakeDigits(digits) ; 
276 }
277
278 //____________________________________________________________________________
279 void AliTOFQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
280 {
281   //
282   // makes data from SDigits
283   //
284
285   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
286   Int_t in[5];
287   Int_t out[5];
288
289   Int_t nentries=sdigits->GetEntriesFast();
290   if(nentries<=0){
291     GetSDigitsData(0)->Fill(-1.) ; 
292   }else{
293     GetSDigitsData(0)->Fill(TMath::Log10(nentries)) ; 
294   } 
295
296   TIter next(sdigits) ; 
297   AliTOFSDigit * sdigit ; 
298   while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
299     
300     for(Int_t i=0;i<sdigit->GetNDigits();i++){
301       GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
302     }
303
304     in[0] = sdigit->GetSector();
305     in[1] = sdigit->GetPlate();
306     in[2] = sdigit->GetStrip();
307     in[3] = sdigit->GetPadx();
308     in[4]= sdigit->GetPadz();
309     GetMapIndeces(in,out);
310     GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
311   }
312 }
313
314 //____________________________________________________________________________
315 void AliTOFQADataMakerSim::MakeSDigits(TTree * sdigitTree)
316 {
317   //
318   // makes data from SDigit Tree
319   //
320   TClonesArray * sdigits = new TClonesArray("AliTOFSDigit", 1000) ; 
321   
322   TBranch * branch = sdigitTree->GetBranch("TOF") ;
323   if ( ! branch ) {
324     AliError("TOF branch in SDigit Tree not found") ; 
325     return;
326   }
327   branch->SetAddress(&sdigits) ;
328   branch->GetEntry(0) ; 
329   MakeSDigits(sdigits) ; 
330 }
331
332 //____________________________________________________________________________ 
333 void AliTOFQADataMakerSim::StartOfDetectorCycle()
334 {
335   //
336   //Detector specific actions at start of cycle
337   //to be implemented  
338 }
339
340 //____________________________________________________________________________ 
341 void AliTOFQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
342 {
343   //Detector specific actions at end of cycle
344   // do the QA checking
345
346   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
347 }
348 //____________________________________________________________________________
349 void AliTOFQADataMakerSim::GetMapIndeces(Int_t* in , Int_t* out)
350 {
351   //
352   //return appropriate indeces for the theta-phi map
353   //
354
355   Int_t npadX = AliTOFGeometry::NpadX();
356   Int_t npadZ = AliTOFGeometry::NpadZ();
357   Int_t nStripA = AliTOFGeometry::NStripA();
358   Int_t nStripB = AliTOFGeometry::NStripB();
359   Int_t nStripC = AliTOFGeometry::NStripC();
360
361   Int_t isector = in[0];
362   Int_t iplate = in[1];
363   Int_t istrip = in[2];
364   Int_t ipadX = in[3]; 
365   Int_t ipadZ = in[4]; 
366   
367   Int_t stripOffset = 0;
368   switch (iplate) {
369   case 0:
370     stripOffset = 0;
371       break;
372   case 1:
373     stripOffset = nStripC;
374     break;
375   case 2:
376     stripOffset = nStripC+nStripB;
377     break;
378   case 3:
379     stripOffset = nStripC+nStripB+nStripA;
380     break;
381   case 4:
382     stripOffset = nStripC+nStripB+nStripA+nStripB;
383     break;
384   default:
385     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
386     break;
387   };
388   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
389   Int_t phiindex=npadX*isector+ipadX+1;
390   out[0]=zindex;  
391   out[1]=phiindex;  
392   
393 }