]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALQADataMakerRec.cxx
updating histogram names (Theo)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMakerRec.cxx
CommitLineData
9e47432c 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/*
16Based on the QA code for PHOS written by Yves Schutz July 2007
17
18Authors: J.Klay (Cal Poly) May 2008
19 S. Salur LBL April 2008
20
21Created one histogram for QA shifter;
22The idea:average counts for all the towers should be flat
23Change all existing histograms as experts
24 --By Yaxian Mao
25
26*/
27
28// --- ROOT system ---
29#include <TClonesArray.h>
30#include <TFile.h>
31#include <TH1F.h>
32#include <TH1I.h>
33#include <TH2F.h>
34#include <TProfile.h>
35
36// --- Standard library ---
37
38
39// --- AliRoot header files ---
40#include "AliESDCaloCluster.h"
41#include "AliESDCaloCells.h"
42#include "AliESDEvent.h"
43#include "AliLog.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"
54
55ClassImp(AliEMCALQADataMakerRec)
56
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),
63 fMinSignalHG(0),
64 fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax)
65 {
66 // ctor
67}
68
69//____________________________________________________________________________
70AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
71 AliQADataMakerRec(),
72 fSuperModules(qadm.GetSuperModules()),
73 fFirstPedestalSample(qadm.GetFirstPedestalSample()),
74 fLastPedestalSample(qadm.GetLastPedestalSample()),
75 fMinSignalHG(qadm.GetMinSignalHG()),
76 fMaxSignalHG(qadm.GetMaxSignalHG())
77{
78 //copy ctor
79 SetName((const char*)qadm.GetName()) ;
80 SetTitle((const char*)qadm.GetTitle());
81}
82
83//__________________________________________________________________
84AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
85{
86 // Equal operator.
87 this->~AliEMCALQADataMakerRec();
88 new(this) AliEMCALQADataMakerRec(qadm);
89 return *this;
90}
91
92//____________________________________________________________________________
93void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
94{
95 //Detector specific actions at end of cycle
96
97// if(fCycleCounter)
98// GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
99
100 // do the QA checking
101 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
102}
103
104//____________________________________________________________________________
105void AliEMCALQADataMakerRec::InitESDs()
106{
107 //Create histograms to controll ESD
108 const Bool_t expert = kTRUE ;
109 const Bool_t image = kTRUE ;
110
111 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts", 200, 0., 100.) ;
112 h1->Sumw2() ;
113 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
114
115 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
116 h2->Sumw2() ;
117 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
118
119 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts", 500, 0., 50.) ;
120 h3->Sumw2() ;
121 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
122
123 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
124 h4->Sumw2() ;
125 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
126
127}
128
129//____________________________________________________________________________
130void AliEMCALQADataMakerRec::InitDigits()
131{
132 // create Digits histograms in Digits subdir
133 const Bool_t expert = kTRUE ;
134 const Bool_t image = kTRUE ;
135
136 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
137 h0->Sumw2() ;
138 Add2DigitsList(h0, 0, !expert, image) ;
139 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
140 h1->Sumw2() ;
141 Add2DigitsList(h1, 1, !expert, image) ;
142}
143
144//____________________________________________________________________________
145void AliEMCALQADataMakerRec::InitRecPoints()
146{
147 // create Reconstructed Points histograms in RecPoints subdir
148 const Bool_t expert = kTRUE ;
149 const Bool_t image = kTRUE ;
150
151 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
152 h0->Sumw2();
153 Add2RecPointsList(h0,kRecPE, !expert, image);
154
155 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
156 h1->Sumw2();
157 Add2RecPointsList(h1,kRecPM, !expert, image);
158
159 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
160 h2->Sumw2();
161 Add2RecPointsList(h2,kRecPDigM, !expert, image);
162
163}
164
165//____________________________________________________________________________
166void AliEMCALQADataMakerRec::InitRaws()
167{
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 ;
172
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
175
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) ;
183
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) ;
191
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) ;
199
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) ;
202 h6->Sumw2() ;
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) ;
205 h7->Sumw2() ;
206 Add2RawsList(h7, kNtotHG, expert, image, !saveCorr) ;
207
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) ;
215
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) ;
223
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) ;
226 h12->Sumw2() ;
227 Add2RawsList(h12, kTowerHG, !expert, image, !saveCorr) ;
228 TH1I * h13 = new TH1I("hTowerLG", "Low Gains on the Tower;Tower", nTot,0, nTot) ;
229 h13->Sumw2() ;
230 Add2RawsList(h13, kTowerLG, !expert, image, !saveCorr) ;
231
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
234
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) ;
239
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) ;
244
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) ;
249
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) ;
252 hT3->Sumw2() ;
253 Add2RawsList(hT3, kNtotTRU, expert, image, !saveCorr) ;
254
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) ;
259
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) ;
264
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
268
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) ;
276
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) ;
284
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) ;
292
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) ;
295 hL6->Sumw2() ;
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) ;
298 hL7->Sumw2() ;
299 Add2RawsList(hL7, kNtotHGLEDMon, expert, image, !saveCorr) ;
300
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) ;
308
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) ;
316
317}
318
319//____________________________________________________________________________
320void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
321{
322 // make QA data from ESDs
323
324 Int_t nTot = 0 ;
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()) ;
329 nTot++ ;
330 }
331 }
332 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
333
334 //fill calo cells
335 AliESDCaloCells* cells = esd->GetEMCALCells();
336 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
337
338 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
339 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
340 }
341
342}
343
344//____________________________________________________________________________
345void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
346{
347 //Fill prepared histograms with Raw digit properties
348
349 //Raw histogram filling not yet implemented
350 //
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
356
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
360
361 rawReader->Reset() ;
362 AliCaloRawStreamV3 in(rawReader,"EMCAL");
30aa89b0 363 rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's
9e47432c 364
365 // setup
366 int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
367 int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
368 int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
369 int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
370
371 int sampleMin = 0;
372 int sampleMax = AliEMCALGeoParams::fgkSampleMax; // 0x3ff = 1023 = 10-bit range
373
374 // for the pedestal calculation
375 Bool_t selectPedestalSamples = kTRUE;
376
377 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
378 int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
379 int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
380 int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
381 int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
382 int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
383
384 // indices for the reading
385 int iSM = 0;
386 int sample = 0;
387 int time = 0;
388 // counters, on sample level
389 int i = 0; // the sample number in current event.
390 int maxTime = 0;
391 int startBin = 0;
392
393 // calc. quantities
394 double meanPed = 0, squaredMean = 0, rmsPed = 0;
395
396 // start loop over input stream
397 while (in.NextDDL()) {
398 int iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
399 while (in.NextChannel()) {
400
401 // counters
402 int max = sampleMin, min = sampleMax; // min and max sample values
403 int nsamples = 0;
404
405 // for the pedestal calculation
406 int sampleSum = 0; // sum of samples
407 int squaredSampleSum = 0; // sum of samples squared
408 int nSum = 0; // number of samples in sum
409
410 while (in.NextBunch()) {
411 const UShort_t *sig = in.GetSignals();
412 startBin = in.GetStartTimeBin();
413 nsamples += in.GetBunchLength();
414 for (i = 0; i < in.GetBunchLength(); i++) {
415 sample = sig[i];
416 time = startBin--;
417
418 // check if it's a min or max value
419 if (sample < min) min = sample;
420 if (sample > max) {
421 max = sample;
422 maxTime = time;
423 }
424
425 // should we add it for the pedestal calculation?
426 if ( (fFirstPedestalSample<=time && time<=fLastPedestalSample) || // sample time in range
427 !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all
428 sampleSum += sample;
429 squaredSampleSum += sample*sample;
430 nSum++;
431 }
432
433 } // loop over samples in bunch
434 } // loop over bunches
435
436 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
437
438 // calculate pedesstal estimate: mean of possibly selected samples
439 if (nSum > 0) {
440 meanPed = sampleSum / (1.0 * nSum);
441 squaredMean = squaredSampleSum / (1.0 * nSum);
442 // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
443 rmsPed = sqrt(squaredMean - meanPed*meanPed);
444 }
445 else {
446 meanPed = 0;
447 squaredMean = 0;
448 rmsPed = 0;
449 }
450
451 // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
452 iSM = in.GetModule(); //The modules are numbered starting from 0
453
454 if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
455
456 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
457 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
458
459
460 if ( in.IsLowGain() ) {
461 //fill the low gain histograms, and counters
462 nTotalSMLG[iSM]++; // one more channel found
463 GetRawsData(kSigLG)->Fill(towerId, max - min);
464 GetRawsData(kTimeLG)->Fill(towerId, maxTime);
465 GetRawsData(kTowerLG)->Fill(towerId);
466 if (nSum>0) { // only fill pedestal info in case it could be calculated
467 GetRawsData(kPedLG)->Fill(towerId, meanPed);
468 GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
469 }
470 } // gain==0
471 else if ( in.IsHighGain() ) {
472 //fill the high gain ones
473 nTotalSMHG[iSM]++; // one more channel found
474 int signal = max - min;
475 // only fill the max-min signal info and maxTime, if the
476 // signal was in the selected range
477 if ( (signal > fMinSignalHG) && (signal < fMaxSignalHG) ) {
478 GetRawsData(kSigHG)->Fill(towerId, signal);
479 GetRawsData(kTimeHG)->Fill(towerId, maxTime);
480 GetRawsData(kTowerHG)->Fill(towerId);
481 } // signal
482 if (nSum>0) { // only fill pedestal info in case it could be calculated
483 GetRawsData(kPedHG)->Fill(towerId, meanPed);
484 GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
485 }
486 }
487 } // low or high gain
488 // TRU
489 else if ( in.IsTRUData() ) {
490 // for TRU data, the mapping class holds the TRU internal 2x2 number (0..95) in the Column var..
491 int iTRU = iRCU; //TRU0 is from RCU0, TRU1 from RCU1
492 if (iRCU>0 && in.GetBranch()>0) iTRU=2; // TRU2 is from branch B on RCU1
493 int TRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU
494 + in.GetColumn();
495
496 //fill the low gain histograms, and counters
497 nTotalSMTRU[iSM]++; // one more channel found
498 GetRawsData(kSigTRU)->Fill(TRU2x2Id, max - min);
499 GetRawsData(kTimeTRU)->Fill(TRU2x2Id, maxTime);
500 if (nSum>0) { // only fill pedestal info in case it could be calculated
501 GetRawsData(kPedTRU)->Fill(TRU2x2Id, meanPed);
502 GetRawsData(kPedRMSTRU)->Fill(TRU2x2Id, rmsPed);
503 }
504 }
505 // LED Mon
506 else if ( in.IsLEDMonData() ) {
507 // for LED Mon data, the mapping class holds the gain info in the Row variable
508 // and the Strip number in the Column..
509 int gain = in.GetRow();
510 int stripId = iSM*nStripsPerSM + in.GetColumn();
511
512 if ( gain == 0 ) {
513 //fill the low gain histograms, and counters
514 nTotalSMLGLEDMon[iSM]++; // one more channel found
515 GetRawsData(kSigLGLEDMon)->Fill(stripId, max - min);
516 GetRawsData(kTimeLGLEDMon)->Fill(stripId, maxTime);
517 if (nSum>0) { // only fill pedestal info in case it could be calculated
518 GetRawsData(kPedLGLEDMon)->Fill(stripId, meanPed);
519 GetRawsData(kPedRMSLGLEDMon)->Fill(stripId, rmsPed);
520 }
521 } // gain==0
522 else if ( gain == 1 ) {
523 //fill the high gain ones
524 nTotalSMHGLEDMon[iSM]++; // one more channel found
525 GetRawsData(kSigHGLEDMon)->Fill(stripId, max - min);
526 GetRawsData(kTimeHGLEDMon)->Fill(stripId, maxTime);
527 if (nSum>0) { // only fill pedestal info in case it could be calculated
528 GetRawsData(kPedHGLEDMon)->Fill(stripId, meanPed);
529 GetRawsData(kPedRMSHGLEDMon)->Fill(stripId, rmsPed);
530 }
531 } // low or high gain
532 } // LEDMon
533
534 } // SM index OK
535
536 } // nsamples>0 check, some data found for this channel; not only trailer/header
537 }// end while over channel
538
539 }//end while over DDL's, of input stream
540
541 // let's also fill the SM and event counter histograms
542 int nTotalHG = 0;
543 int nTotalLG = 0;
544 int nTotalTRU = 0;
545 int nTotalHGLEDMon = 0;
546 int nTotalLGLEDMon = 0;
547 for (iSM=0; iSM<fSuperModules; iSM++) {
548 nTotalLG += nTotalSMLG[iSM];
549 nTotalHG += nTotalSMHG[iSM];
550 nTotalTRU += nTotalSMTRU[iSM];
551 nTotalLG += nTotalSMLGLEDMon[iSM];
552 nTotalHG += nTotalSMHGLEDMon[iSM];
553 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
554 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
555 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
556 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
557 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
558 }
559 GetRawsData(kNtotLG)->Fill(nTotalLG);
560 GetRawsData(kNtotHG)->Fill(nTotalHG);
561 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
562 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
563 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
564
565 // just in case the next rawreader consumer forgets to reset; let's do it here again..
566 rawReader->Reset() ;
567
568 return;
569}
570
571//____________________________________________________________________________
572void AliEMCALQADataMakerRec::MakeDigits()
573{
574 // makes data from Digits
575
576 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
577 TIter next(fDigitsArray) ;
578 AliEMCALDigit * digit ;
579 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
580 GetDigitsData(0)->Fill( digit->GetAmp()) ;
581 }
582
583}
584
585//____________________________________________________________________________
586void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
587{
588 // makes data from Digit Tree
589 if (fDigitsArray)
590 fDigitsArray->Clear() ;
591 else
592 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
593
594 TBranch * branch = digitTree->GetBranch("EMCAL") ;
595 if ( ! branch ) {
596 AliWarning("EMCAL branch in Digit Tree not found") ;
597 } else {
598 branch->SetAddress(&fDigitsArray) ;
599 branch->GetEntry(0) ;
600 MakeDigits() ;
601 }
602
603}
604
605//____________________________________________________________________________
606void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
607{
608 // makes data from RecPoints
609 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
610 if (!emcbranch) {
611 AliError("can't get the branch with the EMCAL clusters !");
612 return;
613 }
614
615 TObjArray * emcrecpoints = new TObjArray(100) ;
616 emcbranch->SetAddress(&emcrecpoints);
617 emcbranch->GetEntry(0);
618
619 GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ;
620 TIter next(emcrecpoints) ;
621 AliEMCALRecPoint * rp ;
622 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
623 GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
624 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
625 }
626 emcrecpoints->Delete();
627 delete emcrecpoints;
628
629}
630
631//____________________________________________________________________________
632void AliEMCALQADataMakerRec::StartOfDetectorCycle()
633{
634 //Detector specific actions at start of cycle
635
636}
637