]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerSim.cxx
online hot channels introduced
[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
93 //____________________________________________________________________________ 
94 void AliPMDQADataMakerSim::InitSDigits()
95 {
96     // create SDigits histograms in SDigits subdir
97   const Bool_t expert   = kTRUE ; 
98   const Bool_t image    = kTRUE ; 
99   
100   TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
101   h0->Sumw2();
102   Add2SDigitsList(h0, 0, !expert, image);
103   
104   TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
105   h1->Sumw2();
106   Add2SDigitsList(h1, 1, !expert, image);
107   
108   TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD);# of SDigits;Entries", 500, 0., 1000.);
109   h2->Sumw2();
110   Add2SDigitsList(h2, 2, !expert, image);
111   
112   TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD);# of SDigits;Entries", 500, 0., 1000.);
113   h3->Sumw2();
114   Add2SDigitsList(h3, 3, !expert, image);
115   
116 }
117
118 //____________________________________________________________________________
119 void AliPMDQADataMakerSim::InitDigits()
120 {
121     // create Digits histograms in Digits subdir
122   const Bool_t expert   = kTRUE ; 
123   const Bool_t image    = kTRUE ; 
124   
125   TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
126   h0->Sumw2();
127   Add2DigitsList(h0, 0, !expert, image);
128   
129   TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.); 
130   h1->Sumw2();
131   Add2DigitsList(h1, 1, !expert, image);
132
133   TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD);# of Digits;Entries", 500, 0, 1000) ; 
134   h2->Sumw2();
135   Add2DigitsList(h2, 2, !expert, image);
136   
137   TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD);# of Digits;Entries", 500, 0, 1000);
138   h3->Sumw2();
139   Add2DigitsList(h3, 3, !expert, image);
140   
141 }
142
143 //____________________________________________________________________________ 
144 void AliPMDQADataMakerSim::MakeHits(TClonesArray *hits)
145 {
146     //make QA data from Hits
147
148   // Check id histograms already created for this Event Specie
149   if ( ! GetHitsData(0) )
150     InitHits() ;
151
152   Int_t premul = 0, cpvmul = 0;
153     Float_t edepkev = 0.;
154     TIter next(hits); 
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             GetHitsData(0)->Fill(edepkev);
163             premul++;
164           }
165         else if (hit->Z() < 361.5)
166           {
167             edepkev = (hit->GetEnergy())/1000.;
168             GetHitsData(1)->Fill(edepkev);
169             cpvmul++;
170           }
171     }
172
173     if(premul <= 0)
174     {
175         GetHitsData(2)->Fill(-1.); 
176     }
177     else
178     {
179         GetHitsData(2)->Fill(premul); 
180     }
181
182     if(cpvmul <= 0)
183     {
184         GetHitsData(3)->Fill(-1.); 
185     }
186     else
187     {
188         GetHitsData(3)->Fill(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     static TClonesArray statichits("AliPMDhit", 1000);
206     statichits.Clear();
207     TClonesArray *hits = &statichits;
208     static TClonesArray staticdummy("AliPMDhit", 1000);
209     staticdummy.Clear();
210     TClonesArray *dummy = &staticdummy;
211     branch->SetAddress(&dummy);
212     Int_t index = 0 ;  
213
214     for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
215         branch->GetEntry(ientry) ; 
216         for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
217             AliPMDhit * hit = dynamic_cast<AliPMDhit *> (dummy->At(ihit)) ; 
218             new((*hits)[index]) AliPMDhit(*hit) ; 
219
220             index++ ;
221         } 
222     }   
223
224     MakeHits(hits);
225
226 }
227 //____________________________________________________________________________
228 void AliPMDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
229 {
230     // makes data from SDigits
231
232   // Check id histograms already created for this Event Specie
233   if ( ! GetSDigitsData(0) )
234     InitSDigits() ;
235
236   Int_t cpvmul = 0, premul = 0;
237     Float_t edepkev = 0.;
238
239     TIter next(sdigits) ; 
240     AliPMDsdigit * sdigit ; 
241     while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
242     {
243         if(sdigit->GetDetector() == 0)
244         {
245           edepkev = (sdigit->GetCellEdep())/1000.;
246           GetSDigitsData(0)->Fill(edepkev);
247           premul++;
248         }
249         if(sdigit->GetDetector() == 1)
250         {
251           edepkev = (sdigit->GetCellEdep())/1000.;
252           GetSDigitsData(1)->Fill(edepkev);
253           cpvmul++;
254         }
255         
256     } 
257     if (premul > 0) GetSDigitsData(2)->Fill(premul);
258     if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
259     
260 }
261
262 //____________________________________________________________________________
263 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
264 {
265     // makes data from SDigit Tree
266
267     TClonesArray * sdigits = new TClonesArray("AliPMDsdigit", 1000) ; 
268     
269     TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
270     branch->SetAddress(&sdigits) ;
271
272     if ( ! branch )
273     {
274         AliWarning("PMD branch in SDigit Tree not found") ; 
275     }
276     else
277     {
278         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
279         {
280             branch->GetEntry(ient) ;
281             MakeSDigits(sdigits) ; 
282         }
283     }
284 }
285
286 //____________________________________________________________________________
287 void AliPMDQADataMakerSim::MakeDigits(TClonesArray * digits)
288 {
289     // makes data from Digits
290     
291   // Check id histograms already created for this Event Specie
292   if ( ! GetDigitsData(0) )
293     InitDigits() ;
294
295   Int_t cpvmul = 0, premul = 0;
296
297     TIter next(digits) ; 
298     AliPMDdigit * digit ; 
299     while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
300     {
301         if(digit->GetDetector() == 0)
302         {
303             GetDigitsData(0)->Fill( digit->GetADC()) ;
304             premul++;
305         }
306         if(digit->GetDetector() == 1)
307         {
308             GetDigitsData(1)->Fill( digit->GetADC());
309             cpvmul++;
310         }
311     }  
312
313     if (premul > 0) GetDigitsData(2)->Fill(premul);
314     if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
315
316
317 }
318
319 //____________________________________________________________________________
320 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
321 {
322     // makes data from Digit Tree
323
324     TClonesArray * digits = new TClonesArray("AliPMDdigit", 1000) ; 
325     
326     TBranch * branch = digitTree->GetBranch("PMDDigit") ;
327     branch->SetAddress(&digits) ;
328
329     if ( ! branch )
330     {
331         AliWarning("PMD branch in Digit Tree not found") ; 
332     }
333     else
334     {
335         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
336         {
337             branch->GetEntry(ient) ; 
338             MakeDigits(digits) ; 
339         }
340         
341     }
342 }
343
344
345 //____________________________________________________________________________ 
346
347 void AliPMDQADataMakerSim::StartOfDetectorCycle()
348 {
349     //Detector specific actions at start of cycle
350     
351 }
352 //____________________________________________________________________________ 
353
354 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
355 {
356     //Detector specific actions at end of cycle
357     // do the QA checking
358     AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
359 }