]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerSim.cxx
c19afe643ff26a23ffa89e96e758ffca61264ae5
[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 in (keV)PRE(PMD)", 500, 0., 500.); 
77   h0->Sumw2() ;
78   Add2HitsList(h0, 0, !expert, image) ;
79   
80   TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution in (keV)CPV(PMD)", 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)", 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)", 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)", 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)", 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)", 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)", 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)", 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)", 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)", 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)", 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     Int_t premul = 0, cpvmul = 0;
149     Float_t edepkev = 0.;
150     TIter next(hits); 
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     static TClonesArray statichits("AliPMDhit", 1000);
202     statichits.Clear();
203     TClonesArray *hits = &statichits;
204     static TClonesArray staticdummy("AliPMDhit", 1000);
205     staticdummy.Clear();
206     TClonesArray *dummy = &staticdummy;
207     branch->SetAddress(&dummy);
208     Int_t index = 0 ;  
209
210     for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
211         branch->GetEntry(ientry) ; 
212         for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
213             AliPMDhit * hit = dynamic_cast<AliPMDhit *> (dummy->At(ihit)) ; 
214             new((*hits)[index]) AliPMDhit(*hit) ; 
215
216             index++ ;
217         } 
218     }   
219
220     MakeHits(hits);
221
222 }
223 //____________________________________________________________________________
224 void AliPMDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
225 {
226     // makes data from SDigits
227
228     Int_t cpvmul = 0, premul = 0;
229     Float_t edepkev = 0.;
230
231     TIter next(sdigits) ; 
232     AliPMDsdigit * sdigit ; 
233     while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
234     {
235         if(sdigit->GetDetector() == 0)
236         {
237           edepkev = (sdigit->GetCellEdep())/1000.;
238           GetSDigitsData(0)->Fill(edepkev);
239           premul++;
240         }
241         if(sdigit->GetDetector() == 1)
242         {
243           edepkev = (sdigit->GetCellEdep())/1000.;
244           GetSDigitsData(1)->Fill(edepkev);
245           cpvmul++;
246         }
247         
248     } 
249     if (premul > 0) GetSDigitsData(2)->Fill(premul);
250     if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
251     
252 }
253
254 //____________________________________________________________________________
255 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
256 {
257     // makes data from SDigit Tree
258
259     TClonesArray * sdigits = new TClonesArray("AliPMDsdigit", 1000) ; 
260     
261     TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
262     branch->SetAddress(&sdigits) ;
263
264     if ( ! branch )
265     {
266         AliWarning("PMD branch in SDigit Tree not found") ; 
267     }
268     else
269     {
270         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
271         {
272             branch->GetEntry(ient) ;
273             MakeSDigits(sdigits) ; 
274         }
275     }
276 }
277
278 //____________________________________________________________________________
279 void AliPMDQADataMakerSim::MakeDigits(TClonesArray * digits)
280 {
281     // makes data from Digits
282     
283     Int_t cpvmul = 0, premul = 0;
284
285     TIter next(digits) ; 
286     AliPMDdigit * digit ; 
287     while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
288     {
289         if(digit->GetDetector() == 0)
290         {
291             GetDigitsData(0)->Fill( digit->GetADC()) ;
292             premul++;
293         }
294         if(digit->GetDetector() == 1)
295         {
296             GetDigitsData(1)->Fill( digit->GetADC());
297             cpvmul++;
298         }
299     }  
300
301     if (premul > 0) GetDigitsData(2)->Fill(premul);
302     if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
303
304
305 }
306
307 //____________________________________________________________________________
308 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
309 {
310     // makes data from Digit Tree
311
312     TClonesArray * digits = new TClonesArray("AliPMDdigit", 1000) ; 
313     
314     TBranch * branch = digitTree->GetBranch("PMDDigit") ;
315     branch->SetAddress(&digits) ;
316
317     if ( ! branch )
318     {
319         AliWarning("PMD branch in Digit Tree not found") ; 
320     }
321     else
322     {
323         for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
324         {
325             branch->GetEntry(ient) ; 
326             MakeDigits(digits) ; 
327         }
328         
329     }
330 }
331
332
333 //____________________________________________________________________________ 
334
335 void AliPMDQADataMakerSim::StartOfDetectorCycle()
336 {
337     //Detector specific actions at start of cycle
338     
339 }
340 //____________________________________________________________________________ 
341
342 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
343 {
344     //Detector specific actions at end of cycle
345     // do the QA checking
346     AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
347 }