]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQualAssDataMaker.cxx
767d8fe31e755850a193ce48bad99b022feffef6
[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   fhHits(0x0), 
51   fhHitsMul(0x0), 
52   fhDigits(0x0),
53   fhDigitsMul(0x0),
54   fhSDigits(0x0),
55   fhSDigitsMul(0x0),
56   fhEmcRecPoints(0x0),
57   fhEmcRecPointsMul(0x0),
58   fhCpvRecPoints(0x0),
59   fhCpvRecPointsMul(0x0),
60   fhTrackSegments(0x0),
61   fhTrackSegmentsMul(0x0),
62   fhRecParticles(0x0),
63   fhRecParticlesMul(0x0),
64   fhESDs(0x0),
65   fhESDsMul(0x0) 
66 {
67   // ctor
68   fDetectorDir = fOutput->GetDirectory(GetName()) ;  
69   if (!fDetectorDir) 
70     fDetectorDir = fOutput->mkdir(GetName()) ;  
71 }
72
73 //____________________________________________________________________________ 
74 AliPHOSQualAssDataMaker::AliPHOSQualAssDataMaker(const AliPHOSQualAssDataMaker& qadm) :
75   AliQualAssDataMaker(), 
76   fhHits(qadm.fhHits), 
77   fhHitsMul(qadm.fhHitsMul), 
78   fhDigits(qadm.fhDigits),
79   fhDigitsMul(qadm.fhDigitsMul),
80   fhSDigits(qadm.fhSDigits),
81   fhSDigitsMul(qadm.fhSDigitsMul), 
82   fhEmcRecPoints(qadm.fhEmcRecPoints),
83   fhEmcRecPointsMul(qadm.fhEmcRecPointsMul), 
84   fhCpvRecPoints(qadm.fhCpvRecPoints),
85   fhCpvRecPointsMul(qadm.fhCpvRecPointsMul), 
86   fhTrackSegments(qadm.fhTrackSegments),
87   fhTrackSegmentsMul(qadm.fhTrackSegmentsMul), 
88   fhRecParticles(qadm.fhRecParticles),
89   fhRecParticlesMul(qadm.fhRecParticlesMul), 
90   fhESDs(qadm.fhESDs), 
91   fhESDsMul(qadm.fhESDsMul) 
92 {
93   //copy ctor 
94   SetName((const char*)qadm.GetName()) ; 
95   SetTitle((const char*)qadm.GetTitle()); 
96 }
97
98 //__________________________________________________________________
99 AliPHOSQualAssDataMaker& AliPHOSQualAssDataMaker::operator = (const AliPHOSQualAssDataMaker& qadm )
100 {
101   // Equal operator.
102   this->~AliPHOSQualAssDataMaker();
103   new(this) AliPHOSQualAssDataMaker(qadm);
104   return *this;
105 }
106  
107 //____________________________________________________________________________ 
108 void AliPHOSQualAssDataMaker::InitESDs()
109 {
110   //create ESDs histograms in ESDs subdir
111   fhESDs     = new TH1F("hPhosESDs",    "ESDs energy distribution in PHOS",       100, 0., 100.) ; 
112   fhESDs->Sumw2() ; 
113   fhESDsMul  = new TH1I("hPhosESDsMul", "ESDs multiplicity distribution in PHOS", 100, 0., 100) ; 
114   fhESDsMul->Sumw2() ;
115 }
116
117 //____________________________________________________________________________ 
118 void AliPHOSQualAssDataMaker::InitHits()
119 {
120   // create Hits histograms in Hits subdir
121   fhHits     = new TH1F("hPhosHits",    "Hits energy distribution in PHOS",       100, 0., 100.) ; 
122   fhHits->Sumw2() ;
123   fhHitsMul  = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS", 500, 0., 10000) ; 
124   fhHitsMul->Sumw2() ;
125 }
126
127 //____________________________________________________________________________ 
128 void AliPHOSQualAssDataMaker::InitDigits()
129 {
130   // create Digits histograms in Digits subdir
131   fhDigits     = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS",    500, 0, 5000) ; 
132   fhDigits->Sumw2() ;
133   fhDigitsMul  = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 500, 0, 1000) ; 
134   fhDigitsMul->Sumw2() ;
135 }
136
137 //____________________________________________________________________________ 
138 void AliPHOSQualAssDataMaker::InitRecParticles()
139 {
140   // create Reconstructed particles histograms in RecParticles subdir
141   fhRecParticles     = new TH1F("hPhosRecParticles",    "RecParticles energy distribution in PHOS",       100, 0., 100.) ; 
142   fhRecParticles->Sumw2() ;
143   fhRecParticlesMul  = new TH1I("hPhosRecParticlesMul", "RecParticles multiplicity distribution in PHOS", 100, 0,  100) ; 
144   fhRecParticlesMul->Sumw2() ;
145 }
146
147 //____________________________________________________________________________ 
148 void AliPHOSQualAssDataMaker::InitRecPoints()
149 {
150   // create Reconstructed Points histograms in RecPoints subdir
151   fhEmcRecPoints     = new TH1F("hEmcPhosRecPoints",    "EMCA RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
152   fhEmcRecPoints->Sumw2() ;
153   fhEmcRecPointsMul  = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
154   fhEmcRecPointsMul->Sumw2() ;
155
156   fhCpvRecPoints     = new TH1F("hCpvPhosRecPoints",    "CPV RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
157   fhCpvRecPoints->Sumw2() ;
158   fhCpvRecPointsMul  = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
159   fhCpvRecPointsMul->Sumw2() ;
160 }
161
162 //____________________________________________________________________________ 
163 void AliPHOSQualAssDataMaker::InitSDigits()
164 {
165   // create SDigits histograms in SDigits subdir
166   fhSDigits     = new TH1F("hPhosSDigits",    "SDigits energy distribution in PHOS",       100, 0., 100.) ; 
167   fhSDigits->Sumw2() ;
168   fhSDigitsMul  = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS", 500, 0,  10000) ; 
169   fhSDigitsMul->Sumw2() ;
170 }
171
172 //____________________________________________________________________________ 
173 void AliPHOSQualAssDataMaker::InitTrackSegments()
174 {
175   // create Track Segments histograms in TrackSegments subdir
176   fhTrackSegments     = new TH1F("hPhosTrackSegments",    "TrackSegments EMC-CPV distance in PHOS",       500, 0., 5000.) ; 
177   fhTrackSegments->Sumw2() ;
178   fhTrackSegmentsMul  = new TH1I("hPhosTrackSegmentsMul", "TrackSegments multiplicity distribution in PHOS", 100, 0,  100) ; 
179   fhTrackSegmentsMul->Sumw2() ;
180 }
181
182 //____________________________________________________________________________
183 void AliPHOSQualAssDataMaker::MakeESDs()
184 {
185   // make QA data from ESDs
186   
187   AliESDEvent * esd = dynamic_cast<AliESDEvent*>(fData) ; 
188   Int_t maxClu = esd->GetNumberOfPHOSClusters() ; 
189   Int_t index = 0, count = 0 ; 
190   for ( index = 0 ; index < maxClu; index++ ) {
191     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
192     fhESDs->Fill(clu->E()) ;
193     count++ ; 
194   }
195   fhESDsMul->Fill(count) ;
196 }
197
198 //____________________________________________________________________________
199 void AliPHOSQualAssDataMaker::MakeHits()
200 {
201   //make QA data from Hits
202
203   TClonesArray * hits = dynamic_cast<TClonesArray*>(fData) ; 
204   fhHitsMul->Fill(hits->GetEntriesFast()) ; 
205   TIter next(hits) ; 
206   AliPHOSHit * hit ; 
207   while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
208     fhHits->Fill( hit->GetEnergy()) ;
209   }
210
211  
212 //____________________________________________________________________________
213 void AliPHOSQualAssDataMaker::MakeDigits()
214 {
215   // makes data from Digits
216
217   TClonesArray * digits = dynamic_cast<TClonesArray*>(fData) ; 
218   fhDigitsMul->Fill(digits->GetEntriesFast()) ; 
219   TIter next(digits) ; 
220   AliPHOSDigit * digit ; 
221   while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
222     fhDigits->Fill( digit->GetEnergy()) ;
223   }  
224 }
225
226 //____________________________________________________________________________
227 void AliPHOSQualAssDataMaker::MakeRecParticles()
228 {
229   // makes data from RecParticles
230
231   TClonesArray * recparticles = dynamic_cast<TClonesArray*>(fData) ; 
232   fhRecParticlesMul->Fill(recparticles->GetEntriesFast()) ; 
233   TIter next(recparticles) ; 
234   AliPHOSRecParticle * recparticle ; 
235   while ( (recparticle = dynamic_cast<AliPHOSRecParticle *>(next())) ) {
236     fhRecParticles->Fill( recparticle->Energy()) ;
237   }
238 }
239
240 //____________________________________________________________________________
241 void AliPHOSQualAssDataMaker::MakeRecPoints()
242 {
243   // makes data from RecPoints
244   TObjArray * recpoints = dynamic_cast<TObjArray*>(fData) ;  
245   TIter next(recpoints) ; 
246   
247   if ( strcmp(fData->GetName(), "EMCRECPOINTS") == 0 ) {
248     fhEmcRecPointsMul->Fill(recpoints->GetEntriesFast()) ; 
249     AliPHOSEmcRecPoint * rp ; 
250     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
251       fhEmcRecPoints->Fill( rp->GetEnergy()) ;
252     }
253   } 
254   else if  ( strcmp(fData->GetName(), "CPVRECPOINTS") == 0 ) {
255     fhCpvRecPointsMul->Fill(recpoints->GetEntriesFast()) ; 
256     AliPHOSCpvRecPoint * rp ; 
257     while ( (rp = dynamic_cast<AliPHOSCpvRecPoint *>(next())) ) {
258       fhCpvRecPoints->Fill( rp->GetEnergy()) ;
259     }
260   }  
261 }
262
263 //____________________________________________________________________________
264 void AliPHOSQualAssDataMaker::MakeSDigits()
265 {
266   // makes data from SDigits
267
268   TClonesArray * sdigits = dynamic_cast<TClonesArray*>(fData) ; 
269   fhSDigitsMul->Fill(sdigits->GetEntriesFast()) ; 
270   TIter next(sdigits) ; 
271   AliPHOSDigit * sdigit ; 
272   while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
273     fhSDigits->Fill( sdigit->GetEnergy()) ;
274   } 
275 }
276
277 //____________________________________________________________________________
278 void AliPHOSQualAssDataMaker::MakeTrackSegments()
279 {
280   // makes data from TrackSegments
281
282   TClonesArray * tracksegments = dynamic_cast<TClonesArray*>(fData) ;
283
284   fhTrackSegmentsMul->Fill(tracksegments->GetEntriesFast()) ; 
285   TIter next(tracksegments) ; 
286   AliPHOSTrackSegment * ts ; 
287   while ( (ts = dynamic_cast<AliPHOSTrackSegment *>(next())) ) {
288     fhTrackSegments->Fill( ts->GetCpvDistance()) ;
289   } 
290 }