]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AD/ADsim/AliADQADataMakerSim.cxx
Changes for Root6 (Mikolaj)
[u/mrichter/AliRoot.git] / AD / ADsim / AliADQADataMakerSim.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 /* $Id: AliADQADataMakerSim.cxx 23123 2007-12-18 09:08:18Z hristov $ */
18
19 //---
20 //  Produces the data needed to calculate the quality assurance. 
21 //  All data must be mergeable objects.
22 //  Author : BC
23 //---
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TDirectory.h>
30 // --- Standard library ---
31
32 // --- AliRoot header files ---
33 #include "AliESDEvent.h"
34 #include "AliLog.h"
35 #include "AliADdigit.h"
36 #include "AliADSDigit.h" 
37 #include "AliADhit.h"
38 #include "AliADQADataMakerSim.h"
39 #include "AliQAChecker.h"
40
41 ClassImp(AliADQADataMakerSim)
42            
43 //____________________________________________________________________________ 
44   AliADQADataMakerSim::AliADQADataMakerSim() : 
45   AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kAD), "AD Quality Assurance Data Maker")
46
47 {
48   // constructor
49
50   
51 }
52
53 //____________________________________________________________________________ 
54 AliADQADataMakerSim::AliADQADataMakerSim(const AliADQADataMakerSim& qadm) :
55   AliQADataMakerSim() 
56 {
57   //copy constructor 
58
59   SetName((const char*)qadm.GetName()) ; 
60   SetTitle((const char*)qadm.GetTitle()); 
61 }
62
63 //__________________________________________________________________
64 AliADQADataMakerSim& AliADQADataMakerSim::operator = (const AliADQADataMakerSim& qadm )
65 {
66   // Assign operator.
67   this->~AliADQADataMakerSim();
68   new(this) AliADQADataMakerSim(qadm);
69   return *this;
70 }
71 //____________________________________________________________________________
72 void AliADQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
73 {
74   //Detector specific actions at end of cycle
75   // do the QA checking
76   ResetEventTrigClasses();
77   AliQAChecker::Instance()->Run(AliQAv1::kAD, task, list) ;
78 }
79
80  
81 //____________________________________________________________________________ 
82 void AliADQADataMakerSim::InitHits()
83 {
84  
85   // create Hits histograms in Hits subdir
86   const Bool_t expert   = kTRUE ; 
87   const Bool_t image    = kTRUE ; 
88   
89   TH1I * h0 = new TH1I("hHitMultiplicity", "Hit multiplicity distribution in AD;# of Hits;Entries", 300, 0, 299) ; 
90   h0->Sumw2() ;
91   Add2HitsList(h0, 0, !expert, image) ;  
92   
93   TH1I * h1 = new TH1I("hHitCellNumber", "Hit cell distribution in AD;Cell;# of Hits", 16, 0, 16) ; 
94   h1->Sumw2() ;
95   Add2HitsList(h1, 1, !expert, image) ; 
96   
97   TH1I * h2 = new TH1I("hHitNPhotons", "Number of photons per hit in AD;# of Photons;Entries", 100000, 0, 100000) ; 
98   h2->Sumw2() ;
99   Add2HitsList(h2, 2, !expert, image) ;
100   
101   TH2I * h3 = new TH2I("hCellNPhotons", "Number of photons per cell in AD;Cell;# of Photons", 16, 0, 16, 100000, 0, 100000) ; 
102   h2->Sumw2() ;
103   Add2HitsList(h3, 3, !expert, image) ;
104   
105   TH2D * h4 = new TH2D("hCellTof", "Time of flight per cell in AD;Cell;Time of Flight [ns]", 16, 0, 16, 6000,40,100) ; 
106   h2->Sumw2() ;
107   Add2HitsList(h4, 4, !expert, image) ;  
108  
109   //
110   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line    
111 }
112
113
114 //____________________________________________________________________________ 
115 void AliADQADataMakerSim::InitSDigits()
116 {
117   // create Digits histograms in Digits subdir
118   const Bool_t expert   = kTRUE ; 
119   const Bool_t image    = kTRUE ; 
120   
121   TH1I *fhSDigCharge[16]; 
122
123   // create SDigits histograms in SDigits subdir
124   TH1I * h0 = new TH1I("hSDigitMultiplicity", "SDigits multiplicity distribution in AD;# of Digits;Entries", 17,-0.5,16.5) ; 
125   h0->Sumw2() ;
126   Add2SDigitsList(h0, 0, !expert, image) ;
127   
128   TH2D * h1 = new TH2D("hSDigitChargePerPM", "SDigits total amplified charge per PM in AD;PM number;Charge [C]", 16, 0, 16, 1000,1e-13, 1e-9); 
129   h1->Sumw2() ;
130   Add2SDigitsList(h1, 1, !expert, image) ;
131   
132   //
133   ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
134 }
135
136
137
138 //____________________________________________________________________________ 
139 void AliADQADataMakerSim::InitDigits()
140 {
141   // create Digits histograms in Digits subdir
142   const Bool_t expert   = kTRUE ; 
143   const Bool_t image    = kTRUE ; 
144
145   // create Digits histograms in Digits subdir
146   TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in AD;# of Digits;Entries", 17,-0.5,16.5) ; 
147   h0->Sumw2() ;
148   Add2DigitsList(h0, 0, !expert, image) ;
149      
150   TH2D * h1 = new TH2D("hDigitLeadingTimePerPM", "Leading time distribution per PM in AD;PM number;Leading Time [ns]",16,0,16, 1000, 200, 300); 
151   h1->Sumw2() ;
152   Add2DigitsList(h1, 1, !expert, image) ; 
153   
154   TH2D * h2 = new TH2D("hDigitTimeWidthPerPM", "Time width distribution per PM in AD;PM number;Time width [ns]",16,0,16, 1000, 0, 100); 
155   h2->Sumw2() ;
156   Add2DigitsList(h2, 2, !expert, image) ;
157   
158   TH2I * h3 = new TH2I("hDigitChargePerClockPerPM", "Charge array per PM in AD;PM number; Clock",16,0,16,21, -10.5, 10.5);
159   h3->Sumw2();
160   Add2DigitsList(h3, 3, !expert, image) ;
161   
162   TH1I * h4 = new TH1I("hDigitBBflagsAD","Number of BB flags in AD; # of BB flags; Entries",17,-0.5,16.5);
163   h4->Sumw2();
164   Add2DigitsList(h4, 4, !expert, image) ;
165   
166   TH1I * h5 = new TH1I("hDigitBBflagsADA","Number of BB flags in ADA; # of BB flags; Entries",9,-0.5,8.5);
167   h5->Sumw2();
168   Add2DigitsList(h5, 5, !expert, image) ;
169   
170   TH1I * h6 = new TH1I("hDigitBBflagsADC","Number of BB flags in ADC; # of BB flags; Entries",9,-0.5,8.5);
171   h6->Sumw2();
172   Add2DigitsList(h6, 6, !expert, image) ;
173   
174   TH2D * h7 = new TH2D("hDigitTotalChargePerPM", "Total Charge per PM in AD;PM number; Charge [ADC counts]",16,0,16,10000,0,10000);
175   h7->Sumw2();
176   Add2DigitsList(h7, 7, !expert, image) ;
177   
178   TH2I * h8 = new TH2I("hDigitMaxChargeClockPerPM", "Clock with maximum charge per PM in AD;PM number; Clock ",16,0,16,21, -10.5, 10.5);
179   h8->Sumw2();
180   Add2DigitsList(h8, 8, !expert, image) ;
181    
182   //
183   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
184 }
185
186
187 //____________________________________________________________________________
188 void AliADQADataMakerSim::MakeHits()
189 {
190         //make QA data from Hits
191
192   Int_t nhits = fHitsArray->GetEntriesFast();
193   FillHitsData(0,nhits) ;    // fills Hit multiplicity
194   for (Int_t ihit=0;ihit<nhits;ihit++) 
195     {
196            AliADhit  * ADHit   = (AliADhit*) fHitsArray->UncheckedAt(ihit);
197            if (!ADHit) {
198               AliError("The unchecked hit doesn't exist");
199               break;
200            }
201            FillHitsData(1,ADHit->GetCell());
202            FillHitsData(2,ADHit->GetNphot());
203            FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
204            FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());
205         }
206 }
207
208
209 //____________________________________________________________________________
210
211 void AliADQADataMakerSim::MakeHits(TTree *hitTree)
212 {
213   //fills QA histos for Hits
214  if (fHitsArray)
215    fHitsArray->Clear() ; 
216   else 
217     fHitsArray = new TClonesArray("AliADhit", 1000);
218   
219   TBranch * branch = hitTree->GetBranch("AD") ;
220   if ( ! branch ) {
221     AliWarning("AD branch in Hit Tree not found") ;
222   } else {
223
224    if (branch) {
225       branch->SetAddress(&fHitsArray);
226     }else{
227       AliError("Branch AD hit not found");
228       exit(111);
229     } 
230     // Check id histograms already created for this Event Specie
231     if ( ! GetHitsData(0) )
232       InitHits() ;
233     
234     Int_t ntracks    = (Int_t) hitTree->GetEntries();
235     
236     if (ntracks<=0) return;
237     // Start loop on tracks in the hits containers
238     for (Int_t track=0; track<ntracks;track++) {
239       branch->GetEntry(track);
240       Int_t nhits = fHitsArray->GetEntriesFast();
241       FillHitsData(0,nhits) ;    // fills Hit multiplicity
242       for (Int_t ihit=0;ihit<nhits;ihit++) 
243         {
244           AliADhit  * ADHit   = (AliADhit*) fHitsArray->UncheckedAt(ihit);
245           if (!ADHit) {
246             AliError("The unchecked hit doesn't exist");
247             break;
248           }
249           FillHitsData(1,ADHit->GetCell());
250           FillHitsData(2,ADHit->GetNphot());
251           FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
252           FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());      
253         }
254     }
255   }
256   //
257   IncEvCountCycleHits();
258   IncEvCountTotalHits();
259   //
260 }
261
262
263
264 //____________________________________________________________________________
265 void AliADQADataMakerSim::MakeSDigits(TTree *sdigitTree)
266 {
267     // makes data from Digit Tree
268         
269   if (fSDigitsArray)
270     fSDigitsArray->Clear() ; 
271   else 
272     fSDigitsArray = new TClonesArray("AliADSDigit", 1000) ; 
273
274     TBranch * branch = sdigitTree->GetBranch("ADSDigit") ;
275     if ( ! branch ) {
276          AliWarning("AD branch in SDigit Tree not found") ; 
277     } else {
278          branch->SetAddress(&fSDigitsArray) ;
279          branch->GetEntry(0) ; 
280          MakeSDigits() ; 
281     }  
282     //
283     IncEvCountCycleDigits();
284     IncEvCountTotalDigits();
285     //    
286 }
287
288 //____________________________________________________________________________
289 void AliADQADataMakerSim::MakeSDigits()
290 {
291   // makes data from SDigits
292
293   FillSDigitsData(0,fSDigitsArray->GetEntriesFast()) ; 
294   TIter next(fSDigitsArray) ; 
295     AliADSDigit *ADSDigit ; 
296     while ( (ADSDigit = dynamic_cast<AliADSDigit *>(next())) ) {
297          Int_t   PMNumber  = ADSDigit->PMNumber();
298          Int_t   Nbins = ADSDigit->GetNBins();
299          Double_t totCharge = 0;
300         
301          for(Int_t i = 0; i<Nbins; i++)totCharge += ADSDigit->GetCharges()[i];       
302          FillSDigitsData(1, PMNumber, totCharge) ;
303     }  
304 }
305
306 //____________________________________________________________________________
307 void AliADQADataMakerSim::MakeDigits()
308 {
309   // makes data from Digits
310
311   FillDigitsData(0,fDigitsArray->GetEntriesFast()) ; 
312   TIter next(fDigitsArray) ; 
313     AliADdigit *ADDigit ; 
314     Int_t nBBflagsADA = 0;
315     Int_t nBBflagsADC = 0;
316     
317     while ( (ADDigit = dynamic_cast<AliADdigit *>(next())) ) {
318          Int_t totCharge = 0;
319          Int_t   PMNumber  = ADDigit->PMNumber();
320
321          if(PMNumber<8 && ADDigit->GetBBflag()) nBBflagsADC++;
322          if(PMNumber>7 && ADDigit->GetBBflag()) nBBflagsADA++;
323          
324          Short_t adc[21];
325          for(Int_t iClock=0; iClock<21; iClock++) { 
326          adc[iClock]= ADDigit->ChargeADC(iClock);
327          FillDigitsData(3, PMNumber,(float)iClock-10,(float)adc[iClock]);
328          totCharge += adc[iClock];
329          }
330             
331          FillDigitsData(1,PMNumber,ADDigit->Time()); 
332          FillDigitsData(2,PMNumber,ADDigit->Width());
333          FillDigitsData(7,PMNumber,totCharge);
334          FillDigitsData(8,PMNumber,TMath::LocMax(21,adc)-10); 
335          
336     }
337     FillDigitsData(4,nBBflagsADA+nBBflagsADC);
338     FillDigitsData(5,nBBflagsADA);
339     FillDigitsData(6,nBBflagsADC);  
340 }
341
342 //____________________________________________________________________________
343 void AliADQADataMakerSim::MakeDigits(TTree *digitTree)
344 {
345     // makes data from Digit Tree
346         
347   if (fDigitsArray)
348     fDigitsArray->Clear() ; 
349   else 
350     fDigitsArray = new TClonesArray("AliADdigit", 1000) ; 
351
352     TBranch * branch = digitTree->GetBranch("ADDigit") ;
353     if ( ! branch ) {
354          AliWarning("AD branch in Digit Tree not found") ; 
355     } else {
356          branch->SetAddress(&fDigitsArray) ;
357          branch->GetEntry(0) ; 
358          MakeDigits() ; 
359     }  
360     //
361     IncEvCountCycleDigits();
362     IncEvCountTotalDigits();
363     //    
364 }
365
366
367 //____________________________________________________________________________
368 void AliADQADataMakerSim::StartOfDetectorCycle()
369 {
370   //Detector specific actions at start of cycle
371
372 }