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 **************************************************************************/
16 Based on the QA code for PHOS written by Yves Schutz July 2007
18 Authors: J.Klay (Cal Poly) May 2008
19 S. Salur LBL April 2008
23 // --- ROOT system ---
24 #include <TClonesArray.h>
31 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
36 #include "AliESDCaloCells.h"
37 #include "AliESDEvent.h"
39 #include "AliEMCALQADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliEMCALDigit.h"
42 #include "AliEMCALRecPoint.h"
43 #include "AliEMCALRawUtils.h"
44 #include "AliEMCALReconstructor.h"
45 #include "AliEMCALRecParam.h"
46 #include "AliRawReader.h"
47 #include "AliCaloRawStreamV3.h"
48 #include "AliEMCALGeoParams.h"
50 ClassImp(AliEMCALQADataMakerRec)
52 //____________________________________________________________________________
53 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec() :
54 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
55 fSuperModules(4) // FIXME!!! number of SuperModules; 4 for 2009; update default to 12 for later runs..
60 //____________________________________________________________________________
61 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
62 AliQADataMakerRec(), fSuperModules()
65 SetName((const char*)qadm.GetName()) ;
66 SetTitle((const char*)qadm.GetTitle());
67 fSuperModules = qadm.GetSuperModules();
70 //__________________________________________________________________
71 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
74 this->~AliEMCALQADataMakerRec();
75 new(this) AliEMCALQADataMakerRec(qadm);
79 //____________________________________________________________________________
80 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
82 //Detector specific actions at end of cycle
84 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
87 //____________________________________________________________________________
88 void AliEMCALQADataMakerRec::InitESDs()
90 //Create histograms to controll ESD
91 const Bool_t expert = kTRUE ;
92 const Bool_t image = kTRUE ;
94 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [MeV];Counts", 200, 0., 20.) ;
96 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
98 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
100 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
102 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [MeV];Counts", 500, 0., 250.) ;
104 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
106 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
108 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
112 //____________________________________________________________________________
113 void AliEMCALQADataMakerRec::InitDigits()
115 // create Digits histograms in Digits subdir
116 const Bool_t expert = kTRUE ;
117 const Bool_t image = kTRUE ;
119 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
121 Add2DigitsList(h0, 0, !expert, image) ;
122 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
124 Add2DigitsList(h1, 1, !expert, image) ;
127 //____________________________________________________________________________
128 void AliEMCALQADataMakerRec::InitRecPoints()
130 // create Reconstructed Points histograms in RecPoints subdir
131 const Bool_t expert = kTRUE ;
132 const Bool_t image = kTRUE ;
134 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [MeV];Counts",200, 0.,20.); //GeV
136 Add2RecPointsList(h0,kRecPE, !expert, image);
138 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
140 Add2RecPointsList(h1,kRecPM, !expert, image);
142 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
144 Add2RecPointsList(h2,kRecPDigM, !expert, image);
148 //____________________________________________________________________________
149 void AliEMCALQADataMakerRec::InitRaws()
151 // create Raws histograms in Raws subdir
152 const Bool_t expert = kTRUE ;
153 const Bool_t saveCorr = kTRUE ;
154 const Bool_t image = kTRUE ;
156 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
157 int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
159 // counter info: number of channels per event (bins are SM index)
160 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
161 fSuperModules, -0.5, fSuperModules-0.5) ;
162 Add2RawsList(h0, kNsmodLG, !expert, image, !saveCorr) ;
163 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
164 fSuperModules, -0.5, fSuperModules-0.5) ;
165 Add2RawsList(h1, kNsmodHG, !expert, image, !saveCorr) ;
167 // where did max sample occur? (bins are towers)
168 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
169 nTot, -0.5, nTot-0.5) ;
170 Add2RawsList(h2, kTimeLG, !expert, image, !saveCorr) ;
171 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
172 nTot, -0.5, nTot-0.5) ;
173 Add2RawsList(h3, kTimeHG, !expert, image, !saveCorr) ;
175 // how much above pedestal was the max sample? (bins are towers)
176 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
177 nTot, -0.5, nTot-0.5) ;
178 Add2RawsList(h4, kSigLG, !expert, image, !saveCorr) ;
179 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
180 nTot, -0.5, nTot-0.5) ;
181 Add2RawsList(h5, kSigHG, !expert, image, !saveCorr) ;
183 // total counter: channels per event
184 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
186 Add2RawsList(h6, kNtotLG, !expert, image, !saveCorr) ;
187 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
189 Add2RawsList(h7, kNtotHG, !expert, image, !saveCorr) ;
191 // pedestal (bins are towers)
192 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
193 nTot, -0.5, nTot-0.5) ;
194 Add2RawsList(h8, kPedLG, !expert, image, !saveCorr) ;
195 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
196 nTot, -0.5, nTot-0.5) ;
197 Add2RawsList(h9, kPedHG, !expert, image, !saveCorr) ;
199 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are towers)
200 TProfile * h10 = new TProfile("hLowEmcalRawPedRMS", "Low Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
201 nTot, -0.5, nTot-0.5) ;
202 Add2RawsList(h10, kPedRMSLG, !expert, image, !saveCorr) ;
203 TProfile * h11 = new TProfile("hHighEmcalRawPedRMS", "High Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
204 nTot, -0.5, nTot-0.5) ;
205 Add2RawsList(h11, kPedRMSHG, !expert, image, !saveCorr) ;
208 // now repeat the same for TRU and LEDMon data
209 int nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
211 // counter info: number of channels per event (bins are SM index)
212 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
213 fSuperModules, -0.5, fSuperModules-0.5) ;
214 Add2RawsList(hT0, kNsmodTRU, !expert, image, !saveCorr) ;
216 // where did max sample occur? (bins are TRU channels)
217 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
218 nTot2x2, -0.5, nTot2x2-0.5) ;
219 Add2RawsList(hT1, kTimeTRU, !expert, image, !saveCorr) ;
221 // how much above pedestal was the max sample? (bins are TRU channels)
222 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
223 nTot2x2, -0.5, nTot2x2-0.5) ;
224 Add2RawsList(hT2, kSigTRU, !expert, image, !saveCorr) ;
226 // total counter: channels per event
227 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
229 Add2RawsList(hT3, kNtotTRU, !expert, image, !saveCorr) ;
231 // pedestal (bins are TRU channels)
232 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
233 nTot2x2, -0.5, nTot2x2-0.5) ;
234 Add2RawsList(hT4, kPedTRU, !expert, image, !saveCorr) ;
236 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are TRU channels)
237 TProfile * hT5 = new TProfile("hTRUEmcalRawPedRMS", "TRU EMC: Pedestal RMS vs 2x2Id;2x2 Id;Width [ADC counts]",
238 nTot2x2, -0.5, nTot2x2-0.5) ;
239 Add2RawsList(hT5, kPedRMSTRU, !expert, image, !saveCorr) ;
241 // and also LED Mon..
242 // LEDMon has both high and low gain channels, just as regular FEE/towers
243 int nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
245 // counter info: number of channels per event (bins are SM index)
246 TProfile * hL0 = new TProfile("hLowLEDMONEmcalSupermodules", "LowLEDMON Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
247 fSuperModules, -0.5, fSuperModules-0.5) ;
248 Add2RawsList(hL0, kNsmodLGLEDMon, !expert, image, !saveCorr) ;
249 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
250 fSuperModules, -0.5, fSuperModules-0.5) ;
251 Add2RawsList(hL1, kNsmodHGLEDMon, !expert, image, !saveCorr) ;
253 // where did max sample occur? (bins are strips)
254 TProfile * hL2 = new TProfile("hLowLEDMONEmcalRawtime", "LowLEDMON Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
255 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
256 Add2RawsList(hL2, kTimeLGLEDMon, !expert, image, !saveCorr) ;
257 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
258 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
259 Add2RawsList(hL3, kTimeHGLEDMon, !expert, image, !saveCorr) ;
261 // how much above pedestal was the max sample? (bins are strips)
262 TProfile * hL4 = new TProfile("hLowLEDMONEmcalRawMaxMinusMin", "LowLEDMON Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
263 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
264 Add2RawsList(hL4, kSigLGLEDMon, !expert, image, !saveCorr) ;
265 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
266 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
267 Add2RawsList(hL5, kSigHGLEDMon, !expert, image, !saveCorr) ;
269 // total counter: channels per event
270 TH1I * hL6 = new TH1I("hLowLEDMONNtot", "LowLEDMON Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
272 Add2RawsList(hL6, kNtotLGLEDMon, !expert, image, !saveCorr) ;
273 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
275 Add2RawsList(hL7, kNtotHGLEDMon, !expert, image, !saveCorr) ;
277 // pedestal (bins are strips)
278 TProfile * hL8 = new TProfile("hLowLEDMONEmcalRawPed", "LowLEDMON Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
279 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
280 Add2RawsList(hL8, kPedLGLEDMon, !expert, image, !saveCorr) ;
281 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
282 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
283 Add2RawsList(hL9, kPedHGLEDMon, !expert, image, !saveCorr) ;
285 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are strips)
286 TProfile * hL10 = new TProfile("hLowLEDMONEmcalRawPedRMS", "LowLEDMON Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
287 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
288 Add2RawsList(hL10, kPedRMSLGLEDMon, !expert, image, !saveCorr) ;
289 TProfile * hL11 = new TProfile("hHighLEDMonEmcalRawPedRMS", "HighLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
290 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
291 Add2RawsList(hL11, kPedRMSHGLEDMon, !expert, image, !saveCorr) ;
295 //____________________________________________________________________________
296 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
298 // make QA data from ESDs
301 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
302 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
303 if( clu->IsEMCAL() ) {
304 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
308 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
311 AliESDCaloCells* cells = esd->GetEMCALCells();
312 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
314 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
315 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
320 //____________________________________________________________________________
321 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
323 //Fill prepared histograms with Raw digit properties
325 //Raw histogram filling not yet implemented
327 //Need to figure out how to get the info we want without having to
328 //actually run Raw2Digits twice.
329 //I suspect what we actually want is a raw digits method, not a true
330 //emcal raw data method, but this doesn't seem to be allowed in
331 //AliQADataMakerRec.h
333 // For now, to avoid redoing the expensive signal fits we just
334 // look at max vs min of the signal spextra, a la online usage in
335 // AliCaloCalibPedestal
338 AliCaloRawStreamV3 in(rawReader,"EMCAL");
341 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
342 int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
343 int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
344 int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
347 int sampleMax = 0x3ff; // 1023 = 10-bit range
349 // for the pedestal calculation
350 Bool_t selectPedestalSamples = kTRUE;
351 int firstPedestalSample = 0;
352 int lastPedestalSample = 15;
354 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
355 int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
356 int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
357 int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
358 int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
359 int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
361 // indices for the reading
365 // counters, on sample level
366 int i = 0; // the sample number in current event.
371 double meanPed = 0, squaredMean = 0, rmsPed = 0;
373 // start loop over input stream
374 while (in.NextDDL()) {
375 while (in.NextChannel()) {
378 int max = sampleMin, min = sampleMax; // min and max sample values
380 // for the pedestal calculation
381 int sampleSum = 0; // sum of samples
382 int squaredSampleSum = 0; // sum of samples squared
383 int nSum = 0; // number of samples in sum
385 while (in.NextBunch()) {
386 const UShort_t *sig = in.GetSignals();
387 startBin = in.GetStartTimeBin();
388 for (i = 0; i < in.GetBunchLength(); i++) {
392 // check if it's a min or max value
393 if (sample < min) min = sample;
399 // should we add it for the pedestal calculation?
400 if ( (firstPedestalSample<=time && time<=lastPedestalSample) || // sample time in range
401 !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all
403 squaredSampleSum += sample*sample;
407 } // loop over samples in bunch
408 } // loop over bunches
410 // calculate pedesstal estimate: mean of possibly selected samples
412 meanPed = sampleSum / (1.0 * nSum);
413 squaredMean = squaredSampleSum / (1.0 * nSum);
414 // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
415 rmsPed = sqrt(squaredMean - meanPed*meanPed);
423 // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
424 iSM = in.GetModule(); //The modules are numbered starting from 0
426 if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
428 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
429 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
431 if ( in.IsLowGain() ) {
432 //fill the low gain histograms, and counters
433 nTotalSMLG[iSM]++; // one more channel found
434 GetRawsData(kSigLG)->Fill(towerId, max - min);
435 GetRawsData(kTimeLG)->Fill(towerId, maxTime);
436 if (nSum>0) { // only fill pedestal info in case it could be calculated
437 GetRawsData(kPedLG)->Fill(towerId, meanPed);
438 GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
441 else if ( in.IsHighGain() ) {
442 //fill the high gain ones
443 nTotalSMHG[iSM]++; // one more channel found
444 GetRawsData(kSigHG)->Fill(towerId, max - min);
445 GetRawsData(kTimeHG)->Fill(towerId, maxTime);
446 if (nSum>0) { // only fill pedestal info in case it could be calculated
447 GetRawsData(kPedHG)->Fill(towerId, meanPed);
448 GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
451 } // low or high gain
453 else if ( in.IsTRUData() ) {
454 // for TRU data, the mapping class holds the channel info in the Column..
455 int TRU2x2Id = iSM*n2x2PerSM + in.GetColumn();
457 //fill the low gain histograms, and counters
458 nTotalSMTRU[iSM]++; // one more channel found
459 GetRawsData(kSigTRU)->Fill(TRU2x2Id, max - min);
460 GetRawsData(kTimeTRU)->Fill(TRU2x2Id, maxTime);
461 if (nSum>0) { // only fill pedestal info in case it could be calculated
462 GetRawsData(kPedTRU)->Fill(TRU2x2Id, meanPed);
463 GetRawsData(kPedRMSTRU)->Fill(TRU2x2Id, rmsPed);
467 else if ( in.IsLEDMonData() ) {
468 // for LED Mon data, the mapping class holds the gain info in the Row variable
469 // and the Strip number in the Column..
470 int gain = in.GetRow();
471 int stripId = iSM*nStripsPerSM + in.GetColumn();
474 //fill the low gain histograms, and counters
475 nTotalSMLGLEDMon[iSM]++; // one more channel found
476 GetRawsData(kSigLGLEDMon)->Fill(stripId, max - min);
477 GetRawsData(kTimeLGLEDMon)->Fill(stripId, maxTime);
478 if (nSum>0) { // only fill pedestal info in case it could be calculated
479 GetRawsData(kPedLGLEDMon)->Fill(stripId, meanPed);
480 GetRawsData(kPedRMSLGLEDMon)->Fill(stripId, rmsPed);
483 else if ( gain == 1 ) {
484 //fill the high gain ones
485 nTotalSMHGLEDMon[iSM]++; // one more channel found
486 GetRawsData(kSigHGLEDMon)->Fill(stripId, max - min);
487 GetRawsData(kTimeHGLEDMon)->Fill(stripId, maxTime);
488 if (nSum>0) { // only fill pedestal info in case it could be calculated
489 GetRawsData(kPedHGLEDMon)->Fill(stripId, meanPed);
490 GetRawsData(kPedRMSHGLEDMon)->Fill(stripId, rmsPed);
492 } // low or high gain
497 }// end while over channel
499 }//end while over DDL's, of input stream
501 // let's also fill the SM and event counter histograms
505 int nTotalHGLEDMon = 0;
506 int nTotalLGLEDMon = 0;
507 for (iSM=0; iSM<fSuperModules; iSM++) {
508 nTotalLG += nTotalSMLG[iSM];
509 nTotalHG += nTotalSMHG[iSM];
510 nTotalTRU += nTotalSMTRU[iSM];
511 nTotalLG += nTotalSMLGLEDMon[iSM];
512 nTotalHG += nTotalSMHGLEDMon[iSM];
513 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
514 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
515 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
516 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
517 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
519 GetRawsData(kNtotLG)->Fill(nTotalLG);
520 GetRawsData(kNtotHG)->Fill(nTotalHG);
521 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
522 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
523 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
525 // just in case the next rawreader consumer forgets to reset; let's do it here again..
531 //____________________________________________________________________________
532 void AliEMCALQADataMakerRec::MakeDigits()
534 // makes data from Digits
536 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
537 TIter next(fDigitsArray) ;
538 AliEMCALDigit * digit ;
539 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
540 GetDigitsData(0)->Fill( digit->GetAmp()) ;
545 //____________________________________________________________________________
546 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
548 // makes data from Digit Tree
550 fDigitsArray->Clear() ;
552 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
554 TBranch * branch = digitTree->GetBranch("EMCAL") ;
556 AliWarning("EMCAL branch in Digit Tree not found") ;
558 branch->SetAddress(&fDigitsArray) ;
559 branch->GetEntry(0) ;
565 //____________________________________________________________________________
566 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
568 // makes data from RecPoints
569 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
571 AliError("can't get the branch with the EMCAL clusters !");
575 TObjArray * emcrecpoints = new TObjArray(100) ;
576 emcbranch->SetAddress(&emcrecpoints);
577 emcbranch->GetEntry(0);
579 GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ;
580 TIter next(emcrecpoints) ;
581 AliEMCALRecPoint * rp ;
582 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
583 GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
584 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
586 emcrecpoints->Delete();
591 //____________________________________________________________________________
592 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
594 //Detector specific actions at start of cycle