]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerSim.cxx
added comments to new data members (Leticia)
[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 per event;TOF hit number;Counts ",101, -1., 100.) ; 
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", 25000, 0., 610.) ; 
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", 700, 0., 700) ; 
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 //____________________________________________________________________________ 
105 void AliTOFQADataMakerSim::InitDigits()
106 {
107   //
108   // create Digits histograms in Digits subdir
109   //
110
111   const Bool_t expert   = kTRUE ; 
112   const Bool_t image    = kTRUE ; 
113   
114   TH1F * h0 = new TH1F("hTOFDigits",    "Number of TOF Digit per event;TOF digit number;Counts ",101, -1., 100.) ;
115   h0->Sumw2() ;
116   Add2DigitsList(h0, 0, !expert, image) ;
117
118   TH1F * h1  = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 25000, 0., 610.) ; 
119   h1->Sumw2() ;
120   Add2DigitsList(h1, 1, !expert, image) ;
121
122   TH1F * h2  = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns);Digitized ToT time [ns];Counts", 1000, 0., 48.8) ; 
123   h2->Sumw2() ;
124   Add2DigitsList(h2, 2, !expert, image) ;
125
126   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) ; 
127   h3->Sumw2() ;
128   h3->GetYaxis()->SetTitleOffset(1.15);
129   Add2DigitsList(h3, 3, !expert, image) ;
130
131 }
132
133 //____________________________________________________________________________ 
134 void AliTOFQADataMakerSim::InitSDigits()
135 {
136   //
137   // create SDigits histograms in SDigits subdir
138   //
139
140   const Bool_t expert   = kTRUE ; 
141   const Bool_t image    = kTRUE ; 
142   
143   TH1F * h0 = new TH1F("hTOFSDigits",    "Number of TOF SDigits per event;TOF sdigit number;Counts ",101, -1., 100.) ;
144   h0->Sumw2() ;
145   Add2SDigitsList(h0, 0, !expert, image) ;
146
147   TH1F * h1  = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns);SDigitized TOF time [ns];Counts", 25000, 0., 610) ; 
148   h1->Sumw2() ;
149   Add2SDigitsList(h1, 1, !expert, image) ;
150
151   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) ; 
152   h2->Sumw2() ;
153   h2->GetYaxis()->SetTitleOffset(1.15);
154   Add2SDigitsList(h2, 2, !expert, image) ;
155
156 }
157
158 //____________________________________________________________________________
159 void AliTOFQADataMakerSim::MakeHits()
160 {
161   //
162   //make QA data from Hits
163   //
164
165   Int_t in[5];
166   Int_t out[5];
167
168   Int_t nentries= fHitsArray->GetEntriesFast();
169   if(nentries<=0) {
170     GetHitsData(0)->Fill(-1.) ; 
171   } else{
172     GetHitsData(0)->Fill(nentries) ; 
173   }
174   TIter next(fHitsArray) ; 
175   AliTOFhitT0 * hit ; 
176   while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
177
178     GetHitsData(1)->Fill( hit->GetTof()*1.E9) ;//in ns
179     GetHitsData(2)->Fill( hit->GetLen()) ;//in cm
180   
181     in[0] = hit->GetSector();
182     in[1] = hit->GetPlate();
183     in[2]= hit->GetStrip();
184     in[3]= hit->GetPadx();
185     in[4]= hit->GetPadz();
186     GetMapIndeces(in,out);
187     GetHitsData(3)->Fill( out[0],out[1]) ;//hit map
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   if (fHitsArray) 
211     fHitsArray->Clear() ; 
212   else 
213     fHitsArray = new TClonesArray("AliTOFhitT0", 1000) ;
214   
215   branch->SetAddress(&fHitsArray);
216   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
217     branch->GetEntry(ientry) ; 
218     MakeHits() ; 
219     fHitsArray->Clear() ; 
220   }     
221 }
222
223 //____________________________________________________________________________
224 void AliTOFQADataMakerSim::MakeDigits()
225 {
226   //
227   // makes data from Digits
228   //
229   
230   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
231   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
232   Int_t in[5];
233   Int_t out[5];
234
235   Int_t nentries=fDigitsArray->GetEntriesFast();
236   if(nentries<=0){
237     GetDigitsData(0)->Fill(-1.) ; 
238   }else{
239     GetDigitsData(0)->Fill(nentries) ; 
240   } 
241
242   TIter next(fDigitsArray) ; 
243   AliTOFdigit * digit ; 
244   while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
245     
246     GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
247     GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
248
249     in[0] = digit->GetSector();
250     in[1] = digit->GetPlate();
251     in[2] = digit->GetStrip();
252     in[3] = digit->GetPadx();
253     in[4]= digit->GetPadz();
254     GetMapIndeces(in,out);
255     GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
256   }
257
258 }
259
260
261 //____________________________________________________________________________
262 void AliTOFQADataMakerSim::MakeDigits(TTree * digitTree)
263 {
264   //
265   // makes data from Digit Tree
266   //
267   if (fDigitsArray) 
268     fDigitsArray->Clear() ; 
269   else 
270     fDigitsArray = new TClonesArray("AliTOFdigit", 1000) ; 
271   
272   TBranch * branch = digitTree->GetBranch("TOF") ;
273   if ( ! branch ) {
274     AliError("TOF branch in Digit Tree not found") ; 
275     return;
276   }
277   branch->SetAddress(&fDigitsArray) ;
278   branch->GetEntry(0) ; 
279   MakeDigits() ; 
280 }
281
282 //____________________________________________________________________________
283 void AliTOFQADataMakerSim::MakeSDigits()
284 {
285   //
286   // makes data from SDigits
287   //
288   
289   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
290   Int_t in[5];
291   Int_t out[5];
292
293   Int_t nentries=fSDigitsArray->GetEntriesFast();
294   if(nentries<=0){
295     GetSDigitsData(0)->Fill(-1.) ; 
296   }else{
297     GetSDigitsData(0)->Fill(nentries) ; 
298   } 
299
300   TIter next(fSDigitsArray) ; 
301   AliTOFSDigit * sdigit ; 
302   while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
303     
304     for(Int_t i=0;i<sdigit->GetNDigits();i++){
305       GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
306     }
307
308     in[0] = sdigit->GetSector();
309     in[1] = sdigit->GetPlate();
310     in[2] = sdigit->GetStrip();
311     in[3] = sdigit->GetPadx();
312     in[4]= sdigit->GetPadz();
313     GetMapIndeces(in,out);
314     GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
315   }
316 }
317
318 //____________________________________________________________________________
319 void AliTOFQADataMakerSim::MakeSDigits(TTree * sdigitTree)
320 {
321   //
322   // makes data from SDigit Tree
323   //
324   if (fSDigitsArray) 
325     fSDigitsArray->Clear() ; 
326   else 
327     fSDigitsArray = new TClonesArray("AliTOFSDigit", 1000) ; 
328   
329   TBranch * branch = sdigitTree->GetBranch("TOF") ;
330   if ( ! branch ) {
331     AliError("TOF branch in SDigit Tree not found") ; 
332     return;
333   }
334   branch->SetAddress(&fSDigitsArray) ;
335   branch->GetEntry(0) ; 
336   MakeSDigits() ; 
337 }
338
339 //____________________________________________________________________________ 
340 void AliTOFQADataMakerSim::StartOfDetectorCycle()
341 {
342   //
343   //Detector specific actions at start of cycle
344   //to be implemented  
345 }
346
347 //____________________________________________________________________________ 
348 void AliTOFQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
349 {
350   //Detector specific actions at end of cycle
351   // do the QA checking
352
353   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
354 }
355 //____________________________________________________________________________
356 void AliTOFQADataMakerSim::GetMapIndeces(Int_t* in , Int_t* out)
357 {
358   //
359   //return appropriate indeces for the theta-phi map
360   //
361
362   Int_t npadX = AliTOFGeometry::NpadX();
363   Int_t npadZ = AliTOFGeometry::NpadZ();
364   Int_t nStripA = AliTOFGeometry::NStripA();
365   Int_t nStripB = AliTOFGeometry::NStripB();
366   Int_t nStripC = AliTOFGeometry::NStripC();
367
368   Int_t isector = in[0];
369   Int_t iplate = in[1];
370   Int_t istrip = in[2];
371   Int_t ipadX = in[3]; 
372   Int_t ipadZ = in[4]; 
373   
374   Int_t stripOffset = 0;
375   switch (iplate) {
376   case 0:
377     stripOffset = 0;
378       break;
379   case 1:
380     stripOffset = nStripC;
381     break;
382   case 2:
383     stripOffset = nStripC+nStripB;
384     break;
385   case 3:
386     stripOffset = nStripC+nStripB+nStripA;
387     break;
388   case 4:
389     stripOffset = nStripC+nStripB+nStripA+nStripB;
390     break;
391   default:
392     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
393     break;
394   };
395   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
396   Int_t phiindex=npadX*isector+ipadX+1;
397   out[0]=zindex;  
398   out[1]=phiindex;  
399   
400 }