]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALQADataMakerRec.cxx
LEDMON->LEDMon for LG histo titles also + nsamples check needed for zero-suppressed...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMakerRec.cxx
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 /*
16 Based on the QA code for PHOS written by Yves Schutz July 2007
17
18 Authors:  J.Klay (Cal Poly) May 2008
19           S. Salur LBL April 2008
20
21 */
22
23 // --- ROOT system ---
24 #include <TClonesArray.h>
25 #include <TFile.h> 
26 #include <TH1F.h> 
27 #include <TH1I.h> 
28 #include <TH2F.h> 
29 #include <TProfile.h> 
30
31 // --- Standard library ---
32
33
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
36 #include "AliESDCaloCells.h"
37 #include "AliESDEvent.h"
38 #include "AliLog.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"
49
50 ClassImp(AliEMCALQADataMakerRec)
51            
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),
58     fMinSignalHG(0),
59     fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax)
60 {
61   // ctor
62 }
63
64 //____________________________________________________________________________ 
65 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
66   AliQADataMakerRec(), 
67   fSuperModules(qadm.GetSuperModules()), 
68   fFirstPedestalSample(qadm.GetFirstPedestalSample()), 
69   fLastPedestalSample(qadm.GetLastPedestalSample()),  
70   fMinSignalHG(qadm.GetMinSignalHG()),
71   fMaxSignalHG(qadm.GetMaxSignalHG())
72 {
73   //copy ctor 
74   SetName((const char*)qadm.GetName()) ; 
75   SetTitle((const char*)qadm.GetTitle()); 
76 }
77
78 //__________________________________________________________________
79 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
80 {
81   // Equal operator.
82   this->~AliEMCALQADataMakerRec();
83   new(this) AliEMCALQADataMakerRec(qadm);
84   return *this;
85 }
86  
87 //____________________________________________________________________________ 
88 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
89 {
90   //Detector specific actions at end of cycle
91   // do the QA checking
92   AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;  
93 }
94
95 //____________________________________________________________________________ 
96 void AliEMCALQADataMakerRec::InitESDs()
97 {
98   //Create histograms to controll ESD
99   const Bool_t expert   = kTRUE ; 
100   const Bool_t image    = kTRUE ; 
101   
102   TH1F * h1 = new TH1F("hESDCaloClusterE",  "ESDs CaloCluster energy in EMCAL;Energy [MeV];Counts",    200, 0., 20.) ; 
103   h1->Sumw2() ;
104   Add2ESDsList(h1, kESDCaloClusE, !expert, image)  ;                                                     
105
106   TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0,  100) ; 
107   h2->Sumw2() ;
108   Add2ESDsList(h2, kESDCaloClusM, !expert, image)  ;
109
110   TH1F * h3 = new TH1F("hESDCaloCellA",  "ESDs CaloCell amplitude in EMCAL;Energy [MeV];Counts",    500, 0., 250.) ; 
111   h3->Sumw2() ;
112   Add2ESDsList(h3, kESDCaloCellA, !expert, image)  ;  
113  
114   TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0,  1000) ; 
115   h4->Sumw2() ;
116   Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
117         
118 }
119
120 //____________________________________________________________________________ 
121 void AliEMCALQADataMakerRec::InitDigits()
122 {
123   // create Digits histograms in Digits subdir
124   const Bool_t expert   = kTRUE ; 
125   const Bool_t image    = kTRUE ; 
126   
127   TH1I * h0 = new TH1I("hEmcalDigits",    "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts",    500, 0, 500) ; 
128   h0->Sumw2() ;
129   Add2DigitsList(h0, 0, !expert, image) ;
130   TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ; 
131   h1->Sumw2() ;
132   Add2DigitsList(h1, 1, !expert, image) ;
133 }
134
135 //____________________________________________________________________________ 
136 void AliEMCALQADataMakerRec::InitRecPoints()
137 {
138   // create Reconstructed Points histograms in RecPoints subdir
139   const Bool_t expert   = kTRUE ; 
140   const Bool_t image    = kTRUE ; 
141   
142   TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [MeV];Counts",200, 0.,20.); //GeV
143   h0->Sumw2();
144   Add2RecPointsList(h0,kRecPE, !expert, image);
145
146   TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
147   h1->Sumw2();
148   Add2RecPointsList(h1,kRecPM, !expert, image);
149
150   TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
151   h2->Sumw2();
152   Add2RecPointsList(h2,kRecPDigM, !expert, image);
153
154 }
155
156 //____________________________________________________________________________ 
157 void AliEMCALQADataMakerRec::InitRaws()
158 {
159   // create Raws histograms in Raws subdir
160    const Bool_t expert   = kTRUE ; 
161    const Bool_t saveCorr = kTRUE ; 
162    const Bool_t image    = kTRUE ; 
163
164   int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
165   int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
166
167   // counter info: number of channels per event (bins are SM index)
168   TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
169                                fSuperModules, -0.5, fSuperModules-0.5) ;
170   Add2RawsList(h0, kNsmodLG, !expert, image, !saveCorr) ;
171   TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",  
172                                fSuperModules, -0.5, fSuperModules-0.5) ; 
173   Add2RawsList(h1, kNsmodHG, !expert, image, !saveCorr) ;
174
175   // where did max sample occur? (bins are towers)
176   TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
177                                nTot, -0.5, nTot-0.5) ;
178   Add2RawsList(h2, kTimeLG, !expert, image, !saveCorr) ;
179   TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
180                                nTot, -0.5, nTot-0.5) ;
181   Add2RawsList(h3, kTimeHG, !expert, image, !saveCorr) ;
182
183   // how much above pedestal was the max sample?  (bins are towers)
184   TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]", 
185                                nTot, -0.5, nTot-0.5) ;
186   Add2RawsList(h4, kSigLG, !expert, image, !saveCorr) ;
187   TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
188                                nTot, -0.5, nTot-0.5) ;
189   Add2RawsList(h5, kSigHG, !expert, image, !saveCorr) ;
190
191   // total counter: channels per event
192   TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
193   h6->Sumw2() ;
194   Add2RawsList(h6, kNtotLG, !expert, image, !saveCorr) ;
195   TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
196   h7->Sumw2() ;
197   Add2RawsList(h7, kNtotHG, !expert, image, !saveCorr) ;
198
199   // pedestal (bins are towers)
200   TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]", 
201                                nTot, -0.5, nTot-0.5) ;
202   Add2RawsList(h8, kPedLG, !expert, image, !saveCorr) ;
203   TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
204                                nTot, -0.5, nTot-0.5) ;
205   Add2RawsList(h9, kPedHG, !expert, image, !saveCorr) ;
206
207   // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are towers)
208   TProfile * h10 = new TProfile("hLowEmcalRawPedRMS", "Low Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]", 
209                                 nTot, -0.5, nTot-0.5) ;
210   Add2RawsList(h10, kPedRMSLG, !expert, image, !saveCorr) ;
211   TProfile * h11 = new TProfile("hHighEmcalRawPedRMS", "High Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
212                                 nTot, -0.5, nTot-0.5) ;
213   Add2RawsList(h11, kPedRMSHG, !expert, image, !saveCorr) ;
214
215
216   // now repeat the same for TRU and LEDMon data
217   int nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
218
219   // counter info: number of channels per event (bins are SM index)
220   TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
221                                 fSuperModules, -0.5, fSuperModules-0.5) ;
222   Add2RawsList(hT0, kNsmodTRU, !expert, image, !saveCorr) ;
223
224   // where did max sample occur? (bins are TRU channels)
225   TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]", 
226                                 nTot2x2, -0.5, nTot2x2-0.5) ;
227   Add2RawsList(hT1, kTimeTRU, !expert, image, !saveCorr) ;
228
229   // how much above pedestal was the max sample?  (bins are TRU channels)
230   TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]", 
231                                 nTot2x2, -0.5, nTot2x2-0.5) ;
232   Add2RawsList(hT2, kSigTRU, !expert, image, !saveCorr) ;
233
234   // total counter: channels per event
235   TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
236   hT3->Sumw2() ;
237   Add2RawsList(hT3, kNtotTRU, !expert, image, !saveCorr) ;
238
239   // pedestal (bins are TRU channels)
240   TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]", 
241                                 nTot2x2, -0.5, nTot2x2-0.5) ;
242   Add2RawsList(hT4, kPedTRU, !expert, image, !saveCorr) ;
243
244   // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are TRU channels)
245   TProfile * hT5 = new TProfile("hTRUEmcalRawPedRMS", "TRU EMC: Pedestal RMS vs 2x2Id;2x2 Id;Width [ADC counts]", 
246                                 nTot2x2, -0.5, nTot2x2-0.5) ;
247   Add2RawsList(hT5, kPedRMSTRU, !expert, image, !saveCorr) ;
248
249   // and also LED Mon..
250   // LEDMon has both high and low gain channels, just as regular FEE/towers
251   int nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
252
253   // counter info: number of channels per event (bins are SM index)
254   TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
255                                fSuperModules, -0.5, fSuperModules-0.5) ;
256   Add2RawsList(hL0, kNsmodLGLEDMon, !expert, image, !saveCorr) ;
257   TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",  
258                                fSuperModules, -0.5, fSuperModules-0.5) ; 
259   Add2RawsList(hL1, kNsmodHGLEDMon, !expert, image, !saveCorr) ;
260
261   // where did max sample occur? (bins are strips)
262   TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
263                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
264   Add2RawsList(hL2, kTimeLGLEDMon, !expert, image, !saveCorr) ;
265   TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
266                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
267   Add2RawsList(hL3, kTimeHGLEDMon, !expert, image, !saveCorr) ;
268
269   // how much above pedestal was the max sample?  (bins are strips)
270   TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]", 
271                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
272   Add2RawsList(hL4, kSigLGLEDMon, !expert, image, !saveCorr) ;
273   TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
274                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
275   Add2RawsList(hL5, kSigHGLEDMon, !expert, image, !saveCorr) ;
276
277   // total counter: channels per event
278   TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
279   hL6->Sumw2() ;
280   Add2RawsList(hL6, kNtotLGLEDMon, !expert, image, !saveCorr) ;
281   TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
282   hL7->Sumw2() ;
283   Add2RawsList(hL7, kNtotHGLEDMon, !expert, image, !saveCorr) ;
284
285   // pedestal (bins are strips)
286   TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]", 
287                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
288   Add2RawsList(hL8, kPedLGLEDMon, !expert, image, !saveCorr) ;
289   TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
290                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
291   Add2RawsList(hL9, kPedHGLEDMon, !expert, image, !saveCorr) ;
292
293   // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are strips)
294   TProfile * hL10 = new TProfile("hLowLEDMonEmcalRawPedRMS", "LowLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]", 
295                                 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
296   Add2RawsList(hL10, kPedRMSLGLEDMon, !expert, image, !saveCorr) ;
297   TProfile * hL11 = new TProfile("hHighLEDMonEmcalRawPedRMS", "HighLEDMon Gain EMC: Pedestal RMS vs stripId;Strip Id;Width [ADC counts]",
298                                 nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
299   Add2RawsList(hL11, kPedRMSHGLEDMon, !expert, image, !saveCorr) ;
300   
301 }
302
303 //____________________________________________________________________________
304 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
305 {
306   // make QA data from ESDs
307
308   Int_t nTot = 0 ; 
309   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
310     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
311     if( clu->IsEMCAL() ) {
312       GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
313       nTot++ ;
314     } 
315   }
316   GetESDsData(kESDCaloClusM)->Fill(nTot) ;
317
318   //fill calo cells
319   AliESDCaloCells* cells = esd->GetEMCALCells();
320   GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
321
322   for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
323     GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
324   }
325
326 }
327
328 //____________________________________________________________________________
329 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
330 {
331   //Fill prepared histograms with Raw digit properties
332
333   //Raw histogram filling not yet implemented
334   //
335   //Need to figure out how to get the info we want without having to
336   //actually run Raw2Digits twice.
337   //I suspect what we actually want is a raw digits method, not a true
338   //emcal raw data method, but this doesn't seem to be allowed in
339   //AliQADataMakerRec.h
340
341   // For now, to avoid redoing the expensive signal fits we just
342   // look at max vs min of the signal spextra, a la online usage in
343   // AliCaloCalibPedestal
344
345   rawReader->Reset() ;
346   AliCaloRawStreamV3 in(rawReader,"EMCAL"); 
347
348   // setup
349   int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
350   int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
351   int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
352   int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
353
354   int sampleMin = 0; 
355   int sampleMax = AliEMCALGeoParams::fgkSampleMax; // 0x3ff = 1023 = 10-bit range
356
357   // for the pedestal calculation
358   Bool_t selectPedestalSamples = kTRUE;
359
360   // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
361   int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
362   int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
363   int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
364   int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
365   int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
366
367   // indices for the reading
368   int iSM = 0;
369   int sample = 0;
370   int time = 0;
371   // counters, on sample level
372   int i = 0; // the sample number in current event.
373   int maxTime = 0;
374   int startBin = 0;
375
376   // calc. quantities
377   double meanPed = 0, squaredMean = 0, rmsPed = 0;
378
379   // start loop over input stream  
380   while (in.NextDDL()) {
381     int iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
382     while (in.NextChannel()) {
383
384       // counters
385       int max = sampleMin, min = sampleMax; // min and max sample values
386       int nsamples = 0;
387
388       // for the pedestal calculation
389       int sampleSum = 0; // sum of samples
390       int squaredSampleSum = 0; // sum of samples squared
391       int nSum = 0; // number of samples in sum
392       
393       while (in.NextBunch()) {
394         const UShort_t *sig = in.GetSignals();
395         startBin = in.GetStartTimeBin();
396         nsamples += in.GetBunchLength();
397         for (i = 0; i < in.GetBunchLength(); i++) {
398           sample = sig[i];
399           time = startBin--;
400
401           // check if it's a min or max value
402           if (sample < min) min = sample;
403           if (sample > max) {
404             max = sample;
405             maxTime = time;
406           }
407
408           // should we add it for the pedestal calculation?
409           if ( (fFirstPedestalSample<=time && time<=fLastPedestalSample) || // sample time in range
410                !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all 
411             sampleSum += sample;
412             squaredSampleSum += sample*sample;
413             nSum++;
414           }
415           
416         } // loop over samples in bunch
417       } // loop over bunches
418     
419       if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
420
421       // calculate pedesstal estimate: mean of possibly selected samples
422       if (nSum > 0) {
423         meanPed = sampleSum / (1.0 * nSum);
424         squaredMean = squaredSampleSum / (1.0 * nSum);
425         // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
426         rmsPed = sqrt(squaredMean - meanPed*meanPed); 
427       }
428       else {
429         meanPed = 0;
430         squaredMean = 0;
431         rmsPed  = 0;
432       }
433
434       // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
435       iSM = in.GetModule(); //The modules are numbered starting from 0
436
437       if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
438
439         if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
440           int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
441
442           if ( in.IsLowGain() ) { 
443             //fill the low gain histograms, and counters
444             nTotalSMLG[iSM]++; // one more channel found
445             GetRawsData(kSigLG)->Fill(towerId, max - min);
446             GetRawsData(kTimeLG)->Fill(towerId, maxTime);
447             if (nSum>0) { // only fill pedestal info in case it could be calculated
448               GetRawsData(kPedLG)->Fill(towerId, meanPed);
449               GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
450             }
451           } // gain==0
452           else if ( in.IsHighGain() ) {         
453             //fill the high gain ones
454             nTotalSMHG[iSM]++; // one more channel found
455             int signal = max - min;
456             // only fill the max-min signal info and maxTime, if the
457             // signal was in the selected range 
458             if ( (signal > fMinSignalHG) && (signal < fMaxSignalHG) ) { 
459               GetRawsData(kSigHG)->Fill(towerId, signal);
460               GetRawsData(kTimeHG)->Fill(towerId, maxTime);
461             } // signal
462             if (nSum>0) { // only fill pedestal info in case it could be calculated
463               GetRawsData(kPedHG)->Fill(towerId, meanPed);
464               GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
465             }
466           }
467         } // low or high gain
468         // TRU
469         else if ( in.IsTRUData() ) {
470           // for TRU data, the mapping class holds the TRU internal 2x2 number (0..95) in the Column var..
471           int iTRU = iRCU; //TRU0 is from RCU0, TRU1 from RCU1
472           if (iRCU>0 && in.GetBranch()>0) iTRU=2; // TRU2 is from branch B on RCU1
473           int TRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU 
474             + in.GetColumn();
475
476           //fill the low gain histograms, and counters
477           nTotalSMTRU[iSM]++; // one more channel found
478           GetRawsData(kSigTRU)->Fill(TRU2x2Id, max - min);
479           GetRawsData(kTimeTRU)->Fill(TRU2x2Id, maxTime);
480           if (nSum>0) { // only fill pedestal info in case it could be calculated
481             GetRawsData(kPedTRU)->Fill(TRU2x2Id, meanPed);
482             GetRawsData(kPedRMSTRU)->Fill(TRU2x2Id, rmsPed);
483           }
484         }
485         // LED Mon
486         else if ( in.IsLEDMonData() ) {
487           // for LED Mon data, the mapping class holds the gain info in the Row variable
488           // and the Strip number in the Column..
489           int gain = in.GetRow(); 
490           int stripId = iSM*nStripsPerSM + in.GetColumn();
491           
492           if ( gain == 0 ) { 
493             //fill the low gain histograms, and counters
494             nTotalSMLGLEDMon[iSM]++; // one more channel found
495             GetRawsData(kSigLGLEDMon)->Fill(stripId, max - min);
496             GetRawsData(kTimeLGLEDMon)->Fill(stripId, maxTime);
497             if (nSum>0) { // only fill pedestal info in case it could be calculated
498               GetRawsData(kPedLGLEDMon)->Fill(stripId, meanPed);
499               GetRawsData(kPedRMSLGLEDMon)->Fill(stripId, rmsPed);
500             }
501           } // gain==0
502           else if ( gain == 1 ) {               
503             //fill the high gain ones
504             nTotalSMHGLEDMon[iSM]++; // one more channel found
505             GetRawsData(kSigHGLEDMon)->Fill(stripId, max - min);
506             GetRawsData(kTimeHGLEDMon)->Fill(stripId, maxTime);
507             if (nSum>0) { // only fill pedestal info in case it could be calculated
508               GetRawsData(kPedHGLEDMon)->Fill(stripId, meanPed);
509               GetRawsData(kPedRMSHGLEDMon)->Fill(stripId, rmsPed);
510             }
511           } // low or high gain
512         } // LEDMon
513
514       } // SM index OK
515
516       } // nsamples>0 check, some data found for this channel; not only trailer/header
517     }// end while over channel 
518    
519   }//end while over DDL's, of input stream 
520
521   // let's also fill the SM and event counter histograms
522   int nTotalHG = 0;
523   int nTotalLG = 0;
524   int nTotalTRU = 0;
525   int nTotalHGLEDMon = 0;
526   int nTotalLGLEDMon = 0;
527   for (iSM=0; iSM<fSuperModules; iSM++) {  
528     nTotalLG += nTotalSMLG[iSM]; 
529     nTotalHG += nTotalSMHG[iSM]; 
530     nTotalTRU += nTotalSMTRU[iSM]; 
531     nTotalLG += nTotalSMLGLEDMon[iSM]; 
532     nTotalHG += nTotalSMHGLEDMon[iSM]; 
533     GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]); 
534     GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]); 
535     GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]); 
536     GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]); 
537     GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]); 
538   }
539   GetRawsData(kNtotLG)->Fill(nTotalLG);
540   GetRawsData(kNtotHG)->Fill(nTotalHG);
541   GetRawsData(kNtotTRU)->Fill(nTotalTRU);
542   GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
543   GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
544
545   // just in case the next rawreader consumer forgets to reset; let's do it here again..
546   rawReader->Reset() ;
547
548   return;
549 }
550
551 //____________________________________________________________________________
552 void AliEMCALQADataMakerRec::MakeDigits()
553 {
554   // makes data from Digits
555
556   GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ; 
557   TIter next(fDigitsArray) ; 
558   AliEMCALDigit * digit ; 
559   while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
560     GetDigitsData(0)->Fill( digit->GetAmp()) ;
561   }  
562   
563 }
564
565 //____________________________________________________________________________
566 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
567 {
568   // makes data from Digit Tree
569   if (fDigitsArray) 
570     fDigitsArray->Clear() ; 
571   else
572     fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; 
573   
574   TBranch * branch = digitTree->GetBranch("EMCAL") ;
575   if ( ! branch ) {
576     AliWarning("EMCAL branch in Digit Tree not found") ; 
577   } else {
578     branch->SetAddress(&fDigitsArray) ;
579     branch->GetEntry(0) ; 
580     MakeDigits() ; 
581   }
582   
583 }
584
585 //____________________________________________________________________________
586 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
587 {
588   // makes data from RecPoints
589   TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
590   if (!emcbranch) { 
591     AliError("can't get the branch with the EMCAL clusters !");
592     return;
593   }
594   
595   TObjArray * emcrecpoints = new TObjArray(100) ;
596   emcbranch->SetAddress(&emcrecpoints);
597   emcbranch->GetEntry(0);
598   
599   GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ; 
600   TIter next(emcrecpoints) ; 
601   AliEMCALRecPoint * rp ; 
602   while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
603     GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
604     GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
605   }
606   emcrecpoints->Delete();
607   delete emcrecpoints;
608   
609 }
610
611 //____________________________________________________________________________ 
612 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
613 {
614   //Detector specific actions at start of cycle
615   
616 }
617