1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 Based on the QA code for PHOS written by Yves Schutz July 2007
18 Authors: J.Klay (Cal Poly) May 2008
19 S. Salur LBL April 2008
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
28 // --- ROOT system ---
29 #include <TClonesArray.h>
36 // --- Standard library ---
39 // --- AliRoot header files ---
41 #include "AliESDCaloCluster.h"
42 #include "AliESDCaloCells.h"
43 #include "AliESDEvent.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"
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"
67 ClassImp(AliEMCALQADataMakerRec)
69 //____________________________________________________________________________
70 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) :
71 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
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),
81 fMaxSignalLG(AliEMCALGeoParams::fgkSampleMax),
83 fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax),
85 fMaxSignalTRU(AliEMCALGeoParams::fgkSampleMax),
86 fMinSignalLGLEDMon(0),
87 fMaxSignalLGLEDMon(AliEMCALGeoParams::fgkSampleMax),
88 fMinSignalHGLEDMon(0),
89 fMaxSignalHGLEDMon(AliEMCALGeoParams::fgkSampleMax)
92 SetFittingAlgorithm(fitAlgo);
93 fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
94 fRawAnalyzerTRU->SetFixTau(kTRUE);
95 fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
98 //____________________________________________________________________________
99 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
101 fFittingAlgorithm(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())
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
129 //__________________________________________________________________
130 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
133 this->~AliEMCALQADataMakerRec();
134 new(this) AliEMCALQADataMakerRec(qadm);
138 //____________________________________________________________________________
139 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
141 //Detector specific actions at end of cycle
144 // GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
146 // do the QA checking
147 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
150 //____________________________________________________________________________
151 void AliEMCALQADataMakerRec::InitESDs()
153 //Create histograms to controll ESD
154 const Bool_t expert = kTRUE ;
155 const Bool_t image = kTRUE ;
157 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts", 200, 0., 100.) ;
159 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
161 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
163 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
165 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts", 500, 0., 50.) ;
167 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
169 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
171 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
175 //____________________________________________________________________________
176 void AliEMCALQADataMakerRec::InitDigits()
178 // create Digits histograms in Digits subdir
179 const Bool_t expert = kTRUE ;
180 const Bool_t image = kTRUE ;
182 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
184 Add2DigitsList(h0, 0, !expert, image) ;
185 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
187 Add2DigitsList(h1, 1, !expert, image) ;
190 //____________________________________________________________________________
191 void AliEMCALQADataMakerRec::InitRecPoints()
193 // create Reconstructed Points histograms in RecPoints subdir
194 const Bool_t expert = kTRUE ;
195 const Bool_t image = kTRUE ;
197 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
199 Add2RecPointsList(h0,kRecPE, !expert, image);
201 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
203 Add2RecPointsList(h1,kRecPM, !expert, image);
205 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
207 Add2RecPointsList(h2,kRecPDigM, !expert, image);
211 //____________________________________________________________________________
212 void AliEMCALQADataMakerRec::InitRaws()
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 ;
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
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) ;
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) ;
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) ;
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) ;
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) ;
252 Add2RawsList(h7, kNtotHG, expert, image, !saveCorr) ;
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) ;
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) ;
265 Add2RawsList(h12, kTowerHG, !expert, image, !saveCorr) ;
266 TH1I * h13 = new TH1I("hTowerLG", "Low Gains on the Tower;Tower", nTot,0, nTot) ;
268 Add2RawsList(h13, kTowerLG, !expert, image, !saveCorr) ;
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
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) ;
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) ;
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) ;
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) ;
291 Add2RawsList(hT3, kNtotTRU, expert, image, !saveCorr) ;
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) ;
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);
301 Add2RawsList(hT5, kNL0TRU, expert, image, !saveCorr);
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);
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
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) ;
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) ;
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) ;
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) ;
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) ;
341 Add2RawsList(hL7, kNtotHGLEDMon, expert, image, !saveCorr) ;
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) ;
353 //____________________________________________________________________________
354 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
356 // make QA data from ESDs
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()) ;
366 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
369 AliESDCaloCells* cells = esd->GetEMCALCells();
370 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
372 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
373 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
378 //____________________________________________________________________________
379 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
381 //Fill prepared histograms with Raw digit properties
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
391 AliCaloRawStreamV3 in(rawReader,"EMCAL");
392 rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's
394 AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
396 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
397 SetEventSpecie(AliRecoParam::kCalib) ;
400 fRawAnalyzer->SetIsZeroSuppressed(true); // TMP - should use stream->IsZeroSuppressed(), or altro cfg registers later
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;
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};
415 int nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
417 // start loop over input stream
419 while (in.NextDDL()) {
420 int iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
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
428 vector<AliCaloBunchInfo> bunchlist;
429 while (in.NextBunch()) {
430 nsamples += in.GetBunchLength();
431 bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
434 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
437 // indices for pedestal calc.
438 int firstPedSample = 0;
439 int lastPedSample = 0;
440 bool isTRUL0IdData = false;
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;
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;
462 vector<int> pedSamples;
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();
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();
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] );
489 // printf("nPed %d\n", nPed);
491 else { // TRU L0 Id Data
492 // which TRU the channel belongs to?
493 int iTRUId = in.GetModule()*3 + (iRCU*in.GetBranch() + iRCU);
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);
513 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
514 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
515 if ( in.IsLowGain() ) {
517 GetRawsData(kTowerLG)->Fill(towerId);
518 if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) {
519 GetRawsData(kSigLG)->Fill(towerId, amp);
520 GetRawsData(kTimeLG)->Fill(towerId, time);
523 for (int i=0; i<nPed; i++) {
524 GetRawsData(kPedLG)->Fill(towerId, pedSamples[i]);
528 else if ( in.IsHighGain() ) {
530 GetRawsData(kTowerHG)->Fill(towerId);
531 if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) {
532 GetRawsData(kSigHG)->Fill(towerId, amp);
533 GetRawsData(kTimeHG)->Fill(towerId, time);
536 for (int i=0; i<nPed; i++) {
537 GetRawsData(kPedHG)->Fill(towerId, pedSamples[i]);
541 } // low or high gain
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
549 if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) {
550 GetRawsData(kSigTRU)->Fill(iTRU2x2Id, amp);
551 GetRawsData(kTimeTRU)->Fill(iTRU2x2Id, time);
554 for (int i=0; i<nPed; i++) {
555 GetRawsData(kPedTRU)->Fill(iTRU2x2Id, pedSamples[i]);
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();
567 nTotalSMLGLEDMon[iSM]++;
568 if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) {
569 GetRawsData(kSigLGLEDMon)->Fill(stripId, amp);
570 GetRawsData(kTimeLGLEDMon)->Fill(stripId, time);
573 for (int i=0; i<nPed; i++) {
574 GetRawsData(kPedLGLEDMon)->Fill(stripId, pedSamples[i]);
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);
585 for (int i=0; i<nPed; i++) {
586 GetRawsData(kPedHGLEDMon)->Fill(stripId, pedSamples[i]);
589 } // low or high gain
594 } // nsamples>0 check, some data found for this channel; not only trailer/header
595 }// end while over channel
597 }//end while over DDL's, of input stream
599 // let's also fill the SM and event counter histograms
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]);
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);
625 SetEventSpecie(saveSpecie) ;
626 // just in case the next rawreader consumer forgets to reset; let's do it here again..
632 //____________________________________________________________________________
633 void AliEMCALQADataMakerRec::MakeDigits()
635 // makes data from Digits
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()) ;
646 //____________________________________________________________________________
647 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
649 // makes data from Digit Tree
651 fDigitsArray->Clear() ;
653 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
655 TBranch * branch = digitTree->GetBranch("EMCAL") ;
657 AliWarning("EMCAL branch in Digit Tree not found") ;
659 branch->SetAddress(&fDigitsArray) ;
660 branch->GetEntry(0) ;
666 //____________________________________________________________________________
667 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
669 // makes data from RecPoints
670 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
672 AliError("can't get the branch with the EMCAL clusters !");
676 TObjArray * emcrecpoints = new TObjArray(100) ;
677 emcbranch->SetAddress(&emcrecpoints);
678 emcbranch->GetEntry(0);
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());
687 emcrecpoints->Delete();
692 //____________________________________________________________________________
693 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
695 //Detector specific actions at start of cycle
699 //____________________________________________________________________________
700 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)
702 //Set fitting algorithm and initialize it if this same algorithm was not set before.
703 //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
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());
710 //Initialize the requested algorithm
711 if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
712 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
714 fFittingAlgorithm = fitAlgo;
715 if (fRawAnalyzer) delete fRawAnalyzer; // delete prev. analyzer if existed.
717 if (fitAlgo == kFastFit) {
718 fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
720 else if (fitAlgo == kNeuralNet) {
721 fRawAnalyzer = new AliCaloRawAnalyzerNN();
723 else if (fitAlgo == kLMS) {
724 fRawAnalyzer = new AliCaloRawAnalyzerLMS();
726 else if (fitAlgo == kPeakFinder) {
727 fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
729 else if (fitAlgo == kCrude) {
730 fRawAnalyzer = new AliCaloRawAnalyzerCrude();
733 AliWarning("EMCAL QA invalid fit algorithm choice") ;