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