]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerSim.cxx
resolving new library dependency libAOD
[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()
145 {
146     //make QA data from Hits
147
148   Int_t premul = 0, cpvmul = 0;
149   Float_t edepkev = 0.;
150   TIter next(fHitsArray); 
151   AliPMDhit * hit; 
152     
153   while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
154     {
155       if (hit->Z() > 361.5)
156         {
157           edepkev = (hit->GetEnergy())/1000.;
158           GetHitsData(0)->Fill(edepkev);
159           premul++;
160         }
161       else if (hit->Z() < 361.5)
162         {
163           edepkev = (hit->GetEnergy())/1000.;
164           GetHitsData(1)->Fill(edepkev);
165           cpvmul++;
166         }
167     }
168   
169   if(premul <= 0)
170     {
171       GetHitsData(2)->Fill(-1.); 
172     }
173   else
174     {
175       GetHitsData(2)->Fill(premul); 
176     }
177   
178   if(cpvmul <= 0)
179     {
180       GetHitsData(3)->Fill(-1.); 
181     }
182   else
183     {
184       GetHitsData(3)->Fill(cpvmul); 
185     }
186   
187 }
188
189 //____________________________________________________________________________
190 void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
191 {
192   // make QA data from Hit Tree
193   
194   TBranch * branch = hitTree->GetBranch("PMD") ;
195   if ( ! branch )
196     {
197       AliWarning("PMD branch in Hit Tree not found") ;
198       return;
199     }
200
201   if (fHitsArray) 
202     fHitsArray->Clear() ; 
203   else 
204     fHitsArray = new TClonesArray("AliPMDhit", 1000);
205
206   branch->SetAddress(&fHitsArray);
207
208   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
209     branch->GetEntry(ientry) ; 
210     MakeHits();
211     fHitsArray->Clear() ; 
212   }     
213 }
214 //____________________________________________________________________________
215 void AliPMDQADataMakerSim::MakeSDigits()
216 {
217     // makes data from SDigits
218
219   Int_t cpvmul = 0, premul = 0;
220   Float_t edepkev = 0.;
221   
222   TIter next(fSDigitsArray) ; 
223   AliPMDsdigit * sdigit ; 
224   while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
225     {
226       if(sdigit->GetDetector() == 0)
227         {
228           edepkev = (sdigit->GetCellEdep())/1000.;
229           GetSDigitsData(0)->Fill(edepkev);
230           premul++;
231         }
232       if(sdigit->GetDetector() == 1)
233         {
234           edepkev = (sdigit->GetCellEdep())/1000.;
235           GetSDigitsData(1)->Fill(edepkev);
236           cpvmul++;
237         }
238         
239     } 
240   if (premul > 0) GetSDigitsData(2)->Fill(premul);
241   if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
242   
243 }
244
245 //____________________________________________________________________________
246 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
247 {
248     // makes data from SDigit Tree
249   
250   if (fSDigitsArray) 
251     fSDigitsArray->Clear() ; 
252   else 
253     fSDigitsArray = new TClonesArray("AliPMDsdigit", 1000) ; 
254     
255   TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
256   if ( ! branch )
257     {
258       AliWarning("PMD branch in SDigit Tree not found") ; 
259     }
260   else
261     {
262       branch->SetAddress(&fSDigitsArray) ;
263       branch->GetEntry(0) ;
264       MakeSDigits() ; 
265     }
266 }
267
268 //____________________________________________________________________________
269 void AliPMDQADataMakerSim::MakeDigits()
270 {
271     // makes data from Digits
272
273   Int_t cpvmul = 0, premul = 0;
274   
275   TIter next(fDigitsArray) ; 
276   AliPMDdigit * digit ; 
277   while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
278     {
279       if(digit->GetDetector() == 0)
280         {
281           GetDigitsData(0)->Fill( digit->GetADC()) ;
282           premul++;
283         }
284       if(digit->GetDetector() == 1)
285         {
286           GetDigitsData(1)->Fill( digit->GetADC());
287           cpvmul++;
288         }
289     }  
290   
291   if (premul > 0) GetDigitsData(2)->Fill(premul);
292   if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
293   
294   
295 }
296
297 //____________________________________________________________________________
298 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
299 {
300     // makes data from Digit Tree
301
302   if (fDigitsArray) 
303     fDigitsArray->Clear() ; 
304   else 
305     fDigitsArray = new TClonesArray("AliPMDdigit", 1000) ; 
306   
307   TBranch * branch = digitTree->GetBranch("PMDDigit") ;
308
309   if ( ! branch )
310     {
311       AliWarning("PMD branch in Digit Tree not found") ; 
312     }
313   else
314     {
315       branch->SetAddress(&fDigitsArray) ;
316       for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
317         {
318           branch->GetEntry(ient) ; 
319           MakeDigits() ; 
320           fDigitsArray->Clear() ; 
321           
322         }
323       
324     }
325 }
326
327
328 //____________________________________________________________________________ 
329
330 void AliPMDQADataMakerSim::StartOfDetectorCycle()
331 {
332     //Detector specific actions at start of cycle
333     
334 }
335 //____________________________________________________________________________ 
336
337 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
338 {
339     //Detector specific actions at end of cycle
340     // do the QA checking
341     AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
342 }