]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerSim.cxx
Fix bug in building local list of valid files.
[u/mrichter/AliRoot.git] / PMD / AliPMDQADataMakerSim.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 quality assurance. 
19   All data must be mergeable objects.
20   B.K. Nandi
21 */
22
23 // --- ROOT system ---
24 #include <TClonesArray.h>
25 #include <TFile.h> 
26 #include <TH1F.h> 
27 #include <TH1I.h> 
28 #include <TH2F.h> 
29 #include <TTree.h>
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34
35 #include "AliLog.h"
36 #include "AliPMDhit.h"
37 #include "AliPMDsdigit.h"
38 #include "AliPMDdigit.h"
39 #include "AliPMDQADataMakerSim.h"
40 #include "AliQAChecker.h"
41
42 ClassImp(AliPMDQADataMakerSim)
43            
44 //____________________________________________________________________________ 
45 AliPMDQADataMakerSim::AliPMDQADataMakerSim() : 
46     AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kPMD), "PMD Quality Assurance Data Maker")
47 {
48     // ctor
49 }
50
51 //____________________________________________________________________________ 
52 AliPMDQADataMakerSim::AliPMDQADataMakerSim(const AliPMDQADataMakerSim& qadm) :
53     AliQADataMakerSim()
54 {
55     //copy ctor 
56     SetName((const char*)qadm.GetName()) ; 
57     SetTitle((const char*)qadm.GetTitle()); 
58 }
59
60 //__________________________________________________________________
61 AliPMDQADataMakerSim& AliPMDQADataMakerSim::operator = (const AliPMDQADataMakerSim& qadm )
62 {
63     // Assign operator.
64     this->~AliPMDQADataMakerSim();
65     new(this) AliPMDQADataMakerSim(qadm);
66     return *this;
67 }
68  
69 //____________________________________________________________________________ 
70 void AliPMDQADataMakerSim::InitHits()
71 {
72     // create Hits histograms in Hits subdir
73   const Bool_t expert   = kTRUE ; 
74   const Bool_t image    = kTRUE ; 
75   
76   TH1F *h0 = new TH1F("hPreHitsEdep","Hits energy distribution PRE(PMD);Energy [keV];Counts", 500, 0., 500.); 
77   h0->Sumw2() ;
78   Add2HitsList(h0, 0, !expert, image) ;
79   
80   TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution CPV(PMD);Energy [keV];Counts", 500, 0., 500.); 
81   h1->Sumw2() ;
82   Add2HitsList(h1, 1, !expert, image) ;
83   
84   TH1I *h2 = new TH1I("hPreHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ; 
85   h2->Sumw2() ;
86   Add2HitsList(h2, 2, !expert, image) ;
87   
88   TH1I *h3 = new TH1I("hCpvHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ; 
89   h2->Sumw2() ;
90   Add2HitsList(h3, 3, !expert, image) ;
91   //
92   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
93 }
94
95 //____________________________________________________________________________ 
96 void AliPMDQADataMakerSim::InitSDigits()
97 {
98     // create SDigits histograms in SDigits subdir
99   const Bool_t expert   = kTRUE ; 
100   const Bool_t image    = kTRUE ; 
101   
102   TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
103   h0->Sumw2();
104   Add2SDigitsList(h0, 0, !expert, image);
105   
106   TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
107   h1->Sumw2();
108   Add2SDigitsList(h1, 1, !expert, image);
109   
110   TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD);# of SDigits;Entries", 500, 0., 1000.);
111   h2->Sumw2();
112   Add2SDigitsList(h2, 2, !expert, image);
113   
114   TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD);# of SDigits;Entries", 500, 0., 1000.);
115   h3->Sumw2();
116   Add2SDigitsList(h3, 3, !expert, image);
117   //
118   ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line  
119 }
120
121 //____________________________________________________________________________
122 void AliPMDQADataMakerSim::InitDigits()
123 {
124     // create Digits histograms in Digits subdir
125   const Bool_t expert   = kTRUE ; 
126   const Bool_t image    = kTRUE ; 
127   
128   TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
129   h0->Sumw2();
130   Add2DigitsList(h0, 0, !expert, image);
131   
132   TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.); 
133   h1->Sumw2();
134   Add2DigitsList(h1, 1, !expert, image);
135
136   TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD);# of Digits;Entries", 500, 0, 1000) ; 
137   h2->Sumw2();
138   Add2DigitsList(h2, 2, !expert, image);
139   
140   TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD);# of Digits;Entries", 500, 0, 1000);
141   h3->Sumw2();
142   Add2DigitsList(h3, 3, !expert, image);
143   //
144   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line  
145 }
146
147 //____________________________________________________________________________ 
148 void AliPMDQADataMakerSim::MakeHits()
149 {
150     //make QA data from Hits
151
152   Int_t premul = 0, cpvmul = 0;
153   Float_t edepkev = 0.;
154   TIter next(fHitsArray); 
155   AliPMDhit * hit; 
156     
157   while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
158     {
159       if (hit->Z() > 361.5)
160         {
161           edepkev = (hit->GetEnergy())/1000.;
162           FillHitsData(0,edepkev);
163           premul++;
164         }
165       else if (hit->Z() < 361.5)
166         {
167           edepkev = (hit->GetEnergy())/1000.;
168           FillHitsData(1,edepkev);
169           cpvmul++;
170         }
171     }
172   
173   if(premul <= 0)
174     {
175       FillHitsData(2,-1.); 
176     }
177   else
178     {
179       FillHitsData(2,premul); 
180     }
181   
182   if(cpvmul <= 0)
183     {
184       FillHitsData(3,-1.); 
185     }
186   else
187     {
188       FillHitsData(3,cpvmul); 
189     }
190   
191 }
192
193 //____________________________________________________________________________
194 void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
195 {
196   // make QA data from Hit Tree
197   
198   TBranch * branch = hitTree->GetBranch("PMD") ;
199   if ( ! branch )
200     {
201       AliWarning("PMD branch in Hit Tree not found") ;
202       return;
203     }
204
205   if (fHitsArray) 
206     fHitsArray->Clear() ; 
207   else 
208     fHitsArray = new TClonesArray("AliPMDhit", 1000);
209
210   branch->SetAddress(&fHitsArray);
211
212   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
213     branch->GetEntry(ientry) ; 
214     MakeHits();
215     fHitsArray->Clear() ; 
216   }     
217   //
218   IncEvCountCycleHits();
219   IncEvCountTotalHits();
220   //
221 }
222 //____________________________________________________________________________
223 void AliPMDQADataMakerSim::MakeSDigits()
224 {
225     // makes data from SDigits
226
227   Int_t cpvmul = 0, premul = 0;
228   Float_t edepkev = 0.;
229   
230   TIter next(fSDigitsArray) ; 
231   AliPMDsdigit * sdigit ; 
232   while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
233     {
234       if(sdigit->GetDetector() == 0)
235         {
236           edepkev = (sdigit->GetCellEdep())/1000.;
237           FillSDigitsData(0,edepkev);
238           premul++;
239         }
240       if(sdigit->GetDetector() == 1)
241         {
242           edepkev = (sdigit->GetCellEdep())/1000.;
243           FillSDigitsData(1,edepkev);
244           cpvmul++;
245         }
246         
247     } 
248   if (premul > 0) FillSDigitsData(2,premul);
249   if (cpvmul > 0) FillSDigitsData(3,cpvmul);
250   
251 }
252
253 //____________________________________________________________________________
254 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
255 {
256     // makes data from SDigit Tree
257   
258   if (fSDigitsArray) 
259     fSDigitsArray->Clear() ; 
260   else 
261     fSDigitsArray = new TClonesArray("AliPMDsdigit", 1000) ; 
262     
263   TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
264   if ( ! branch )
265     {
266       AliWarning("PMD branch in SDigit Tree not found") ; 
267     }
268   else
269     {
270       branch->SetAddress(&fSDigitsArray) ;
271       branch->GetEntry(0) ;
272       MakeSDigits() ; 
273     }
274   //
275   IncEvCountCycleSDigits();
276   IncEvCountTotalSDigits();
277   //
278 }
279
280 //____________________________________________________________________________
281 void AliPMDQADataMakerSim::MakeDigits()
282 {
283     // makes data from Digits
284
285   Int_t cpvmul = 0, premul = 0;
286   
287   TIter next(fDigitsArray) ; 
288   AliPMDdigit * digit ; 
289   while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
290     {
291       if(digit->GetDetector() == 0)
292         {
293           FillDigitsData(0, digit->GetADC()) ;
294           premul++;
295         }
296       if(digit->GetDetector() == 1)
297         {
298           FillDigitsData(1, digit->GetADC());
299           cpvmul++;
300         }
301     }  
302   
303   if (premul > 0) FillDigitsData(2,premul);
304   if (cpvmul > 0) FillDigitsData(3,cpvmul);
305   
306   
307 }
308
309 //____________________________________________________________________________
310 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
311 {
312     // makes data from Digit Tree
313
314   if (fDigitsArray) 
315     fDigitsArray->Clear() ; 
316   else 
317     fDigitsArray = new TClonesArray("AliPMDdigit", 1000) ; 
318   
319   TBranch * branch = digitTree->GetBranch("PMDDigit") ;
320
321   if ( ! branch )
322     {
323       AliWarning("PMD branch in Digit Tree not found") ; 
324     }
325   else
326     {
327       branch->SetAddress(&fDigitsArray) ;
328       for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
329         {
330           branch->GetEntry(ient) ; 
331           MakeDigits() ; 
332           fDigitsArray->Clear() ; 
333           
334         }
335       
336     }
337   //
338   IncEvCountCycleDigits();
339   IncEvCountTotalDigits();
340   //
341 }
342
343
344 //____________________________________________________________________________ 
345
346 void AliPMDQADataMakerSim::StartOfDetectorCycle()
347 {
348     //Detector specific actions at start of cycle
349     
350 }
351 //____________________________________________________________________________ 
352
353 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
354 {
355     //Detector specific actions at end of cycle
356     // do the QA checking
357     AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
358 }