]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerSim.cxx
Fixes for bug #49914: Compilation breaks in trunk, and bug #48629: Trunk cannot read...
[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
74     TH1F *h0 = new TH1F("hPreHitsEdep","Hits energy distribution in (keV)PRE(PMD)", 500, 0., 500.); 
75     h0->Sumw2() ;
76     Add2HitsList(h0, 0) ;
77
78     TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution in (keV)CPV(PMD)", 500, 0., 500.); 
79     h1->Sumw2() ;
80     Add2HitsList(h1, 1) ;
81
82     TH1I *h2 = new TH1I("hPreHitsMult","Hits multiplicity distribution in PRE(PMD)", 500, 0, 3000) ; 
83     h2->Sumw2() ;
84     Add2HitsList(h2, 2) ;
85
86     TH1I *h3 = new TH1I("hCpvHitsMult","Hits multiplicity distribution in PRE(PMD)", 500, 0, 3000) ; 
87     h2->Sumw2() ;
88     Add2HitsList(h3, 3) ;
89 }
90
91 //____________________________________________________________________________ 
92 void AliPMDQADataMakerSim::InitSDigits()
93 {
94     // create SDigits histograms in SDigits subdir
95
96     TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD)", 500, 0., 500.);
97     h0->Sumw2();
98     Add2SDigitsList(h0, 0);
99
100     TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD)", 500, 0., 500.);
101     h1->Sumw2();
102     Add2SDigitsList(h1, 1);
103
104     TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD)", 500, 0., 1000.);
105     h2->Sumw2();
106     Add2SDigitsList(h2, 2);
107
108     TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD)", 500, 0., 1000.);
109     h3->Sumw2();
110     Add2SDigitsList(h3, 3);
111
112 }
113
114 //____________________________________________________________________________
115 void AliPMDQADataMakerSim::InitDigits()
116 {
117     // create Digits histograms in Digits subdir
118
119     TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD)", 100, 0., 2000.);
120     h0->Sumw2();
121     Add2DigitsList(h0, 0);
122
123     TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD)", 100, 0., 2000.); 
124     h1->Sumw2();
125     Add2DigitsList(h1, 1);
126
127     TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD)", 500, 0, 1000) ; 
128     h2->Sumw2();
129     Add2DigitsList(h2, 2);
130
131     TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD)", 500, 0, 1000);
132     h3->Sumw2();
133     Add2DigitsList(h3, 3);
134
135 }
136
137 //____________________________________________________________________________ 
138 void AliPMDQADataMakerSim::MakeHits(TClonesArray *hits)
139 {
140     //make QA data from Hits
141
142     Int_t premul = 0, cpvmul = 0;
143     Float_t edepkev = 0.;
144     TIter next(hits); 
145     AliPMDhit * hit; 
146     
147     while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
148       {
149         if (hit->Z() > 361.5)
150           {
151             edepkev = (hit->GetEnergy())/1000.;
152             GetHitsData(0)->Fill(edepkev);
153             premul++;
154           }
155         else if (hit->Z() < 361.5)
156           {
157             edepkev = (hit->GetEnergy())/1000.;
158             GetHitsData(1)->Fill(edepkev);
159             cpvmul++;
160           }
161     }
162
163     if(premul <= 0)
164     {
165         GetHitsData(2)->Fill(-1.); 
166     }
167     else
168     {
169         GetHitsData(2)->Fill(premul); 
170     }
171
172     if(cpvmul <= 0)
173     {
174         GetHitsData(3)->Fill(-1.); 
175     }
176     else
177     {
178         GetHitsData(3)->Fill(cpvmul); 
179     }
180
181 }
182
183 //____________________________________________________________________________
184 void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
185 {
186     // make QA data from Hit Tree
187
188     TBranch * branch = hitTree->GetBranch("PMD") ;
189     if ( ! branch )
190     {
191         AliWarning("PMD branch in Hit Tree not found") ;
192         return;
193     }
194
195     static TClonesArray statichits("AliPMDhit", 1000);
196     statichits.Clear();
197     TClonesArray *hits = &statichits;
198     static TClonesArray staticdummy("AliPMDhit", 1000);
199     staticdummy.Clear();
200     TClonesArray *dummy = &staticdummy;
201     branch->SetAddress(&dummy);
202     Int_t index = 0 ;  
203
204     for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
205         branch->GetEntry(ientry) ; 
206         for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
207             AliPMDhit * hit = dynamic_cast<AliPMDhit *> (dummy->At(ihit)) ; 
208             new((*hits)[index]) AliPMDhit(*hit) ; 
209
210             index++ ;
211         } 
212     }   
213
214     MakeHits(hits);
215
216 }
217 //____________________________________________________________________________
218 void AliPMDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
219 {
220     // makes data from SDigits
221
222     Int_t cpvmul = 0, premul = 0;
223     Float_t edepkev = 0.;
224
225     TIter next(sdigits) ; 
226     AliPMDsdigit * sdigit ; 
227     while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
228     {
229         if(sdigit->GetDetector() == 0)
230         {
231           edepkev = (sdigit->GetCellEdep())/1000.;
232           GetSDigitsData(0)->Fill(edepkev);
233           premul++;
234         }
235         if(sdigit->GetDetector() == 1)
236         {
237           edepkev = (sdigit->GetCellEdep())/1000.;
238           GetSDigitsData(1)->Fill(edepkev);
239           cpvmul++;
240         }
241         
242     } 
243     if (premul > 0) GetSDigitsData(2)->Fill(premul);
244     if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
245     
246 }
247
248 //____________________________________________________________________________
249 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
250 {
251     // makes data from SDigit Tree
252
253     TClonesArray * sdigits = new TClonesArray("AliPMDsdigit", 1000) ; 
254     
255     TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
256     branch->SetAddress(&sdigits) ;
257
258     if ( ! branch )
259     {
260         AliWarning("PMD branch in SDigit Tree not found") ; 
261     }
262     else
263     {
264         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
265         {
266             branch->GetEntry(ient) ;
267             MakeSDigits(sdigits) ; 
268         }
269     }
270 }
271
272 //____________________________________________________________________________
273 void AliPMDQADataMakerSim::MakeDigits(TClonesArray * digits)
274 {
275     // makes data from Digits
276     
277     Int_t cpvmul = 0, premul = 0;
278
279     TIter next(digits) ; 
280     AliPMDdigit * digit ; 
281     while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
282     {
283         if(digit->GetDetector() == 0)
284         {
285             GetDigitsData(0)->Fill( digit->GetADC()) ;
286             premul++;
287         }
288         if(digit->GetDetector() == 1)
289         {
290             GetDigitsData(1)->Fill( digit->GetADC());
291             cpvmul++;
292         }
293     }  
294
295     if (premul > 0) GetDigitsData(2)->Fill(premul);
296     if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
297
298
299 }
300
301 //____________________________________________________________________________
302 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
303 {
304     // makes data from Digit Tree
305
306     TClonesArray * digits = new TClonesArray("AliPMDdigit", 1000) ; 
307     
308     TBranch * branch = digitTree->GetBranch("PMDDigit") ;
309     branch->SetAddress(&digits) ;
310
311     if ( ! branch )
312     {
313         AliWarning("PMD branch in Digit Tree not found") ; 
314     }
315     else
316     {
317         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
318         {
319             branch->GetEntry(ient) ; 
320             MakeDigits(digits) ; 
321         }
322         
323     }
324 }
325
326
327 //____________________________________________________________________________ 
328
329 void AliPMDQADataMakerSim::StartOfDetectorCycle()
330 {
331     //Detector specific actions at start of cycle
332     
333 }
334 //____________________________________________________________________________ 
335
336 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
337 {
338     //Detector specific actions at end of cycle
339     // do the QA checking
340     AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
341 }