1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 Produces the data needed to calculate the quality assurance.
21 All data must be mergeable objects.
22 Y. Schutz CERN July 2007
25 // --- ROOT system ---
26 #include <TClonesArray.h>
32 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
36 #include "AliESDEvent.h"
38 #include "AliPHOSDigit.h"
39 #include "AliPHOSHit.h"
40 #include "AliPHOSQADataMaker.h"
41 #include "AliQAChecker.h"
42 #include "AliPHOSCpvRecPoint.h"
43 #include "AliPHOSEmcRecPoint.h"
44 #include "AliPHOSRecParticle.h"
45 #include "AliPHOSTrackSegment.h"
46 #include "AliPHOSRawDecoder.h"
48 ClassImp(AliPHOSQADataMaker)
50 //____________________________________________________________________________
51 AliPHOSQADataMaker::AliPHOSQADataMaker() :
52 AliQADataMaker(AliQA::GetDetName(AliQA::kPHOS), "PHOS Quality Assurance Data Maker")
57 //____________________________________________________________________________
58 AliPHOSQADataMaker::AliPHOSQADataMaker(const AliPHOSQADataMaker& qadm) :
62 SetName((const char*)qadm.GetName()) ;
63 SetTitle((const char*)qadm.GetTitle());
66 //__________________________________________________________________
67 AliPHOSQADataMaker& AliPHOSQADataMaker::operator = (const AliPHOSQADataMaker& qadm )
70 this->~AliPHOSQADataMaker();
71 new(this) AliPHOSQADataMaker(qadm);
75 //____________________________________________________________________________
76 void AliPHOSQADataMaker::EndOfDetectorCycle(AliQA::TASKINDEX task, TList * list)
78 //Detector specific actions at end of cycle
80 AliQAChecker::Instance()->Run(AliQA::kPHOS, task, list) ;
83 //____________________________________________________________________________
84 void AliPHOSQADataMaker::InitESDs()
86 //create ESDs histograms in ESDs subdir
87 TH1F * h0 = new TH1F("hPhosESDs", "ESDs energy distribution in PHOS", 100, 0., 100.) ;
90 TH1I * h1 = new TH1I("hPhosESDsMul", "ESDs multiplicity distribution in PHOS", 100, 0., 100) ;
95 //____________________________________________________________________________
96 void AliPHOSQADataMaker::InitHits()
98 // create Hits histograms in Hits subdir
99 TH1F * h0 = new TH1F("hPhosHits", "Hits energy distribution in PHOS", 100, 0., 100.) ;
101 Add2HitsList(h0, 0) ;
102 TH1I * h1 = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS", 500, 0., 10000) ;
104 Add2HitsList(h1, 1) ;
107 //____________________________________________________________________________
108 void AliPHOSQADataMaker::InitDigits()
110 // create Digits histograms in Digits subdir
111 TH1I * h0 = new TH1I("hPhosDigits", "Digits amplitude distribution in PHOS", 500, 0, 5000) ;
113 Add2DigitsList(h0, 0) ;
114 TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 500, 0, 1000) ;
116 Add2DigitsList(h1, 0) ;
119 //____________________________________________________________________________
120 //void AliPHOSQADataMaker::InitRecParticles()
122 // // create Reconstructed particles histograms in RecParticles subdir
123 // fhRecParticles = new TH1F("hPhosRecParticles", "RecParticles energy distribution in PHOS", 100, 0., 100.) ;
124 // fhRecParticles->Sumw2() ;
125 // fhRecParticlesMul = new TH1I("hPhosRecParticlesMul", "RecParticles multiplicity distribution in PHOS", 100, 0, 100) ;
126 // fhRecParticlesMul->Sumw2() ;
129 //____________________________________________________________________________
130 void AliPHOSQADataMaker::InitRecPoints()
132 // create Reconstructed Points histograms in RecPoints subdir
133 TH1F * h0 = new TH1F("hEmcPhosRecPoints", "EMCA RecPoints energy distribution in PHOS", 100, 0., 100.) ;
135 Add2RecPointsList(h0, 0) ;
136 TH1I * h1 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0, 100) ;
138 Add2RecPointsList(h1, 1) ;
140 TH1F * h2 = new TH1F("hCpvPhosRecPoints", "CPV RecPoints energy distribution in PHOS", 100, 0., 100.) ;
142 Add2RecPointsList(h2, 2) ;
143 TH1I * h3 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0, 100) ;
145 Add2RecPointsList(h3, 3) ;
148 //____________________________________________________________________________
149 void AliPHOSQADataMaker::InitRaws()
151 // create Raws histograms in Raws subdir
152 const Int_t modMax = 5 ;
153 TH2I * h0[modMax*2] ;
156 for (Int_t mod = 0; mod < modMax; mod++) {
157 sprintf(title, "Low Gain Rows x Columns for PHOS module %d", mod) ;
158 sprintf(name, "hLowPHOSxyMod%d", mod) ;
159 h0[mod] = new TH2I(name, title, 64, 1, 65, 56, 1, 57) ;
160 Add2RawsList(h0[mod], mod) ;
161 sprintf(title, "High Gain Rows x Columns for PHOS module %d", mod) ;
162 sprintf(name, "hHighPHOSxyMod%d", mod) ;
163 h0[mod+modMax] = new TH2I(name, title, 64, 1, 65, 56, 1, 57) ;
164 Add2RawsList(h0[mod+modMax], mod+modMax) ;
166 TH1I * h10 = new TH1I("hLowPhosModules", "Low Gain Hits in EMCA PHOS modules", 6, 0, 6) ;
168 Add2RawsList(h10, 10) ;
169 TH1I * h11 = new TH1I("hHighPhosModules", "High Gain Hits in EMCA PHOS modules", 6, 0, 6) ;
171 Add2RawsList(h11, 11) ;
172 TH1F * h12 = new TH1F("hLowPhosRawtime", "Low Gain Time of raw hits in PHOS", 100, 0, 100.) ;
174 Add2RawsList(h12, 12) ;
175 TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 100, 0, 100.) ;
177 Add2RawsList(h13, 13) ;
178 TH1F * h14 = new TH1F("hLowPhosRawEnergy", "Low Gain Energy of raw hits in PHOS", 100, 0., 100.) ;
180 Add2RawsList(h14, 14) ;
181 TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS", 100, 0., 100.) ;
183 Add2RawsList(h15, 15) ;
187 //____________________________________________________________________________
188 void AliPHOSQADataMaker::InitSDigits()
190 // create SDigits histograms in SDigits subdir
191 TH1F * h0 = new TH1F("hPhosSDigits", "SDigits energy distribution in PHOS", 100, 0., 100.) ;
193 Add2SDigitsList(h0, 0) ;
194 TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS", 500, 0, 10000) ;
196 Add2SDigitsList(h1, 1) ;
199 //____________________________________________________________________________
200 //void AliPHOSQADataMaker::InitTrackSegments()
202 // // create Track Segments histograms in TrackSegments subdir
203 // fhTrackSegments = new TH1F("hPhosTrackSegments", "TrackSegments EMC-CPV distance in PHOS", 500, 0., 5000.) ;
204 // fhTrackSegments->Sumw2() ;
205 // fhTrackSegmentsMul = new TH1I("hPhosTrackSegmentsMul", "TrackSegments multiplicity distribution in PHOS", 100, 0, 100) ;
206 // fhTrackSegmentsMul->Sumw2() ;
209 //____________________________________________________________________________
210 void AliPHOSQADataMaker::MakeESDs(AliESDEvent * esd)
212 // make QA data from ESDs
214 Int_t maxClu = esd->GetNumberOfPHOSClusters() ;
215 Int_t index = 0, count = 0 ;
216 for ( index = 0 ; index < maxClu; index++ ) {
217 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
218 GetESDsData(0)->Fill(clu->E()) ;
221 GetESDsData(1)->Fill(count) ;
224 //____________________________________________________________________________
225 void AliPHOSQADataMaker::MakeHits(TClonesArray * hits)
227 //make QA data from Hits
229 GetHitsData(1)->Fill(hits->GetEntriesFast()) ;
232 while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
233 GetHitsData(0)->Fill( hit->GetEnergy()) ;
236 //____________________________________________________________________________
237 void AliPHOSQADataMaker::MakeDigits(TClonesArray * digits)
239 // makes data from Digits
241 GetDigitsData(1)->Fill(digits->GetEntriesFast()) ;
243 AliPHOSDigit * digit ;
244 while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
245 GetDigitsData(0)->Fill( digit->GetEnergy()) ;
249 //____________________________________________________________________________
250 // void AliPHOSQADataMaker::MakeRecParticles(TTree * recpar)
252 // // makes data from RecParticles
254 // TClonesArray * recparticles = dynamic_cast<TClonesArray*>(fData) ;
255 // fhRecParticlesMul->Fill(recparticles->GetEntriesFast()) ;
256 // TIter next(recparticles) ;
257 // AliPHOSRecParticle * recparticle ;
258 // while ( (recparticle = dynamic_cast<AliPHOSRecParticle *>(next())) ) {
259 // fhRecParticles->Fill( recparticle->Energy()) ;
263 //____________________________________________________________________________
264 void AliPHOSQADataMaker::MakeRaws(AliRawReader* rawReader)
266 const Int_t modMax = 5 ;
268 AliPHOSRawDecoder decoder(rawReader);
269 decoder.SetOldRCUFormat(kTRUE);
270 decoder.SubtractPedestals(kTRUE);
273 while (decoder.NextDigit()) {
274 Int_t module = decoder.GetModule() ;
275 Int_t row = decoder.GetRow() ;
276 Int_t col = decoder.GetColumn() ;
277 Double_t time = decoder.GetTime() ;
278 Double_t energy = decoder.GetEnergy() ;
279 Bool_t lowGain = decoder.IsLowGain();
281 GetRawsData(module)->Fill(row, col) ;
282 GetRawsData(10)->Fill(module) ;
283 GetRawsData(12)->Fill(time) ;
284 GetRawsData(14)->Fill(energy) ;
286 GetRawsData(module+modMax)->Fill(row, col) ;
287 GetRawsData(11)->Fill(module) ;
288 GetRawsData(13)->Fill(time) ;
289 GetRawsData(15)->Fill(energy) ;
291 //AliInfo(Form(" %d %d %d %d %f %f\n", count, module, row, col, time, energy)) ;
296 //____________________________________________________________________________
297 void AliPHOSQADataMaker::MakeRecPoints(TTree * clustersTree)
300 // makes data from RecPoints
301 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
303 AliError("can't get the branch with the PHOS EMC clusters !");
306 TObjArray * emcrecpoints = new TObjArray(100) ;
307 emcbranch->SetAddress(&emcrecpoints);
308 emcbranch->GetEntry(0);
310 GetRecPointsData(1)->Fill(emcrecpoints->GetEntriesFast()) ;
311 TIter next(emcrecpoints) ;
312 AliPHOSEmcRecPoint * rp ;
313 while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
314 GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
316 emcrecpoints->Delete();
320 TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
322 AliError("can't get the branch with the PHOS CPV clusters !");
325 TObjArray *cpvrecpoints = new TObjArray(100) ;
326 cpvbranch->SetAddress(&cpvrecpoints);
327 cpvbranch->GetEntry(0);
329 GetRecPointsData(1)->Fill(cpvrecpoints->GetEntriesFast()) ;
330 TIter next(cpvrecpoints) ;
331 AliPHOSCpvRecPoint * rp ;
332 while ( (rp = dynamic_cast<AliPHOSCpvRecPoint *>(next())) ) {
333 GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
335 cpvrecpoints->Delete();
340 //____________________________________________________________________________
341 void AliPHOSQADataMaker::MakeSDigits(TClonesArray * sdigits)
343 // makes data from SDigits
345 GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ;
346 TIter next(sdigits) ;
347 AliPHOSDigit * sdigit ;
348 while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
349 GetSDigitsData(0)->Fill( sdigit->GetEnergy()) ;
353 //____________________________________________________________________________
354 // void AliPHOSQADataMaker::MakeTrackSegments(TTree * ts)
356 // // makes data from TrackSegments
358 // TClonesArray * tracksegments = dynamic_cast<TClonesArray*>(fData) ;
360 // fhTrackSegmentsMul->Fill(tracksegments->GetEntriesFast()) ;
361 // TIter next(tracksegments) ;
362 // AliPHOSTrackSegment * ts ;
363 // while ( (ts = dynamic_cast<AliPHOSTrackSegment *>(next())) ) {
364 // fhTrackSegments->Fill( ts->GetCpvDistance()) ;
368 //____________________________________________________________________________
369 void AliPHOSQADataMaker::StartOfDetectorCycle()
371 //Detector specific actions at start of cycle