]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AD/AliADQADataMakerSim.cxx
Moving required CMake version from 2.8.4 to 2.8.8
[u/mrichter/AliRoot.git] / AD / 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", 1000, 0, 50000) ; 
98   h2->Sumw2() ;
99   Add2HitsList(h2, 2, expert, image) ;  
100  
101   //
102   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line    
103 }
104
105
106 //____________________________________________________________________________ 
107 void AliADQADataMakerSim::InitSDigits()
108 {
109   // create Digits histograms in Digits subdir
110   const Bool_t expert   = kTRUE ; 
111   const Bool_t image    = kTRUE ; 
112   
113   TH1I *fhSDigCharge[16]; 
114
115   // create SDigits histograms in SDigits subdir
116   TH1I * h0 = new TH1I("hSDigitMultiplicity", "SDigits multiplicity distribution in AD;# of Digits;Entries", 100, 0, 99) ; 
117   h0->Sumw2() ;
118   Add2DigitsList(h0, 0, !expert, image) ;
119   
120   for (Int_t i=0; i<16; i++)
121     {
122        fhSDigCharge[i] = new TH1I(Form("hSDigitCharge%d", i),Form("SDigit charges in cell %d; Time;Entries",i),1700,0.,1700);
123        
124        Add2SDigitsList(fhSDigCharge[i],i+1, !expert, image);
125        
126      }  
127   //
128   ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
129 }
130
131
132
133 //____________________________________________________________________________ 
134 void AliADQADataMakerSim::InitDigits()
135 {
136   // create Digits histograms in Digits subdir
137   const Bool_t expert   = kTRUE ; 
138   const Bool_t image    = kTRUE ; 
139   
140   TH1I *fhDigTDC[16]; 
141   TH1I *fhDigADC[16]; 
142
143   // create Digits histograms in Digits subdir
144   TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in AD;# of Digits;Entries", 100, 0, 99) ; 
145   h0->Sumw2() ;
146   Add2DigitsList(h0, 0, !expert, image) ;
147   
148   for (Int_t i=0; i<16; i++)
149     {
150        fhDigTDC[i] = new TH1I(Form("hDigitTDC%d", i),Form("Digit TDC in cell %d; TDC value;Entries",i),300,0.,149.);
151        fhDigADC[i]= new TH1I(Form("hDigitADC%d", i),Form("Digit ADC in cell %d;ADC value;Entries",i),1024,0.,1023.);
152        
153        Add2DigitsList(fhDigTDC[i],i+1, !expert, image);
154        Add2DigitsList(fhDigADC[i],i+1+16, !expert, image);  
155      }  
156   //
157   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
158 }
159
160
161 //____________________________________________________________________________
162 void AliADQADataMakerSim::MakeHits()
163 {
164         //make QA data from Hits
165
166   Int_t nhits = fHitsArray->GetEntriesFast();
167   FillHitsData(0,nhits) ;    // fills Hit multiplicity
168   for (Int_t ihit=0;ihit<nhits;ihit++) 
169     {
170            AliADhit  * ADHit   = (AliADhit*) fHitsArray->UncheckedAt(ihit);
171            if (!ADHit) {
172               AliError("The unchecked hit doesn't exist");
173               break;
174            }
175            FillHitsData(1,ADHit->GetCell());
176            FillHitsData(2,ADHit->GetNphot());
177         }
178 }
179
180
181 //____________________________________________________________________________
182
183 void AliADQADataMakerSim::MakeHits(TTree *hitTree)
184 {
185   //fills QA histos for Hits
186  if (fHitsArray)
187    fHitsArray->Clear() ; 
188   else 
189     fHitsArray = new TClonesArray("AliADhit", 1000);
190   
191   TBranch * branch = hitTree->GetBranch("AD") ;
192   if ( ! branch ) {
193     AliWarning("AD branch in Hit Tree not found") ;
194   } else {
195
196    if (branch) {
197       branch->SetAddress(&fHitsArray);
198     }else{
199       AliError("Branch AD hit not found");
200       exit(111);
201     } 
202     // Check id histograms already created for this Event Specie
203     if ( ! GetHitsData(0) )
204       InitHits() ;
205     
206     Int_t ntracks    = (Int_t) hitTree->GetEntries();
207     
208     if (ntracks<=0) return;
209     // Start loop on tracks in the hits containers
210     for (Int_t track=0; track<ntracks;track++) {
211       branch->GetEntry(track);
212       Int_t nhits = fHitsArray->GetEntriesFast();
213       FillHitsData(0,nhits) ;    // fills Hit multiplicity
214       for (Int_t ihit=0;ihit<nhits;ihit++) 
215         {
216           AliADhit  * ADHit   = (AliADhit*) fHitsArray->UncheckedAt(ihit);
217           if (!ADHit) {
218             AliError("The unchecked hit doesn't exist");
219             break;
220           }
221           FillHitsData(1,ADHit->GetCell());
222           FillHitsData(2,ADHit->GetNphot());     
223         }
224     }
225   }
226   //
227   IncEvCountCycleHits();
228   IncEvCountTotalHits();
229   //
230 }
231
232
233
234 //____________________________________________________________________________
235 void AliADQADataMakerSim::MakeSDigits(TTree *sdigitTree)
236 {
237     // makes data from Digit Tree
238         
239   if (fSDigitsArray)
240     fSDigitsArray->Clear() ; 
241   else 
242     fSDigitsArray = new TClonesArray("AliADSDigit", 1000) ; 
243
244     TBranch * branch = sdigitTree->GetBranch("ADSDigit") ;
245     if ( ! branch ) {
246          AliWarning("AD branch in SDigit Tree not found") ; 
247     } else {
248          branch->SetAddress(&fSDigitsArray) ;
249          branch->GetEntry(0) ; 
250          MakeSDigits() ; 
251     }  
252     //
253     IncEvCountCycleDigits();
254     IncEvCountTotalDigits();
255     //    
256 }
257
258 //____________________________________________________________________________
259 void AliADQADataMakerSim::MakeSDigits()
260 {
261   // makes data from SDigits
262
263   FillSDigitsData(0,fSDigitsArray->GetEntriesFast()) ; 
264   TIter next(fSDigitsArray) ; 
265     AliADSDigit *ADSDigit ; 
266     while ( (ADSDigit = dynamic_cast<AliADSDigit *>(next())) ) {
267          Int_t   PMNumber  = ADSDigit->PMNumber();       
268          FillSDigitsData(PMNumber +1, ADSDigit->GetNBins()) ; 
269          
270     }  
271 }
272
273 //____________________________________________________________________________
274 void AliADQADataMakerSim::MakeDigits()
275 {
276   // makes data from Digits
277
278   FillDigitsData(0,fDigitsArray->GetEntriesFast()) ; 
279   TIter next(fDigitsArray) ; 
280     AliADdigit *ADDigit ; 
281     while ( (ADDigit = dynamic_cast<AliADdigit *>(next())) ) {
282          Int_t   PMNumber  = ADDigit->PMNumber();    
283          FillDigitsData(PMNumber +1, ADDigit->Time()) ;    // in 100 of picoseconds
284          FillDigitsData(PMNumber +1+16, ADDigit->ADC()) ;
285     }  
286 }
287
288 //____________________________________________________________________________
289 void AliADQADataMakerSim::MakeDigits(TTree *digitTree)
290 {
291     // makes data from Digit Tree
292         
293   if (fDigitsArray)
294     fDigitsArray->Clear() ; 
295   else 
296     fDigitsArray = new TClonesArray("AliADdigit", 1000) ; 
297
298     TBranch * branch = digitTree->GetBranch("ADDigit") ;
299     if ( ! branch ) {
300          AliWarning("AD branch in Digit Tree not found") ; 
301     } else {
302          branch->SetAddress(&fDigitsArray) ;
303          branch->GetEntry(0) ; 
304          MakeDigits() ; 
305     }  
306     //
307     IncEvCountCycleDigits();
308     IncEvCountTotalDigits();
309     //    
310 }
311
312
313 //____________________________________________________________________________
314 void AliADQADataMakerSim::StartOfDetectorCycle()
315 {
316   //Detector specific actions at start of cycle
317
318 }