]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALQADataMakerRec.cxx
a3d151771e8ae0966aa2ca5c2e2df13e57f2d7ab
[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 #include "AliRawEventHeaderBase.h"
55
56 #include "AliCaloBunchInfo.h"
57 #include "AliCaloFitResults.h"
58 #include "AliCaloRawAnalyzerFastFit.h"
59 #include "AliCaloRawAnalyzerNN.h"
60 #include "AliCaloRawAnalyzerLMS.h"
61 #include "AliCaloRawAnalyzerPeakFinder.h"
62 #include "AliCaloRawAnalyzerCrude.h"
63
64 ClassImp(AliEMCALQADataMakerRec)
65            
66 //____________________________________________________________________________ 
67 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) : 
68   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
69   fFittingAlgorithm(0),
70   fRawAnalyzer(0),
71   fRawAnalyzerTRU(0),
72   fSuperModules(4), // FIXME!!! number of SuperModules; 4 for 2009; update default to 12 for later runs..
73   fFirstPedestalSample(0),
74   fLastPedestalSample(3),
75   fFirstPedestalSampleTRU(0),
76   fLastPedestalSampleTRU(3),
77   fMinSignalLG(0),
78   fMaxSignalLG(AliEMCALGeoParams::fgkSampleMax),
79   fMinSignalHG(0),
80   fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax),
81   fMinSignalTRU(0),
82   fMaxSignalTRU(AliEMCALGeoParams::fgkSampleMax),
83   fMinSignalLGLEDMon(0),
84   fMaxSignalLGLEDMon(AliEMCALGeoParams::fgkSampleMax),
85   fMinSignalHGLEDMon(0),
86   fMaxSignalHGLEDMon(AliEMCALGeoParams::fgkSampleMax)
87 {
88   // ctor
89   SetFittingAlgorithm(fitAlgo);
90   fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
91   fRawAnalyzerTRU->SetFixTau(kTRUE); 
92   fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
93 }
94
95 //____________________________________________________________________________ 
96 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
97   AliQADataMakerRec(), 
98   fFittingAlgorithm(0),
99   fRawAnalyzer(0),
100   fRawAnalyzerTRU(0),
101   fSuperModules(qadm.GetSuperModules()), 
102   fFirstPedestalSample(qadm.GetFirstPedestalSample()), 
103   fLastPedestalSample(qadm.GetLastPedestalSample()),  
104   fFirstPedestalSampleTRU(qadm.GetFirstPedestalSampleTRU()), 
105   fLastPedestalSampleTRU(qadm.GetLastPedestalSampleTRU()),  
106   fMinSignalLG(qadm.GetMinSignalLG()),
107   fMaxSignalLG(qadm.GetMaxSignalLG()),
108   fMinSignalHG(qadm.GetMinSignalHG()),
109   fMaxSignalHG(qadm.GetMaxSignalHG()),
110   fMinSignalTRU(qadm.GetMinSignalTRU()),
111   fMaxSignalTRU(qadm.GetMaxSignalTRU()),
112   fMinSignalLGLEDMon(qadm.GetMinSignalLGLEDMon()),
113   fMaxSignalLGLEDMon(qadm.GetMaxSignalLGLEDMon()),
114   fMinSignalHGLEDMon(qadm.GetMinSignalHGLEDMon()),
115   fMaxSignalHGLEDMon(qadm.GetMaxSignalHGLEDMon())
116 {
117   //copy ctor 
118   SetName((const char*)qadm.GetName()) ; 
119   SetTitle((const char*)qadm.GetTitle()); 
120   SetFittingAlgorithm(qadm.GetFittingAlgorithm());
121   fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
122   fRawAnalyzerTRU->SetFixTau(kTRUE); 
123   fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
124 }
125
126 //__________________________________________________________________
127 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
128 {
129   // Equal operator.
130   this->~AliEMCALQADataMakerRec();
131   new(this) AliEMCALQADataMakerRec(qadm);
132   return *this;
133 }
134  
135 //____________________________________________________________________________ 
136 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
137 {
138   //Detector specific actions at end of cycle
139         
140 //  if(fCycleCounter)
141 //        GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
142
143   // do the QA checking
144   AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;  
145 }
146
147 //____________________________________________________________________________ 
148 void AliEMCALQADataMakerRec::InitESDs()
149 {
150   //Create histograms to controll ESD
151   const Bool_t expert   = kTRUE ; 
152   const Bool_t image    = kTRUE ; 
153   
154   TH1F * h1 = new TH1F("hESDCaloClusterE",  "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts",    200, 0., 100.) ; 
155   h1->Sumw2() ;
156   Add2ESDsList(h1, kESDCaloClusE, !expert, image)  ;                                                     
157
158   TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0,  100) ; 
159   h2->Sumw2() ;
160   Add2ESDsList(h2, kESDCaloClusM, !expert, image)  ;
161
162   TH1F * h3 = new TH1F("hESDCaloCellA",  "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts",    500, 0., 50.) ; 
163   h3->Sumw2() ;
164   Add2ESDsList(h3, kESDCaloCellA, !expert, image)  ;  
165  
166   TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0,  1000) ; 
167   h4->Sumw2() ;
168   Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
169         
170 }
171
172 //____________________________________________________________________________ 
173 void AliEMCALQADataMakerRec::InitDigits()
174 {
175   // create Digits histograms in Digits subdir
176   const Bool_t expert   = kTRUE ; 
177   const Bool_t image    = kTRUE ; 
178   
179   TH1I * h0 = new TH1I("hEmcalDigits",    "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts",    500, 0, 500) ; 
180   h0->Sumw2() ;
181   Add2DigitsList(h0, 0, !expert, image) ;
182   TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ; 
183   h1->Sumw2() ;
184   Add2DigitsList(h1, 1, !expert, image) ;
185 }
186
187 //____________________________________________________________________________ 
188 void AliEMCALQADataMakerRec::InitRecPoints()
189 {
190   // create Reconstructed Points histograms in RecPoints subdir
191   const Bool_t expert   = kTRUE ; 
192   const Bool_t image    = kTRUE ; 
193   
194   TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
195   h0->Sumw2();
196   Add2RecPointsList(h0,kRecPE, !expert, image);
197
198   TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
199   h1->Sumw2();
200   Add2RecPointsList(h1,kRecPM, !expert, image);
201
202   TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
203   h2->Sumw2();
204   Add2RecPointsList(h2,kRecPDigM, !expert, image);
205
206 }
207
208 //____________________________________________________________________________ 
209 void AliEMCALQADataMakerRec::InitRaws()
210 {
211   // create Raws histograms in Raws subdir
212    const Bool_t expert   = kTRUE ; 
213    const Bool_t saveCorr = kTRUE ; 
214    const Bool_t image    = kTRUE ; 
215
216   int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
217   int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
218
219   // counter info: number of channels per event (bins are SM index)
220   TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
221                                fSuperModules, -0.5, fSuperModules-0.5) ;
222   Add2RawsList(h0, kNsmodLG, expert, image, !saveCorr) ;
223   TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",  
224                                fSuperModules, -0.5, fSuperModules-0.5) ; 
225   Add2RawsList(h1, kNsmodHG, expert, image, !saveCorr) ;
226
227   // where did max sample occur? (bins are towers)
228   TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
229                                nTot, -0.5, nTot-0.5) ;
230   Add2RawsList(h2, kTimeLG, expert, image, !saveCorr) ;
231   TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
232                                nTot, -0.5, nTot-0.5) ;
233   Add2RawsList(h3, kTimeHG, expert, image, !saveCorr) ;
234
235   // how much above pedestal was the max sample?  (bins are towers)
236   TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]", 
237                                nTot, -0.5, nTot-0.5) ;
238   Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
239   TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
240                                nTot, -0.5, nTot-0.5) ;
241   Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
242
243   // total counter: channels per event
244   TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
245   h6->Sumw2() ;
246   Add2RawsList(h6, kNtotLG, expert, image, !saveCorr) ;
247   TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
248   h7->Sumw2() ;
249   Add2RawsList(h7, kNtotHG, expert, image, !saveCorr) ;
250
251   // pedestal (bins are towers)
252   TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]", 
253                                nTot, -0.5, nTot-0.5) ;
254   Add2RawsList(h8, kPedLG, expert, image, !saveCorr) ;
255   TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
256                                nTot, -0.5, nTot-0.5) ;
257   Add2RawsList(h9, kPedHG, expert, image, !saveCorr) ;
258         
259  //number of events per tower, for shifter fast check   
260   TH1I * h12 = new TH1I("hTowerHG", "High Gains on the Tower;Tower", nTot,0, nTot) ;
261   h12->Sumw2() ;
262   Add2RawsList(h12, kTowerHG, !expert, image, !saveCorr) ;
263   TH1I * h13 = new TH1I("hTowerLG", "Low Gains on the Tower;Tower", nTot,0, nTot) ;
264   h13->Sumw2() ;
265   Add2RawsList(h13, kTowerLG, !expert, image, !saveCorr) ;              
266
267   // now repeat the same for TRU and LEDMon data
268   int nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
269
270   // counter info: number of channels per event (bins are SM index)
271   TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
272                                 fSuperModules, -0.5, fSuperModules-0.5) ;
273   Add2RawsList(hT0, kNsmodTRU, expert, image, !saveCorr) ;
274
275   // where did max sample occur? (bins are TRU channels)
276   TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]", 
277                                 nTot2x2, -0.5, nTot2x2-0.5) ;
278   Add2RawsList(hT1, kTimeTRU, expert, image, !saveCorr) ;
279
280   // how much above pedestal was the max sample?  (bins are TRU channels)
281   TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]", 
282                                 nTot2x2, -0.5, nTot2x2-0.5) ;
283   Add2RawsList(hT2, kSigTRU, expert, image, !saveCorr) ;
284
285   // total counter: channels per event
286   TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
287   hT3->Sumw2() ;
288   Add2RawsList(hT3, kNtotTRU, expert, image, !saveCorr) ;
289
290   // pedestal (bins are TRU channels)
291   TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]", 
292                                 nTot2x2, -0.5, nTot2x2-0.5) ;
293   Add2RawsList(hT4, kPedTRU, expert, image, !saveCorr) ;
294
295   // and also LED Mon..
296   // LEDMon has both high and low gain channels, just as regular FEE/towers
297   int nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
298
299   // counter info: number of channels per event (bins are SM index)
300   TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
301                                fSuperModules, -0.5, fSuperModules-0.5) ;
302   Add2RawsList(hL0, kNsmodLGLEDMon, expert, image, !saveCorr) ;
303   TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",  
304                                fSuperModules, -0.5, fSuperModules-0.5) ; 
305   Add2RawsList(hL1, kNsmodHGLEDMon, expert, image, !saveCorr) ;
306
307   // where did max sample occur? (bins are strips)
308   TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
309                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
310   Add2RawsList(hL2, kTimeLGLEDMon, expert, image, !saveCorr) ;
311   TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
312                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
313   Add2RawsList(hL3, kTimeHGLEDMon, expert, image, !saveCorr) ;
314
315   // how much above pedestal was the max sample?  (bins are strips)
316   TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]", 
317                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
318   Add2RawsList(hL4, kSigLGLEDMon, expert, image, !saveCorr) ;
319   TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
320                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
321   Add2RawsList(hL5, kSigHGLEDMon, expert, image, !saveCorr) ;
322
323   // total counter: channels per event
324   TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
325   hL6->Sumw2() ;
326   Add2RawsList(hL6, kNtotLGLEDMon, expert, image, !saveCorr) ;
327   TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
328   hL7->Sumw2() ;
329   Add2RawsList(hL7, kNtotHGLEDMon, expert, image, !saveCorr) ;
330
331   // pedestal (bins are strips)
332   TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]", 
333                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
334   Add2RawsList(hL8, kPedLGLEDMon, expert, image, !saveCorr) ;
335   TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
336                                nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
337   Add2RawsList(hL9, kPedHGLEDMon, expert, image, !saveCorr) ;
338
339 }
340
341 //____________________________________________________________________________
342 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
343 {
344   // make QA data from ESDs
345
346   Int_t nTot = 0 ; 
347   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
348     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
349     if( clu->IsEMCAL() ) {
350       GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
351       nTot++ ;
352     } 
353   }
354   GetESDsData(kESDCaloClusM)->Fill(nTot) ;
355
356   //fill calo cells
357   AliESDCaloCells* cells = esd->GetEMCALCells();
358   GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
359
360   for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
361     GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
362   }
363
364 }
365
366 //____________________________________________________________________________
367 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
368 {
369   //Fill prepared histograms with Raw digit properties
370
371   // setup
372   rawReader->Reset() ;
373   AliCaloRawStreamV3 in(rawReader,"EMCAL"); 
374   rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's 
375
376   AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ; 
377   if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) { 
378     SetEventSpecie(AliRecoParam::kCalib) ;
379   }
380   
381   fRawAnalyzer->SetIsZeroSuppressed(true); // TMP - should use stream->IsZeroSuppressed(), or altro cfg registers later
382
383   int nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
384   int nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
385   int nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
386   int n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
387
388   // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
389   int nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
390   int nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
391   int nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
392   int nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
393   int nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
394
395   // start loop over input stream  
396   int iSM = 0;
397   while (in.NextDDL()) {
398     iSM = in.GetModule(); // SuperModule
399     if (iSM>=0 && iSM<fSuperModules) { // valid module reading
400
401       int iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
402
403       while (in.NextChannel()) {
404
405         int nsamples = 0;
406         vector<AliCaloBunchInfo> bunchlist; 
407         while (in.NextBunch()) {
408           nsamples += in.GetBunchLength();
409           bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
410         } 
411         
412         if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
413           Float_t time = 0; 
414           Float_t amp = 0; 
415           // indices for pedestal calc.
416           int firstPedSample = 0;
417           int lastPedSample = 0;
418
419           if (! in.IsTRUData() ) { // high gain, low gain, LED Mon data - all have the same shaper/sampling 
420             AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
421             amp = fitResults.GetAmp();
422             time = fitResults.GetTof(); 
423             firstPedSample = fFirstPedestalSample;
424             lastPedSample = fLastPedestalSample;
425           }
426           else { // TRU data is special, needs its own analyzer
427             AliCaloFitResults fitResults = fRawAnalyzerTRU->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
428             amp = fitResults.GetAmp();
429             time = fitResults.GetTof(); 
430             firstPedSample = fFirstPedestalSampleTRU;
431             lastPedSample = fLastPedestalSampleTRU;
432           }
433   
434           // pedestal samples
435           int nPed = 0;
436           UShort_t *pedSamples = 0;
437         
438           // select earliest bunch 
439           unsigned int bunchIndex = 0;
440           unsigned int startBin = bunchlist.at(0).GetStartBin();
441           if (bunchlist.size() > 0) {
442             for(unsigned int ui=1; ui < bunchlist.size(); ui++ ) {
443               if (startBin > bunchlist.at(ui).GetStartBin() ) {
444                 startBin = bunchlist.at(ui).GetStartBin();
445                 bunchIndex = ui;
446               }
447             }
448           }
449
450           // check bunch for entries in the pedestal sample range
451           int bunchLength = bunchlist.at(bunchIndex).GetLength(); 
452           const UShort_t *sig = bunchlist.at(bunchIndex).GetData();
453           int timebin = 0;
454           for (int i = 0; i<bunchLength; i++) {
455             timebin = startBin--;
456             if ( firstPedSample<=timebin && timebin<=lastPedSample ) {
457               pedSamples[nPed] = sig[i];
458               nPed++;
459             }       
460           } // i
461
462           // fill histograms
463           if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
464             int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
465           
466             if ( in.IsLowGain() ) { 
467               nTotalSMLG[iSM]++; 
468               GetRawsData(kTowerLG)->Fill(towerId);
469               if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) { 
470                 GetRawsData(kSigLG)->Fill(towerId, amp);
471                 GetRawsData(kTimeLG)->Fill(towerId, time);
472               }
473               if (nPed > 0) {
474                 for (int i=0; i<nPed; i++) {
475                   GetRawsData(kPedLG)->Fill(towerId, pedSamples[i]);
476                 }
477               }
478             } // gain==0
479             else if ( in.IsHighGain() ) {               
480               nTotalSMHG[iSM]++; 
481               GetRawsData(kTowerHG)->Fill(towerId);
482               if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) { 
483                 GetRawsData(kSigHG)->Fill(towerId, amp);
484                 GetRawsData(kTimeHG)->Fill(towerId, time);
485               } 
486               if (nPed > 0) {
487                 for (int i=0; i<nPed; i++) {
488                   GetRawsData(kPedHG)->Fill(towerId, pedSamples[i]);
489                 }
490               }
491             } // gain==1
492           } // low or high gain
493           // TRU
494           else if ( in.IsTRUData() && in.GetColumn()<AliEMCALGeoParams::fgkEMCAL2x2PerTRU) {
495             // for TRU data, the mapping class holds the TRU internal 2x2 number (0..95) in the Column var..
496             int iTRU = iRCU; //TRU0 is from RCU0, TRU1 from RCU1
497             if (iRCU>0 && in.GetBranch()>0) iTRU=2; // TRU2 is from branch B on RCU1
498             int iTRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU 
499               + in.GetColumn();
500             
501             nTotalSMTRU[iSM]++; 
502             if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) { 
503               GetRawsData(kSigTRU)->Fill(iTRU2x2Id, amp);
504               GetRawsData(kTimeTRU)->Fill(iTRU2x2Id, time);
505             }
506             if (nPed > 0) {
507               for (int i=0; i<nPed; i++) {
508                 GetRawsData(kPedTRU)->Fill(iTRU2x2Id, pedSamples[i]);
509               }
510             }
511           }
512           // LED Mon
513           else if ( in.IsLEDMonData() ) {
514             // for LED Mon data, the mapping class holds the gain info in the Row variable
515             // and the Strip number in the Column..
516             int gain = in.GetRow(); 
517             int stripId = iSM*nStripsPerSM + in.GetColumn();
518           
519             if ( gain == 0 ) { 
520               nTotalSMLGLEDMon[iSM]++; 
521               if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) { 
522                 GetRawsData(kSigLGLEDMon)->Fill(stripId, amp);
523                 GetRawsData(kTimeLGLEDMon)->Fill(stripId, time);
524               }
525               if (nPed > 0) {
526                 for (int i=0; i<nPed; i++) {
527                   GetRawsData(kPedLGLEDMon)->Fill(stripId, pedSamples[i]);
528                 }
529               }
530             } // gain==0
531             else if ( gain == 1 ) {             
532               nTotalSMHGLEDMon[iSM]++; 
533               if ( (amp > fMinSignalHGLEDMon) && (amp < fMaxSignalHGLEDMon) ) { 
534                 GetRawsData(kSigHGLEDMon)->Fill(stripId, amp);
535                 GetRawsData(kTimeHGLEDMon)->Fill(stripId, time);
536               }
537               if (nPed > 0) {
538                 for (int i=0; i<nPed; i++) {
539                   GetRawsData(kPedHGLEDMon)->Fill(stripId, pedSamples[i]);
540                 }
541               }
542             } // low or high gain
543           } // LEDMon
544
545         } // SM index OK
546
547       } // nsamples>0 check, some data found for this channel; not only trailer/header
548     }// end while over channel 
549    
550   }//end while over DDL's, of input stream 
551
552   // let's also fill the SM and event counter histograms
553   int nTotalHG = 0;
554   int nTotalLG = 0;
555   int nTotalTRU = 0;
556   int nTotalHGLEDMon = 0;
557   int nTotalLGLEDMon = 0;
558   for (iSM=0; iSM<fSuperModules; iSM++) {  
559     nTotalLG += nTotalSMLG[iSM]; 
560     nTotalHG += nTotalSMHG[iSM]; 
561     nTotalTRU += nTotalSMTRU[iSM]; 
562     nTotalLG += nTotalSMLGLEDMon[iSM]; 
563     nTotalHG += nTotalSMHGLEDMon[iSM]; 
564     GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]); 
565     GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]); 
566     GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]); 
567     GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]); 
568     GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]); 
569   }
570   
571   GetRawsData(kNtotLG)->Fill(nTotalLG);
572   GetRawsData(kNtotHG)->Fill(nTotalHG);
573   GetRawsData(kNtotTRU)->Fill(nTotalTRU);
574   GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
575   GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
576   
577   SetEventSpecie(saveSpecie) ; 
578   // just in case the next rawreader consumer forgets to reset; let's do it here again..
579   rawReader->Reset() ;
580
581   return;
582 }
583
584 //____________________________________________________________________________
585 void AliEMCALQADataMakerRec::MakeDigits()
586 {
587   // makes data from Digits
588
589   GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ; 
590   TIter next(fDigitsArray) ; 
591   AliEMCALDigit * digit ; 
592   while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
593     GetDigitsData(0)->Fill( digit->GetAmp()) ;
594   }  
595   
596 }
597
598 //____________________________________________________________________________
599 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
600 {
601   // makes data from Digit Tree
602   if (fDigitsArray) 
603     fDigitsArray->Clear() ; 
604   else
605     fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; 
606   
607   TBranch * branch = digitTree->GetBranch("EMCAL") ;
608   if ( ! branch ) {
609     AliWarning("EMCAL branch in Digit Tree not found") ; 
610   } else {
611     branch->SetAddress(&fDigitsArray) ;
612     branch->GetEntry(0) ; 
613     MakeDigits() ; 
614   }
615   
616 }
617
618 //____________________________________________________________________________
619 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
620 {
621   // makes data from RecPoints
622   TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
623   if (!emcbranch) { 
624     AliError("can't get the branch with the EMCAL clusters !");
625     return;
626   }
627   
628   TObjArray * emcrecpoints = new TObjArray(100) ;
629   emcbranch->SetAddress(&emcrecpoints);
630   emcbranch->GetEntry(0);
631   
632   GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ; 
633   TIter next(emcrecpoints) ; 
634   AliEMCALRecPoint * rp ; 
635   while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
636     GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
637     GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
638   }
639   emcrecpoints->Delete();
640   delete emcrecpoints;
641   
642 }
643
644 //____________________________________________________________________________ 
645 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
646 {
647   //Detector specific actions at start of cycle
648   
649 }
650
651 //____________________________________________________________________________ 
652 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)              
653 {
654   //Set fitting algorithm and initialize it if this same algorithm was not set before.
655   //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
656   
657   if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
658     //Do nothing, this same algorithm already set before.
659     //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
660     return;
661   }
662   //Initialize the requested algorithm
663   if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
664     //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
665                 
666     fFittingAlgorithm = fitAlgo; 
667     if (fRawAnalyzer) delete fRawAnalyzer;  // delete prev. analyzer if existed.
668                 
669     if (fitAlgo == kFastFit) {
670       fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
671     }
672     else if (fitAlgo == kNeuralNet) {
673       fRawAnalyzer = new AliCaloRawAnalyzerNN();
674     }
675     else if (fitAlgo == kLMS) {
676       fRawAnalyzer = new AliCaloRawAnalyzerLMS();
677     }
678     else if (fitAlgo == kPeakFinder) {
679       fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
680     }
681     else if (fitAlgo == kCrude) {
682       fRawAnalyzer = new AliCaloRawAnalyzerCrude();
683     }
684     else {
685       AliWarning("EMCAL QA invalid fit algorithm choice") ; 
686     }
687
688   }
689   return;
690 }
691