Line data Source code
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 12 : ClassImp(AliPMDQADataMakerSim)
43 :
44 : //____________________________________________________________________________
45 : AliPMDQADataMakerSim::AliPMDQADataMakerSim() :
46 3 : AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kPMD), "PMD Quality Assurance Data Maker")
47 5 : {
48 : // ctor
49 2 : }
50 :
51 : //____________________________________________________________________________
52 : AliPMDQADataMakerSim::AliPMDQADataMakerSim(const AliPMDQADataMakerSim& qadm) :
53 0 : AliQADataMakerSim()
54 0 : {
55 : //copy ctor
56 0 : SetName((const char*)qadm.GetName()) ;
57 0 : SetTitle((const char*)qadm.GetTitle());
58 0 : }
59 :
60 : //__________________________________________________________________
61 : AliPMDQADataMakerSim& AliPMDQADataMakerSim::operator = (const AliPMDQADataMakerSim& qadm )
62 : {
63 : // Assign operator.
64 0 : this->~AliPMDQADataMakerSim();
65 0 : new(this) AliPMDQADataMakerSim(qadm);
66 0 : return *this;
67 0 : }
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 2 : TH1F *h0 = new TH1F("hPreHitsEdep","Hits energy distribution PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
77 1 : h0->Sumw2() ;
78 1 : Add2HitsList(h0, 0, !expert, image) ;
79 :
80 1 : TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
81 1 : h1->Sumw2() ;
82 1 : Add2HitsList(h1, 1, !expert, image) ;
83 :
84 1 : TH1I *h2 = new TH1I("hPreHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ;
85 1 : h2->Sumw2() ;
86 1 : Add2HitsList(h2, 2, !expert, image) ;
87 :
88 1 : TH1I *h3 = new TH1I("hCpvHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ;
89 1 : h2->Sumw2() ;
90 1 : Add2HitsList(h3, 3, !expert, image) ;
91 : //
92 1 : ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
93 1 : }
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 2 : TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
103 1 : h0->Sumw2();
104 1 : Add2SDigitsList(h0, 0, !expert, image);
105 :
106 1 : TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
107 1 : h1->Sumw2();
108 1 : Add2SDigitsList(h1, 1, !expert, image);
109 :
110 1 : TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD);# of SDigits;Entries", 500, 0., 1000.);
111 1 : h2->Sumw2();
112 1 : Add2SDigitsList(h2, 2, !expert, image);
113 :
114 1 : TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD);# of SDigits;Entries", 500, 0., 1000.);
115 1 : h3->Sumw2();
116 1 : Add2SDigitsList(h3, 3, !expert, image);
117 : //
118 1 : ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
119 1 : }
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 2 : TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
129 1 : h0->Sumw2();
130 1 : Add2DigitsList(h0, 0, !expert, image);
131 :
132 1 : TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
133 1 : h1->Sumw2();
134 1 : Add2DigitsList(h1, 1, !expert, image);
135 :
136 1 : TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD);# of Digits;Entries", 500, 0, 1000) ;
137 1 : h2->Sumw2();
138 1 : Add2DigitsList(h2, 2, !expert, image);
139 :
140 1 : TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD);# of Digits;Entries", 500, 0, 1000);
141 1 : h3->Sumw2();
142 1 : Add2DigitsList(h3, 3, !expert, image);
143 : //
144 1 : ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
145 1 : }
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 224 : TIter next(fHitsArray);
155 : AliPMDhit * hit;
156 :
157 2336 : while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
158 : {
159 944 : if (hit->Z() > 361.5)
160 : {
161 393 : edepkev = (hit->GetEnergy())/1000.;
162 393 : FillHitsData(0,edepkev);
163 393 : premul++;
164 393 : }
165 158 : else if (hit->Z() < 361.5)
166 : {
167 79 : edepkev = (hit->GetEnergy())/1000.;
168 79 : FillHitsData(1,edepkev);
169 79 : cpvmul++;
170 79 : }
171 : }
172 :
173 224 : if(premul <= 0)
174 : {
175 200 : FillHitsData(2,-1.);
176 : }
177 : else
178 : {
179 24 : FillHitsData(2,premul);
180 : }
181 :
182 224 : if(cpvmul <= 0)
183 : {
184 200 : FillHitsData(3,-1.);
185 : }
186 : else
187 : {
188 24 : FillHitsData(3,cpvmul);
189 : }
190 :
191 112 : }
192 :
193 : //____________________________________________________________________________
194 : void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
195 : {
196 : // make QA data from Hit Tree
197 :
198 8 : TBranch * branch = hitTree->GetBranch("PMD") ;
199 4 : if ( ! branch )
200 : {
201 0 : AliWarning("PMD branch in Hit Tree not found") ;
202 0 : return;
203 : }
204 :
205 4 : if (fHitsArray)
206 3 : fHitsArray->Clear() ;
207 : else
208 2 : fHitsArray = new TClonesArray("AliPMDhit", 1000);
209 :
210 4 : branch->SetAddress(&fHitsArray);
211 :
212 232 : for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
213 112 : branch->GetEntry(ientry) ;
214 112 : MakeHits();
215 112 : fHitsArray->Clear() ;
216 : }
217 : //
218 4 : IncEvCountCycleHits();
219 4 : IncEvCountTotalHits();
220 : //
221 8 : }
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 8 : TIter next(fSDigitsArray) ;
231 : AliPMDsdigit * sdigit ;
232 44 : while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
233 : {
234 14 : if(sdigit->GetDetector() == 0)
235 : {
236 14 : edepkev = (sdigit->GetCellEdep())/1000.;
237 7 : FillSDigitsData(0,edepkev);
238 7 : premul++;
239 7 : }
240 14 : if(sdigit->GetDetector() == 1)
241 : {
242 0 : edepkev = (sdigit->GetCellEdep())/1000.;
243 0 : FillSDigitsData(1,edepkev);
244 0 : cpvmul++;
245 0 : }
246 :
247 : }
248 5 : if (premul > 0) FillSDigitsData(2,premul);
249 4 : if (cpvmul > 0) FillSDigitsData(3,cpvmul);
250 :
251 4 : }
252 :
253 : //____________________________________________________________________________
254 : void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
255 : {
256 : // makes data from SDigit Tree
257 :
258 8 : if (fSDigitsArray)
259 3 : fSDigitsArray->Clear() ;
260 : else
261 2 : fSDigitsArray = new TClonesArray("AliPMDsdigit", 1000) ;
262 :
263 4 : TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
264 4 : if ( ! branch )
265 : {
266 0 : AliWarning("PMD branch in SDigit Tree not found") ;
267 0 : }
268 : else
269 : {
270 4 : branch->SetAddress(&fSDigitsArray) ;
271 4 : branch->GetEntry(0) ;
272 4 : MakeSDigits() ;
273 : }
274 : //
275 4 : IncEvCountCycleSDigits();
276 4 : IncEvCountTotalSDigits();
277 : //
278 4 : }
279 :
280 : //____________________________________________________________________________
281 : void AliPMDQADataMakerSim::MakeDigits()
282 : {
283 : // makes data from Digits
284 :
285 : Int_t cpvmul = 0, premul = 0;
286 :
287 384 : TIter next(fDigitsArray) ;
288 : AliPMDdigit * digit ;
289 1488 : while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
290 : {
291 360 : if(digit->GetDetector() == 0)
292 : {
293 240 : FillDigitsData(0, digit->GetADC()) ;
294 120 : premul++;
295 120 : }
296 360 : if(digit->GetDetector() == 1)
297 : {
298 120 : FillDigitsData(1, digit->GetADC());
299 60 : cpvmul++;
300 60 : }
301 : }
302 :
303 215 : if (premul > 0) FillDigitsData(2,premul);
304 218 : if (cpvmul > 0) FillDigitsData(3,cpvmul);
305 :
306 :
307 192 : }
308 :
309 : //____________________________________________________________________________
310 : void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
311 : {
312 : // makes data from Digit Tree
313 :
314 8 : if (fDigitsArray)
315 3 : fDigitsArray->Clear() ;
316 : else
317 2 : fDigitsArray = new TClonesArray("AliPMDdigit", 1000) ;
318 :
319 4 : TBranch * branch = digitTree->GetBranch("PMDDigit") ;
320 :
321 4 : if ( ! branch )
322 : {
323 0 : AliWarning("PMD branch in Digit Tree not found") ;
324 0 : }
325 : else
326 : {
327 4 : branch->SetAddress(&fDigitsArray) ;
328 392 : for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
329 : {
330 192 : branch->GetEntry(ient) ;
331 192 : MakeDigits() ;
332 192 : fDigitsArray->Clear() ;
333 :
334 : }
335 :
336 : }
337 : //
338 4 : IncEvCountCycleDigits();
339 4 : IncEvCountTotalDigits();
340 : //
341 4 : }
342 :
343 :
344 : //____________________________________________________________________________
345 :
346 : void AliPMDQADataMakerSim::StartOfDetectorCycle()
347 : {
348 : //Detector specific actions at start of cycle
349 :
350 12 : }
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 18 : AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;
358 9 : }
|