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
21 Created one histogram for QA shifter;
22 The idea:average counts for all the towers should be flat
23 Change all existing histograms as experts
28 // --- ROOT system ---
29 #include <TClonesArray.h>
36 // --- Standard library ---
39 // --- AliRoot header files ---
40 #include "AliESDCaloCluster.h"
41 #include "AliESDCaloCells.h"
42 #include "AliESDEvent.h"
44 #include "AliEMCALQADataMakerRec.h"
45 #include "AliQAChecker.h"
46 #include "AliEMCALDigit.h"
47 #include "AliEMCALRecPoint.h"
48 #include "AliEMCALRawUtils.h"
49 #include "AliEMCALReconstructor.h"
50 #include "AliEMCALRecParam.h"
51 #include "AliRawReader.h"
52 #include "AliCaloRawStreamV3.h"
53 #include "AliEMCALGeoParams.h"
55 ClassImp(AliEMCALQADataMakerRec)
57 //____________________________________________________________________________
58 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec() :
59 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
60 fSuperModules(4), // FIXME!!! number of SuperModules; 4 for 2009; update default to 12 for later runs..
61 fFirstPedestalSample(0),
62 fLastPedestalSample(15),
64 fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax)
69 //____________________________________________________________________________
70 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
72 fSuperModules(qadm.GetSuperModules()),
73 fFirstPedestalSample(qadm.GetFirstPedestalSample()),
74 fLastPedestalSample(qadm.GetLastPedestalSample()),
75 fMinSignalHG(qadm.GetMinSignalHG()),
76 fMaxSignalHG(qadm.GetMaxSignalHG())
79 SetName((const char*)qadm.GetName()) ;
80 SetTitle((const char*)qadm.GetTitle());
83 //__________________________________________________________________
84 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
87 this->~AliEMCALQADataMakerRec();
88 new(this) AliEMCALQADataMakerRec(qadm);
92 //____________________________________________________________________________
93 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
95 //Detector specific actions at end of cycle
98 // GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
100 // do the QA checking
101 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
104 //____________________________________________________________________________
105 void AliEMCALQADataMakerRec::InitESDs()
107 //Create histograms to controll ESD
108 const Bool_t expert = kTRUE ;
109 const Bool_t image = kTRUE ;
111 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts", 200, 0., 100.) ;
113 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
115 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
117 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
119 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts", 500, 0., 50.) ;
121 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
123 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
125 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
129 //____________________________________________________________________________
130 void AliEMCALQADataMakerRec::InitDigits()
132 // create Digits histograms in Digits subdir
133 const Bool_t expert = kTRUE ;
134 const Bool_t image = kTRUE ;
136 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
138 Add2DigitsList(h0, 0, !expert, image) ;
139 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
141 Add2DigitsList(h1, 1, !expert, image) ;
144 //____________________________________________________________________________
145 void AliEMCALQADataMakerRec::InitRecPoints()
147 // create Reconstructed Points histograms in RecPoints subdir
148 const Bool_t expert = kTRUE ;
149 const Bool_t image = kTRUE ;
151 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
153 Add2RecPointsList(h0,kRecPE, !expert, image);
155 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
157 Add2RecPointsList(h1,kRecPM, !expert, image);
159 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
161 Add2RecPointsList(h2,kRecPDigM, !expert, image);
165 //____________________________________________________________________________
166 void AliEMCALQADataMakerRec::InitRaws()
168 // create Raws histograms in Raws subdir
169 const Bool_t expert = kTRUE ;
170 const Bool_t saveCorr = kTRUE ;
171 const Bool_t image = kTRUE ;
173 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
174 int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
176 // counter info: number of channels per event (bins are SM index)
177 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
178 fSuperModules, -0.5, fSuperModules-0.5) ;
179 Add2RawsList(h0, kNsmodLG, expert, image, !saveCorr) ;
180 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
181 fSuperModules, -0.5, fSuperModules-0.5) ;
182 Add2RawsList(h1, kNsmodHG, expert, image, !saveCorr) ;
184 // where did max sample occur? (bins are towers)
185 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
186 nTot, -0.5, nTot-0.5) ;
187 Add2RawsList(h2, kTimeLG, expert, image, !saveCorr) ;
188 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
189 nTot, -0.5, nTot-0.5) ;
190 Add2RawsList(h3, kTimeHG, expert, image, !saveCorr) ;
192 // how much above pedestal was the max sample? (bins are towers)
193 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
194 nTot, -0.5, nTot-0.5) ;
195 Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
196 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
197 nTot, -0.5, nTot-0.5) ;
198 Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
200 // total counter: channels per event
201 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
203 Add2RawsList(h6, kNtotLG, expert, image, !saveCorr) ;
204 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
206 Add2RawsList(h7, kNtotHG, expert, image, !saveCorr) ;
208 // pedestal (bins are towers)
209 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
210 nTot, -0.5, nTot-0.5) ;
211 Add2RawsList(h8, kPedLG, expert, image, !saveCorr) ;
212 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
213 nTot, -0.5, nTot-0.5) ;
214 Add2RawsList(h9, kPedHG, expert, image, !saveCorr) ;
216 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are towers)
217 TProfile * h10 = new TProfile("hLowEmcalRawPedRMS", "Low Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
218 nTot, -0.5, nTot-0.5) ;
219 Add2RawsList(h10, kPedRMSLG, expert, image, !saveCorr) ;
220 TProfile * h11 = new TProfile("hHighEmcalRawPedRMS", "High Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
221 nTot, -0.5, nTot-0.5) ;
222 Add2RawsList(h11, kPedRMSHG, expert, image, !saveCorr) ;
224 //number of events per tower, for shifter fast check
225 TH1I * h12 = new TH1I("hTowerHG", "High Gains on the Tower;Tower", nTot,0, nTot) ;
227 Add2RawsList(h12, kTowerHG, !expert, image, !saveCorr) ;
228 TH1I * h13 = new TH1I("hTowerLG", "Low Gains on the Tower;Tower", nTot,0, nTot) ;
230 Add2RawsList(h13, kTowerLG, !expert, image, !saveCorr) ;
232 // now repeat the same for TRU and LEDMon data
233 int nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
235 // counter info: number of channels per event (bins are SM index)
236 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
237 fSuperModules, -0.5, fSuperModules-0.5) ;
238 Add2RawsList(hT0, kNsmodTRU, expert, image, !saveCorr) ;
240 // where did max sample occur? (bins are TRU channels)
241 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
242 nTot2x2, -0.5, nTot2x2-0.5) ;
243 Add2RawsList(hT1, kTimeTRU, expert, image, !saveCorr) ;
245 // how much above pedestal was the max sample? (bins are TRU channels)
246 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
247 nTot2x2, -0.5, nTot2x2-0.5) ;
248 Add2RawsList(hT2, kSigTRU, expert, image, !saveCorr) ;
250 // total counter: channels per event
251 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
253 Add2RawsList(hT3, kNtotTRU, expert, image, !saveCorr) ;
255 // pedestal (bins are TRU channels)
256 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
257 nTot2x2, -0.5, nTot2x2-0.5) ;
258 Add2RawsList(hT4, kPedTRU, expert, image, !saveCorr) ;
260 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are TRU channels)
261 TProfile * hT5 = new TProfile("hTRUEmcalRawPedRMS", "TRU EMC: Pedestal RMS vs 2x2Id;2x2 Id;Width [ADC counts]",
262 nTot2x2, -0.5, nTot2x2-0.5) ;
263 Add2RawsList(hT5, kPedRMSTRU, expert, image, !saveCorr) ;
265 // and also LED Mon..
266 // LEDMon has both high and low gain channels, just as regular FEE/towers
267 int nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
269 // counter info: number of channels per event (bins are SM index)
270 TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
271 fSuperModules, -0.5, fSuperModules-0.5) ;
272 Add2RawsList(hL0, kNsmodLGLEDMon, expert, image, !saveCorr) ;
273 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
274 fSuperModules, -0.5, fSuperModules-0.5) ;
275 Add2RawsList(hL1, kNsmodHGLEDMon, expert, image, !saveCorr) ;
277 // where did max sample occur? (bins are strips)
278 TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
279 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
280 Add2RawsList(hL2, kTimeLGLEDMon, expert, image, !saveCorr) ;
281 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
282 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
283 Add2RawsList(hL3, kTimeHGLEDMon, expert, image, !saveCorr) ;
285 // how much above pedestal was the max sample? (bins are strips)
286 TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
287 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
288 Add2RawsList(hL4, kSigLGLEDMon, expert, image, !saveCorr) ;
289 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
290 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
291 Add2RawsList(hL5, kSigHGLEDMon, expert, image, !saveCorr) ;
293 // total counter: channels per event
294 TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
296 Add2RawsList(hL6, kNtotLGLEDMon, expert, image, !saveCorr) ;
297 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
299 Add2RawsList(hL7, kNtotHGLEDMon, expert, image, !saveCorr) ;
301 // pedestal (bins are strips)
302 TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
303 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
304 Add2RawsList(hL8, kPedLGLEDMon, expert, image, !saveCorr) ;
305 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
306 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
307 Add2RawsList(hL9, kPedHGLEDMon, expert, image, !saveCorr) ;
309 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are strips)
310 TProfile * hL10 = new TProfile("hLowLEDMonEmcalRawPedRMS", "LowLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
311 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
312 Add2RawsList(hL10, kPedRMSLGLEDMon, expert, image, !saveCorr) ;
313 TProfile * hL11 = new TProfile("hHighLEDMonEmcalRawPedRMS", "HighLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
314 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
315 Add2RawsList(hL11, kPedRMSHGLEDMon, expert, image, !saveCorr) ;
319 //____________________________________________________________________________
320 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
322 // make QA data from ESDs
325 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
326 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
327 if( clu->IsEMCAL() ) {
328 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
332 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
335 AliESDCaloCells* cells = esd->GetEMCALCells();
336 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
338 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
339 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
344 //____________________________________________________________________________
345 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
347 //Fill prepared histograms with Raw digit properties
349 //Raw histogram filling not yet implemented
351 //Need to figure out how to get the info we want without having to
352 //actually run Raw2Digits twice.
353 //I suspect what we actually want is a raw digits method, not a true
354 //emcal raw data method, but this doesn't seem to be allowed in
355 //AliQADataMakerRec.h
357 // For now, to avoid redoing the expensive signal fits we just
358 // look at max vs min of the signal spextra, a la online usage in
359 // AliCaloCalibPedestal
362 AliCaloRawStreamV3 in(rawReader,"EMCAL");
365 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
366 int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
367 int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
368 int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
371 int sampleMax = AliEMCALGeoParams::fgkSampleMax; // 0x3ff = 1023 = 10-bit range
373 // for the pedestal calculation
374 Bool_t selectPedestalSamples = kTRUE;
376 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
377 int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
378 int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
379 int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
380 int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
381 int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
383 // indices for the reading
387 // counters, on sample level
388 int i = 0; // the sample number in current event.
393 double meanPed = 0, squaredMean = 0, rmsPed = 0;
395 // start loop over input stream
396 while (in.NextDDL()) {
397 int iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
398 while (in.NextChannel()) {
401 int max = sampleMin, min = sampleMax; // min and max sample values
404 // for the pedestal calculation
405 int sampleSum = 0; // sum of samples
406 int squaredSampleSum = 0; // sum of samples squared
407 int nSum = 0; // number of samples in sum
409 while (in.NextBunch()) {
410 const UShort_t *sig = in.GetSignals();
411 startBin = in.GetStartTimeBin();
412 nsamples += in.GetBunchLength();
413 for (i = 0; i < in.GetBunchLength(); i++) {
417 // check if it's a min or max value
418 if (sample < min) min = sample;
424 // should we add it for the pedestal calculation?
425 if ( (fFirstPedestalSample<=time && time<=fLastPedestalSample) || // sample time in range
426 !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all
428 squaredSampleSum += sample*sample;
432 } // loop over samples in bunch
433 } // loop over bunches
435 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
437 // calculate pedesstal estimate: mean of possibly selected samples
439 meanPed = sampleSum / (1.0 * nSum);
440 squaredMean = squaredSampleSum / (1.0 * nSum);
441 // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
442 rmsPed = sqrt(squaredMean - meanPed*meanPed);
450 // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
451 iSM = in.GetModule(); //The modules are numbered starting from 0
453 if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
455 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
456 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
459 if ( in.IsLowGain() ) {
460 //fill the low gain histograms, and counters
461 nTotalSMLG[iSM]++; // one more channel found
462 GetRawsData(kSigLG)->Fill(towerId, max - min);
463 GetRawsData(kTimeLG)->Fill(towerId, maxTime);
464 GetRawsData(kTowerLG)->Fill(towerId);
465 if (nSum>0) { // only fill pedestal info in case it could be calculated
466 GetRawsData(kPedLG)->Fill(towerId, meanPed);
467 GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
470 else if ( in.IsHighGain() ) {
471 //fill the high gain ones
472 nTotalSMHG[iSM]++; // one more channel found
473 int signal = max - min;
474 // only fill the max-min signal info and maxTime, if the
475 // signal was in the selected range
476 if ( (signal > fMinSignalHG) && (signal < fMaxSignalHG) ) {
477 GetRawsData(kSigHG)->Fill(towerId, signal);
478 GetRawsData(kTimeHG)->Fill(towerId, maxTime);
479 GetRawsData(kTowerHG)->Fill(towerId);
481 if (nSum>0) { // only fill pedestal info in case it could be calculated
482 GetRawsData(kPedHG)->Fill(towerId, meanPed);
483 GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
486 } // low or high gain
488 else if ( in.IsTRUData() ) {
489 // for TRU data, the mapping class holds the TRU internal 2x2 number (0..95) in the Column var..
490 int iTRU = iRCU; //TRU0 is from RCU0, TRU1 from RCU1
491 if (iRCU>0 && in.GetBranch()>0) iTRU=2; // TRU2 is from branch B on RCU1
492 int TRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU
495 //fill the low gain histograms, and counters
496 nTotalSMTRU[iSM]++; // one more channel found
497 GetRawsData(kSigTRU)->Fill(TRU2x2Id, max - min);
498 GetRawsData(kTimeTRU)->Fill(TRU2x2Id, maxTime);
499 if (nSum>0) { // only fill pedestal info in case it could be calculated
500 GetRawsData(kPedTRU)->Fill(TRU2x2Id, meanPed);
501 GetRawsData(kPedRMSTRU)->Fill(TRU2x2Id, rmsPed);
505 else if ( in.IsLEDMonData() ) {
506 // for LED Mon data, the mapping class holds the gain info in the Row variable
507 // and the Strip number in the Column..
508 int gain = in.GetRow();
509 int stripId = iSM*nStripsPerSM + in.GetColumn();
512 //fill the low gain histograms, and counters
513 nTotalSMLGLEDMon[iSM]++; // one more channel found
514 GetRawsData(kSigLGLEDMon)->Fill(stripId, max - min);
515 GetRawsData(kTimeLGLEDMon)->Fill(stripId, maxTime);
516 if (nSum>0) { // only fill pedestal info in case it could be calculated
517 GetRawsData(kPedLGLEDMon)->Fill(stripId, meanPed);
518 GetRawsData(kPedRMSLGLEDMon)->Fill(stripId, rmsPed);
521 else if ( gain == 1 ) {
522 //fill the high gain ones
523 nTotalSMHGLEDMon[iSM]++; // one more channel found
524 GetRawsData(kSigHGLEDMon)->Fill(stripId, max - min);
525 GetRawsData(kTimeHGLEDMon)->Fill(stripId, maxTime);
526 if (nSum>0) { // only fill pedestal info in case it could be calculated
527 GetRawsData(kPedHGLEDMon)->Fill(stripId, meanPed);
528 GetRawsData(kPedRMSHGLEDMon)->Fill(stripId, rmsPed);
530 } // low or high gain
535 } // nsamples>0 check, some data found for this channel; not only trailer/header
536 }// end while over channel
538 }//end while over DDL's, of input stream
540 // let's also fill the SM and event counter histograms
544 int nTotalHGLEDMon = 0;
545 int nTotalLGLEDMon = 0;
546 for (iSM=0; iSM<fSuperModules; iSM++) {
547 nTotalLG += nTotalSMLG[iSM];
548 nTotalHG += nTotalSMHG[iSM];
549 nTotalTRU += nTotalSMTRU[iSM];
550 nTotalLG += nTotalSMLGLEDMon[iSM];
551 nTotalHG += nTotalSMHGLEDMon[iSM];
552 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
553 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
554 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
555 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
556 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
558 GetRawsData(kNtotLG)->Fill(nTotalLG);
559 GetRawsData(kNtotHG)->Fill(nTotalHG);
560 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
561 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
562 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
564 // just in case the next rawreader consumer forgets to reset; let's do it here again..
570 //____________________________________________________________________________
571 void AliEMCALQADataMakerRec::MakeDigits()
573 // makes data from Digits
575 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
576 TIter next(fDigitsArray) ;
577 AliEMCALDigit * digit ;
578 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
579 GetDigitsData(0)->Fill( digit->GetAmp()) ;
584 //____________________________________________________________________________
585 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
587 // makes data from Digit Tree
589 fDigitsArray->Clear() ;
591 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
593 TBranch * branch = digitTree->GetBranch("EMCAL") ;
595 AliWarning("EMCAL branch in Digit Tree not found") ;
597 branch->SetAddress(&fDigitsArray) ;
598 branch->GetEntry(0) ;
604 //____________________________________________________________________________
605 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
607 // makes data from RecPoints
608 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
610 AliError("can't get the branch with the EMCAL clusters !");
614 TObjArray * emcrecpoints = new TObjArray(100) ;
615 emcbranch->SetAddress(&emcrecpoints);
616 emcbranch->GetEntry(0);
618 GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ;
619 TIter next(emcrecpoints) ;
620 AliEMCALRecPoint * rp ;
621 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
622 GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
623 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
625 emcrecpoints->Delete();
630 //____________________________________________________________________________
631 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
633 //Detector specific actions at start of cycle