]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALQADataMakerRec.cxx
Changes required in QA for the amoreQAshifter agent; change units of titles from...
[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 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
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
55 ClassImp(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 //____________________________________________________________________________ 
70 AliEMCALQADataMakerRec::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 //__________________________________________________________________
84 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
85 {
86   // Equal operator.
87   this->~AliEMCALQADataMakerRec();
88   new(this) AliEMCALQADataMakerRec(qadm);
89   return *this;
90 }
91  
92 //____________________________________________________________________________ 
93 void 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 //____________________________________________________________________________ 
105 void 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 //____________________________________________________________________________ 
130 void 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 //____________________________________________________________________________ 
145 void 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 //____________________________________________________________________________ 
166 void 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 //____________________________________________________________________________
320 void 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 //____________________________________________________________________________
345 void 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"); 
363
364   // setup
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
369
370   int sampleMin = 0; 
371   int sampleMax = AliEMCALGeoParams::fgkSampleMax; // 0x3ff = 1023 = 10-bit range
372
373   // for the pedestal calculation
374   Bool_t selectPedestalSamples = kTRUE;
375
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};
382
383   // indices for the reading
384   int iSM = 0;
385   int sample = 0;
386   int time = 0;
387   // counters, on sample level
388   int i = 0; // the sample number in current event.
389   int maxTime = 0;
390   int startBin = 0;
391
392   // calc. quantities
393   double meanPed = 0, squaredMean = 0, rmsPed = 0;
394
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()) {
399
400       // counters
401       int max = sampleMin, min = sampleMax; // min and max sample values
402       int nsamples = 0;
403
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
408       
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++) {
414           sample = sig[i];
415           time = startBin--;
416
417           // check if it's a min or max value
418           if (sample < min) min = sample;
419           if (sample > max) {
420             max = sample;
421             maxTime = time;
422           }
423
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 
427             sampleSum += sample;
428             squaredSampleSum += sample*sample;
429             nSum++;
430           }
431           
432         } // loop over samples in bunch
433       } // loop over bunches
434     
435       if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
436
437       // calculate pedesstal estimate: mean of possibly selected samples
438       if (nSum > 0) {
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); 
443       }
444       else {
445         meanPed = 0;
446         squaredMean = 0;
447         rmsPed  = 0;
448       }
449
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
452
453       if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
454
455         if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
456           int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
457         
458
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);
468             }
469           } // gain==0
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);
480             } // signal
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);
484             }
485           }
486         } // low or high gain
487         // TRU
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 
493             + in.GetColumn();
494
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);
502           }
503         }
504         // LED Mon
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();
510           
511           if ( gain == 0 ) { 
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);
519             }
520           } // gain==0
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);
529             }
530           } // low or high gain
531         } // LEDMon
532
533       } // SM index OK
534
535       } // nsamples>0 check, some data found for this channel; not only trailer/header
536     }// end while over channel 
537    
538   }//end while over DDL's, of input stream 
539
540   // let's also fill the SM and event counter histograms
541   int nTotalHG = 0;
542   int nTotalLG = 0;
543   int nTotalTRU = 0;
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]); 
557   }
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);
563
564   // just in case the next rawreader consumer forgets to reset; let's do it here again..
565   rawReader->Reset() ;
566
567   return;
568 }
569
570 //____________________________________________________________________________
571 void AliEMCALQADataMakerRec::MakeDigits()
572 {
573   // makes data from Digits
574
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()) ;
580   }  
581   
582 }
583
584 //____________________________________________________________________________
585 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
586 {
587   // makes data from Digit Tree
588   if (fDigitsArray) 
589     fDigitsArray->Clear() ; 
590   else
591     fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; 
592   
593   TBranch * branch = digitTree->GetBranch("EMCAL") ;
594   if ( ! branch ) {
595     AliWarning("EMCAL branch in Digit Tree not found") ; 
596   } else {
597     branch->SetAddress(&fDigitsArray) ;
598     branch->GetEntry(0) ; 
599     MakeDigits() ; 
600   }
601   
602 }
603
604 //____________________________________________________________________________
605 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
606 {
607   // makes data from RecPoints
608   TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
609   if (!emcbranch) { 
610     AliError("can't get the branch with the EMCAL clusters !");
611     return;
612   }
613   
614   TObjArray * emcrecpoints = new TObjArray(100) ;
615   emcbranch->SetAddress(&emcrecpoints);
616   emcbranch->GetEntry(0);
617   
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());
624   }
625   emcrecpoints->Delete();
626   delete emcrecpoints;
627   
628 }
629
630 //____________________________________________________________________________ 
631 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
632 {
633   //Detector specific actions at start of cycle
634   
635 }
636