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..
56 fFirstPedestalSample(0),
57 fLastPedestalSample(15)
62 //____________________________________________________________________________
63 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
65 fSuperModules(qadm.GetSuperModules()),
66 fFirstPedestalSample(qadm.GetFirstPedestalSample()),
67 fLastPedestalSample(qadm.GetLastPedestalSample())
70 SetName((const char*)qadm.GetName()) ;
71 SetTitle((const char*)qadm.GetTitle());
74 //__________________________________________________________________
75 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
78 this->~AliEMCALQADataMakerRec();
79 new(this) AliEMCALQADataMakerRec(qadm);
83 //____________________________________________________________________________
84 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
86 //Detector specific actions at end of cycle
88 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
91 //____________________________________________________________________________
92 void AliEMCALQADataMakerRec::InitESDs()
94 //Create histograms to controll ESD
95 const Bool_t expert = kTRUE ;
96 const Bool_t image = kTRUE ;
98 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [MeV];Counts", 200, 0., 20.) ;
100 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
102 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
104 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
106 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [MeV];Counts", 500, 0., 250.) ;
108 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
110 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
112 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
116 //____________________________________________________________________________
117 void AliEMCALQADataMakerRec::InitDigits()
119 // create Digits histograms in Digits subdir
120 const Bool_t expert = kTRUE ;
121 const Bool_t image = kTRUE ;
123 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
125 Add2DigitsList(h0, 0, !expert, image) ;
126 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
128 Add2DigitsList(h1, 1, !expert, image) ;
131 //____________________________________________________________________________
132 void AliEMCALQADataMakerRec::InitRecPoints()
134 // create Reconstructed Points histograms in RecPoints subdir
135 const Bool_t expert = kTRUE ;
136 const Bool_t image = kTRUE ;
138 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [MeV];Counts",200, 0.,20.); //GeV
140 Add2RecPointsList(h0,kRecPE, !expert, image);
142 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
144 Add2RecPointsList(h1,kRecPM, !expert, image);
146 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
148 Add2RecPointsList(h2,kRecPDigM, !expert, image);
152 //____________________________________________________________________________
153 void AliEMCALQADataMakerRec::InitRaws()
155 // create Raws histograms in Raws subdir
156 const Bool_t expert = kTRUE ;
157 const Bool_t saveCorr = kTRUE ;
158 const Bool_t image = kTRUE ;
160 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
161 int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
163 // counter info: number of channels per event (bins are SM index)
164 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
165 fSuperModules, -0.5, fSuperModules-0.5) ;
166 Add2RawsList(h0, kNsmodLG, !expert, image, !saveCorr) ;
167 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
168 fSuperModules, -0.5, fSuperModules-0.5) ;
169 Add2RawsList(h1, kNsmodHG, !expert, image, !saveCorr) ;
171 // where did max sample occur? (bins are towers)
172 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
173 nTot, -0.5, nTot-0.5) ;
174 Add2RawsList(h2, kTimeLG, !expert, image, !saveCorr) ;
175 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
176 nTot, -0.5, nTot-0.5) ;
177 Add2RawsList(h3, kTimeHG, !expert, image, !saveCorr) ;
179 // how much above pedestal was the max sample? (bins are towers)
180 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
181 nTot, -0.5, nTot-0.5) ;
182 Add2RawsList(h4, kSigLG, !expert, image, !saveCorr) ;
183 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
184 nTot, -0.5, nTot-0.5) ;
185 Add2RawsList(h5, kSigHG, !expert, image, !saveCorr) ;
187 // total counter: channels per event
188 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
190 Add2RawsList(h6, kNtotLG, !expert, image, !saveCorr) ;
191 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
193 Add2RawsList(h7, kNtotHG, !expert, image, !saveCorr) ;
195 // pedestal (bins are towers)
196 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
197 nTot, -0.5, nTot-0.5) ;
198 Add2RawsList(h8, kPedLG, !expert, image, !saveCorr) ;
199 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
200 nTot, -0.5, nTot-0.5) ;
201 Add2RawsList(h9, kPedHG, !expert, image, !saveCorr) ;
203 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are towers)
204 TProfile * h10 = new TProfile("hLowEmcalRawPedRMS", "Low Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
205 nTot, -0.5, nTot-0.5) ;
206 Add2RawsList(h10, kPedRMSLG, !expert, image, !saveCorr) ;
207 TProfile * h11 = new TProfile("hHighEmcalRawPedRMS", "High Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
208 nTot, -0.5, nTot-0.5) ;
209 Add2RawsList(h11, kPedRMSHG, !expert, image, !saveCorr) ;
212 // now repeat the same for TRU and LEDMon data
213 int nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
215 // counter info: number of channels per event (bins are SM index)
216 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
217 fSuperModules, -0.5, fSuperModules-0.5) ;
218 Add2RawsList(hT0, kNsmodTRU, !expert, image, !saveCorr) ;
220 // where did max sample occur? (bins are TRU channels)
221 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
222 nTot2x2, -0.5, nTot2x2-0.5) ;
223 Add2RawsList(hT1, kTimeTRU, !expert, image, !saveCorr) ;
225 // how much above pedestal was the max sample? (bins are TRU channels)
226 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
227 nTot2x2, -0.5, nTot2x2-0.5) ;
228 Add2RawsList(hT2, kSigTRU, !expert, image, !saveCorr) ;
230 // total counter: channels per event
231 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
233 Add2RawsList(hT3, kNtotTRU, !expert, image, !saveCorr) ;
235 // pedestal (bins are TRU channels)
236 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
237 nTot2x2, -0.5, nTot2x2-0.5) ;
238 Add2RawsList(hT4, kPedTRU, !expert, image, !saveCorr) ;
240 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are TRU channels)
241 TProfile * hT5 = new TProfile("hTRUEmcalRawPedRMS", "TRU EMC: Pedestal RMS vs 2x2Id;2x2 Id;Width [ADC counts]",
242 nTot2x2, -0.5, nTot2x2-0.5) ;
243 Add2RawsList(hT5, kPedRMSTRU, !expert, image, !saveCorr) ;
245 // and also LED Mon..
246 // LEDMon has both high and low gain channels, just as regular FEE/towers
247 int nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
249 // counter info: number of channels per event (bins are SM index)
250 TProfile * hL0 = new TProfile("hLowLEDMONEmcalSupermodules", "LowLEDMON Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
251 fSuperModules, -0.5, fSuperModules-0.5) ;
252 Add2RawsList(hL0, kNsmodLGLEDMon, !expert, image, !saveCorr) ;
253 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
254 fSuperModules, -0.5, fSuperModules-0.5) ;
255 Add2RawsList(hL1, kNsmodHGLEDMon, !expert, image, !saveCorr) ;
257 // where did max sample occur? (bins are strips)
258 TProfile * hL2 = new TProfile("hLowLEDMONEmcalRawtime", "LowLEDMON Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
259 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
260 Add2RawsList(hL2, kTimeLGLEDMon, !expert, image, !saveCorr) ;
261 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
262 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
263 Add2RawsList(hL3, kTimeHGLEDMon, !expert, image, !saveCorr) ;
265 // how much above pedestal was the max sample? (bins are strips)
266 TProfile * hL4 = new TProfile("hLowLEDMONEmcalRawMaxMinusMin", "LowLEDMON Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
267 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
268 Add2RawsList(hL4, kSigLGLEDMon, !expert, image, !saveCorr) ;
269 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
270 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
271 Add2RawsList(hL5, kSigHGLEDMon, !expert, image, !saveCorr) ;
273 // total counter: channels per event
274 TH1I * hL6 = new TH1I("hLowLEDMONNtot", "LowLEDMON Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
276 Add2RawsList(hL6, kNtotLGLEDMon, !expert, image, !saveCorr) ;
277 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
279 Add2RawsList(hL7, kNtotHGLEDMon, !expert, image, !saveCorr) ;
281 // pedestal (bins are strips)
282 TProfile * hL8 = new TProfile("hLowLEDMONEmcalRawPed", "LowLEDMON Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
283 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
284 Add2RawsList(hL8, kPedLGLEDMon, !expert, image, !saveCorr) ;
285 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
286 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
287 Add2RawsList(hL9, kPedHGLEDMon, !expert, image, !saveCorr) ;
289 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are strips)
290 TProfile * hL10 = new TProfile("hLowLEDMONEmcalRawPedRMS", "LowLEDMON Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
291 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
292 Add2RawsList(hL10, kPedRMSLGLEDMon, !expert, image, !saveCorr) ;
293 TProfile * hL11 = new TProfile("hHighLEDMonEmcalRawPedRMS", "HighLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
294 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
295 Add2RawsList(hL11, kPedRMSHGLEDMon, !expert, image, !saveCorr) ;
299 //____________________________________________________________________________
300 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
302 // make QA data from ESDs
305 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
306 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
307 if( clu->IsEMCAL() ) {
308 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
312 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
315 AliESDCaloCells* cells = esd->GetEMCALCells();
316 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
318 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
319 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
324 //____________________________________________________________________________
325 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
327 //Fill prepared histograms with Raw digit properties
329 //Raw histogram filling not yet implemented
331 //Need to figure out how to get the info we want without having to
332 //actually run Raw2Digits twice.
333 //I suspect what we actually want is a raw digits method, not a true
334 //emcal raw data method, but this doesn't seem to be allowed in
335 //AliQADataMakerRec.h
337 // For now, to avoid redoing the expensive signal fits we just
338 // look at max vs min of the signal spextra, a la online usage in
339 // AliCaloCalibPedestal
342 AliCaloRawStreamV3 in(rawReader,"EMCAL");
345 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
346 int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
347 int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
348 int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
351 int sampleMax = 0x3ff; // 1023 = 10-bit range
353 // for the pedestal calculation
354 Bool_t selectPedestalSamples = kTRUE;
356 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
357 int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
358 int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
359 int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
360 int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
361 int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
363 // indices for the reading
367 // counters, on sample level
368 int i = 0; // the sample number in current event.
373 double meanPed = 0, squaredMean = 0, rmsPed = 0;
375 // start loop over input stream
376 while (in.NextDDL()) {
377 while (in.NextChannel()) {
380 int max = sampleMin, min = sampleMax; // min and max sample values
382 // for the pedestal calculation
383 int sampleSum = 0; // sum of samples
384 int squaredSampleSum = 0; // sum of samples squared
385 int nSum = 0; // number of samples in sum
387 while (in.NextBunch()) {
388 const UShort_t *sig = in.GetSignals();
389 startBin = in.GetStartTimeBin();
390 for (i = 0; i < in.GetBunchLength(); i++) {
394 // check if it's a min or max value
395 if (sample < min) min = sample;
401 // should we add it for the pedestal calculation?
402 if ( (fFirstPedestalSample<=time && time<=fLastPedestalSample) || // sample time in range
403 !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all
405 squaredSampleSum += sample*sample;
409 } // loop over samples in bunch
410 } // loop over bunches
412 // calculate pedesstal estimate: mean of possibly selected samples
414 meanPed = sampleSum / (1.0 * nSum);
415 squaredMean = squaredSampleSum / (1.0 * nSum);
416 // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
417 rmsPed = sqrt(squaredMean - meanPed*meanPed);
425 // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
426 iSM = in.GetModule(); //The modules are numbered starting from 0
428 if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
430 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
431 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
433 if ( in.IsLowGain() ) {
434 //fill the low gain histograms, and counters
435 nTotalSMLG[iSM]++; // one more channel found
436 GetRawsData(kSigLG)->Fill(towerId, max - min);
437 GetRawsData(kTimeLG)->Fill(towerId, maxTime);
438 if (nSum>0) { // only fill pedestal info in case it could be calculated
439 GetRawsData(kPedLG)->Fill(towerId, meanPed);
440 GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
443 else if ( in.IsHighGain() ) {
444 //fill the high gain ones
445 nTotalSMHG[iSM]++; // one more channel found
446 GetRawsData(kSigHG)->Fill(towerId, max - min);
447 GetRawsData(kTimeHG)->Fill(towerId, maxTime);
448 if (nSum>0) { // only fill pedestal info in case it could be calculated
449 GetRawsData(kPedHG)->Fill(towerId, meanPed);
450 GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
453 } // low or high gain
455 else if ( in.IsTRUData() ) {
456 // for TRU data, the mapping class holds the channel info in the Column..
457 int TRU2x2Id = iSM*n2x2PerSM + in.GetColumn();
459 //fill the low gain histograms, and counters
460 nTotalSMTRU[iSM]++; // one more channel found
461 GetRawsData(kSigTRU)->Fill(TRU2x2Id, max - min);
462 GetRawsData(kTimeTRU)->Fill(TRU2x2Id, maxTime);
463 if (nSum>0) { // only fill pedestal info in case it could be calculated
464 GetRawsData(kPedTRU)->Fill(TRU2x2Id, meanPed);
465 GetRawsData(kPedRMSTRU)->Fill(TRU2x2Id, rmsPed);
469 else if ( in.IsLEDMonData() ) {
470 // for LED Mon data, the mapping class holds the gain info in the Row variable
471 // and the Strip number in the Column..
472 int gain = in.GetRow();
473 int stripId = iSM*nStripsPerSM + in.GetColumn();
476 //fill the low gain histograms, and counters
477 nTotalSMLGLEDMon[iSM]++; // one more channel found
478 GetRawsData(kSigLGLEDMon)->Fill(stripId, max - min);
479 GetRawsData(kTimeLGLEDMon)->Fill(stripId, maxTime);
480 if (nSum>0) { // only fill pedestal info in case it could be calculated
481 GetRawsData(kPedLGLEDMon)->Fill(stripId, meanPed);
482 GetRawsData(kPedRMSLGLEDMon)->Fill(stripId, rmsPed);
485 else if ( gain == 1 ) {
486 //fill the high gain ones
487 nTotalSMHGLEDMon[iSM]++; // one more channel found
488 GetRawsData(kSigHGLEDMon)->Fill(stripId, max - min);
489 GetRawsData(kTimeHGLEDMon)->Fill(stripId, maxTime);
490 if (nSum>0) { // only fill pedestal info in case it could be calculated
491 GetRawsData(kPedHGLEDMon)->Fill(stripId, meanPed);
492 GetRawsData(kPedRMSHGLEDMon)->Fill(stripId, rmsPed);
494 } // low or high gain
499 }// end while over channel
501 }//end while over DDL's, of input stream
503 // let's also fill the SM and event counter histograms
507 int nTotalHGLEDMon = 0;
508 int nTotalLGLEDMon = 0;
509 for (iSM=0; iSM<fSuperModules; iSM++) {
510 nTotalLG += nTotalSMLG[iSM];
511 nTotalHG += nTotalSMHG[iSM];
512 nTotalTRU += nTotalSMTRU[iSM];
513 nTotalLG += nTotalSMLGLEDMon[iSM];
514 nTotalHG += nTotalSMHGLEDMon[iSM];
515 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
516 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
517 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
518 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
519 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
521 GetRawsData(kNtotLG)->Fill(nTotalLG);
522 GetRawsData(kNtotHG)->Fill(nTotalHG);
523 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
524 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
525 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
527 // just in case the next rawreader consumer forgets to reset; let's do it here again..
533 //____________________________________________________________________________
534 void AliEMCALQADataMakerRec::MakeDigits()
536 // makes data from Digits
538 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
539 TIter next(fDigitsArray) ;
540 AliEMCALDigit * digit ;
541 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
542 GetDigitsData(0)->Fill( digit->GetAmp()) ;
547 //____________________________________________________________________________
548 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
550 // makes data from Digit Tree
552 fDigitsArray->Clear() ;
554 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
556 TBranch * branch = digitTree->GetBranch("EMCAL") ;
558 AliWarning("EMCAL branch in Digit Tree not found") ;
560 branch->SetAddress(&fDigitsArray) ;
561 branch->GetEntry(0) ;
567 //____________________________________________________________________________
568 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
570 // makes data from RecPoints
571 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
573 AliError("can't get the branch with the EMCAL clusters !");
577 TObjArray * emcrecpoints = new TObjArray(100) ;
578 emcbranch->SetAddress(&emcrecpoints);
579 emcbranch->GetEntry(0);
581 GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ;
582 TIter next(emcrecpoints) ;
583 AliEMCALRecPoint * rp ;
584 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
585 GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
586 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
588 emcrecpoints->Delete();
593 //____________________________________________________________________________
594 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
596 //Detector specific actions at start of cycle