]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQualAssDataMaker.cxx
Oversight for debugging purpose fixed
[u/mrichter/AliRoot.git] / PHOS / AliPHOSQualAssDataMaker.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 /* $Id$ */
18
19 /*
20   Produces the data needed to calculate the quality assurance. 
21   All data must be mergeable objects.
22   Y. Schutz CERN July 2007
23 */
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliESDCaloCluster.h"
35 #include "AliESDEvent.h"
36 #include "AliLog.h"
37 #include "AliPHOSDigit.h"
38 #include "AliPHOSHit.h"
39 #include "AliPHOSQualAssDataMaker.h"
40 #include "AliPHOSCpvRecPoint.h" 
41 #include "AliPHOSEmcRecPoint.h" 
42 #include "AliPHOSRecParticle.h" 
43 #include "AliPHOSTrackSegment.h" 
44
45 ClassImp(AliPHOSQualAssDataMaker)
46            
47 //____________________________________________________________________________ 
48   AliPHOSQualAssDataMaker::AliPHOSQualAssDataMaker() : 
49   AliQualAssDataMaker(AliQualAss::GetDetName(AliQualAss::kPHOS), "PHOS Quality Assurance Data Maker")
50 {
51   // ctor
52 }
53
54 //____________________________________________________________________________ 
55 AliPHOSQualAssDataMaker::AliPHOSQualAssDataMaker(const AliPHOSQualAssDataMaker& qadm) :
56   AliQualAssDataMaker()
57 {
58   //copy ctor 
59   SetName((const char*)qadm.GetName()) ; 
60   SetTitle((const char*)qadm.GetTitle()); 
61 }
62
63 //__________________________________________________________________
64 AliPHOSQualAssDataMaker& AliPHOSQualAssDataMaker::operator = (const AliPHOSQualAssDataMaker& qadm )
65 {
66   // Equal operator.
67   this->~AliPHOSQualAssDataMaker();
68   new(this) AliPHOSQualAssDataMaker(qadm);
69   return *this;
70 }
71  
72 //____________________________________________________________________________ 
73 void AliPHOSQualAssDataMaker::EndOfDetectorCycle()
74 {
75   //Detector specific actions at end of cycle
76 }
77
78 //____________________________________________________________________________ 
79 void AliPHOSQualAssDataMaker::InitESDs()
80 {
81   //create ESDs histograms in ESDs subdir
82   TH1F * h0 = new TH1F("hPhosESDs",    "ESDs energy distribution in PHOS",       100, 0., 100.) ;  
83   h0->Sumw2() ; 
84   Add2ESDsList(h0, 0) ;
85   TH1I * h1  = new TH1I("hPhosESDsMul", "ESDs multiplicity distribution in PHOS", 100, 0., 100) ; 
86   h1->Sumw2() ;
87   Add2ESDsList(h1, 1) ;
88 }
89
90 //____________________________________________________________________________ 
91 void AliPHOSQualAssDataMaker::InitHits()
92 {
93   // create Hits histograms in Hits subdir
94   TH1F * h0 = new TH1F("hPhosHits",    "Hits energy distribution in PHOS",       100, 0., 100.) ; 
95   h0->Sumw2() ;
96   Add2HitsList(h0, 0) ;
97   TH1I * h1  = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS", 500, 0., 10000) ; 
98   h1->Sumw2() ;
99   Add2HitsList(h1, 1) ;
100 }
101
102 //____________________________________________________________________________ 
103 void AliPHOSQualAssDataMaker::InitDigits()
104 {
105   // create Digits histograms in Digits subdir
106   TH1I * h0 = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS",    500, 0, 5000) ; 
107   h0->Sumw2() ;
108   Add2DigitsList(h0, 0) ;
109   TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 500, 0, 1000) ; 
110   h1->Sumw2() ;
111   Add2DigitsList(h1, 0) ;
112 }
113
114 //____________________________________________________________________________ 
115 //void AliPHOSQualAssDataMaker::InitRecParticles()
116 //{
117 //  // create Reconstructed particles histograms in RecParticles subdir
118 //  fhRecParticles     = new TH1F("hPhosRecParticles",    "RecParticles energy distribution in PHOS",       100, 0., 100.) ; 
119 //  fhRecParticles->Sumw2() ;
120 //  fhRecParticlesMul  = new TH1I("hPhosRecParticlesMul", "RecParticles multiplicity distribution in PHOS", 100, 0,  100) ; 
121 //  fhRecParticlesMul->Sumw2() ;
122 //}
123
124 //____________________________________________________________________________ 
125 void AliPHOSQualAssDataMaker::InitRecPoints()
126 {
127   // create Reconstructed Points histograms in RecPoints subdir
128   TH1F * h0 = new TH1F("hEmcPhosRecPoints",    "EMCA RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
129   h0->Sumw2() ;
130   Add2RecPointsList(h0, 0) ;
131   TH1I * h1 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
132   h1->Sumw2() ;
133   Add2RecPointsList(h1, 1) ;
134
135   TH1F * h2 = new TH1F("hCpvPhosRecPoints",    "CPV RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
136   h2->Sumw2() ;
137   Add2RecPointsList(h2, 2) ;
138   TH1I * h3 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
139   h3->Sumw2() ;
140   Add2RecPointsList(h3, 3) ;
141 }
142
143 //____________________________________________________________________________ 
144 void AliPHOSQualAssDataMaker::InitRaws()
145 {
146   // create Raws histograms in Raws subdir
147   TH1F * h0 = new TH1F("hEmcPhosRaws",    "EMCA Raws in PHOS",       100, 0., 100.) ; 
148   h0->Sumw2() ;
149   Add2RawsList(h0, 0) ;
150 }
151
152 //____________________________________________________________________________ 
153 void AliPHOSQualAssDataMaker::InitSDigits()
154 {
155   // create SDigits histograms in SDigits subdir
156   TH1F * h0 = new TH1F("hPhosSDigits",    "SDigits energy distribution in PHOS",       100, 0., 100.) ; 
157   h0->Sumw2() ;
158   Add2SDigitsList(h0, 0) ;
159   TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS", 500, 0,  10000) ; 
160   h1->Sumw2() ;
161   Add2SDigitsList(h1, 1) ;
162 }
163
164 //____________________________________________________________________________ 
165 //void AliPHOSQualAssDataMaker::InitTrackSegments()
166 //{
167 //  // create Track Segments histograms in TrackSegments subdir
168 //  fhTrackSegments     = new TH1F("hPhosTrackSegments",    "TrackSegments EMC-CPV distance in PHOS",       500, 0., 5000.) ; 
169 //  fhTrackSegments->Sumw2() ;
170 //  fhTrackSegmentsMul  = new TH1I("hPhosTrackSegmentsMul", "TrackSegments multiplicity distribution in PHOS", 100, 0,  100) ; 
171 //  fhTrackSegmentsMul->Sumw2() ;
172 //}
173
174 //____________________________________________________________________________
175 void AliPHOSQualAssDataMaker::MakeESDs(AliESDEvent * esd)
176 {
177   // make QA data from ESDs
178   
179   Int_t maxClu = esd->GetNumberOfPHOSClusters() ; 
180   Int_t index = 0, count = 0 ; 
181   for ( index = 0 ; index < maxClu; index++ ) {
182     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
183     GetESDsData(0)->Fill(clu->E()) ;
184     count++ ; 
185   }
186   GetESDsData(1)->Fill(count) ;
187 }
188
189 //____________________________________________________________________________
190 void AliPHOSQualAssDataMaker::MakeHits(TObject * data)
191 {
192   //make QA data from Hits
193
194   TClonesArray * hits = dynamic_cast<TClonesArray *>(data) ; 
195   if (!hits) {
196     AliError("Wrong type of hits container") ; 
197   } else {
198     GetHitsData(1)->Fill(hits->GetEntriesFast()) ; 
199     TIter next(hits) ; 
200     AliPHOSHit * hit ; 
201     while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
202       GetHitsData(0)->Fill( hit->GetEnergy()) ;
203     }
204   } 
205 }
206 //____________________________________________________________________________
207 void AliPHOSQualAssDataMaker::MakeDigits(TObject * data)
208 {
209   // makes data from Digits
210
211   TClonesArray * digits = dynamic_cast<TClonesArray *>(data) ; 
212   if (!digits) {
213     AliError("Wrong type of digits container") ; 
214   } else {
215     GetDigitsData(1)->Fill(digits->GetEntriesFast()) ; 
216     TIter next(digits) ; 
217     AliPHOSDigit * digit ; 
218     while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
219       GetDigitsData(0)->Fill( digit->GetEnergy()) ;
220     }  
221   }
222 }
223
224 //____________________________________________________________________________
225 // void AliPHOSQualAssDataMaker::MakeRecParticles(TTree * recpar)
226 // {
227 //   // makes data from RecParticles
228
229 //   TClonesArray * recparticles = dynamic_cast<TClonesArray*>(fData) ; 
230 //   fhRecParticlesMul->Fill(recparticles->GetEntriesFast()) ; 
231 //   TIter next(recparticles) ; 
232 //   AliPHOSRecParticle * recparticle ; 
233 //   while ( (recparticle = dynamic_cast<AliPHOSRecParticle *>(next())) ) {
234 //     fhRecParticles->Fill( recparticle->Energy()) ;
235 //   }
236 // }
237
238 //____________________________________________________________________________
239 void AliPHOSQualAssDataMaker::MakeRaws(TObject * data)
240 {
241     GetRawsData(1)->Fill(99) ; 
242 }
243
244 //____________________________________________________________________________
245 void AliPHOSQualAssDataMaker::MakeRecPoints(TTree * clustersTree)
246 {
247   {
248     // makes data from RecPoints
249     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
250     if (!emcbranch) { 
251       AliError("can't get the branch with the PHOS EMC clusters !");
252       return;
253     }
254     TObjArray * emcrecpoints = new TObjArray(100) ;
255     emcbranch->SetAddress(&emcrecpoints);
256     emcbranch->GetEntry(0);
257     
258     GetRecPointsData(1)->Fill(emcrecpoints->GetEntriesFast()) ; 
259     TIter next(emcrecpoints) ; 
260     AliPHOSEmcRecPoint * rp ; 
261     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
262       GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
263     }
264     emcrecpoints->Delete();
265     delete emcrecpoints;
266   }
267   {
268     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
269     if (!cpvbranch) { 
270       AliError("can't get the branch with the PHOS CPV clusters !");
271       return;
272     }
273     TObjArray *cpvrecpoints = new TObjArray(100) ;
274     cpvbranch->SetAddress(&cpvrecpoints);
275     cpvbranch->GetEntry(0);
276     
277     GetRecPointsData(1)->Fill(cpvrecpoints->GetEntriesFast()) ; 
278     TIter next(cpvrecpoints) ; 
279     AliPHOSCpvRecPoint * rp ; 
280     while ( (rp = dynamic_cast<AliPHOSCpvRecPoint *>(next())) ) {
281       GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
282     }
283     cpvrecpoints->Delete();
284     delete cpvrecpoints;
285   }
286 }
287
288 //____________________________________________________________________________
289 void AliPHOSQualAssDataMaker::MakeSDigits(TObject * data)
290 {
291   // makes data from SDigits
292   
293   TClonesArray * sdigits = dynamic_cast<TClonesArray *>(data) ; 
294   if (!sdigits) {
295     AliError("Wrong type of sdigits container") ; 
296   } else {
297     GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ; 
298     TIter next(sdigits) ; 
299     AliPHOSDigit * sdigit ; 
300     while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
301       GetSDigitsData(0)->Fill( sdigit->GetEnergy()) ;
302     } 
303   }
304 }
305
306 //____________________________________________________________________________
307 // void AliPHOSQualAssDataMaker::MakeTrackSegments(TTree * ts)
308 // {
309 //   // makes data from TrackSegments
310
311 //   TClonesArray * tracksegments = dynamic_cast<TClonesArray*>(fData) ;
312
313 //   fhTrackSegmentsMul->Fill(tracksegments->GetEntriesFast()) ; 
314 //   TIter next(tracksegments) ; 
315 //   AliPHOSTrackSegment * ts ; 
316 //   while ( (ts = dynamic_cast<AliPHOSTrackSegment *>(next())) ) {
317 //     fhTrackSegments->Fill( ts->GetCpvDistance()) ;
318 //   } 
319 // }
320
321 //____________________________________________________________________________ 
322 void AliPHOSQualAssDataMaker::StartOfDetectorCycle()
323 {
324   //Detector specific actions at start of cycle
325   
326 }