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