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;-- Yaxian Mao: 11/2009
22 The idea:average counts for all the towers should be flat
23 Change all existing histograms as experts
25 Change histograms for DQM shifter: -- Yaxian Mao 04/2010
26 Calcuate the amplitude ratio from current run and the LED reference, for QAChecker use
27 Also calculate the ratio of amplitude from LED Monitor system (current/Reference), to check LED system
31 // --- ROOT system ---
32 #include <TClonesArray.h>
41 // --- Standard library ---
44 // --- AliRoot header files ---
46 #include "AliESDCaloCluster.h"
47 #include "AliESDCaloCells.h"
48 #include "AliESDEvent.h"
50 #include "AliEMCALQADataMakerRec.h"
51 #include "AliQAChecker.h"
52 #include "AliEMCALDigit.h"
53 #include "AliEMCALRecPoint.h"
54 #include "AliEMCALRawUtils.h"
55 #include "AliEMCALReconstructor.h"
56 #include "AliEMCALRecParam.h"
57 #include "AliRawReader.h"
58 #include "AliCaloRawStreamV3.h"
59 #include "AliEMCALGeoParams.h"
60 #include "AliRawEventHeaderBase.h"
61 #include "AliQAManager.h"
62 #include "AliCDBEntry.h"
64 #include "AliCaloBunchInfo.h"
65 #include "AliCaloFitResults.h"
66 #include "AliCaloRawAnalyzerFastFit.h"
67 #include "AliCaloRawAnalyzerNN.h"
68 #include "AliCaloRawAnalyzerLMS.h"
69 #include "AliCaloRawAnalyzerPeakFinder.h"
70 #include "AliCaloRawAnalyzerCrude.h"
72 #include "AliCaloRawAnalyzerFactory.h"
76 ClassImp(AliEMCALQADataMakerRec)
78 //____________________________________________________________________________
79 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) :
80 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
84 fSuperModules(10), // FIXME!!! number of SuperModules; 10 for 2011; update default for later runs
85 fFirstPedestalSample(0),
86 fLastPedestalSample(3),
87 fFirstPedestalSampleTRU(0),
88 fLastPedestalSampleTRU(3),
90 fMaxSignalLG(AliEMCALGeoParams::fgkSampleMax),
92 fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax),
94 fMaxSignalTRU(AliEMCALGeoParams::fgkSampleMax),
95 fMinSignalLGLEDMon(0),
96 fMaxSignalLGLEDMon(AliEMCALGeoParams::fgkSampleMax),
97 fMinSignalHGLEDMon(0),
98 fMaxSignalHGLEDMon(AliEMCALGeoParams::fgkSampleMax),
99 fCalibRefHistoPro(NULL),
100 fCalibRefHistoH2F(NULL),
101 fLEDMonRefHistoPro(NULL),
102 fHighEmcHistoH2F(NULL)
103 // fTextSM(new TText*[fSuperModules]) ,
109 SetFittingAlgorithm(fitAlgo);
111 //fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
113 fRawAnalyzerTRU = ( AliCaloRawAnalyzerLMS*)AliCaloRawAnalyzerFactory::CreateAnalyzer(kLMS);
115 fRawAnalyzerTRU->SetFixTau(kTRUE);
116 fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
117 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
118 // fTextSM[sm] = NULL ;
122 //____________________________________________________________________________
123 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
125 fFittingAlgorithm(0),
128 fSuperModules(qadm.GetSuperModules()),
129 fFirstPedestalSample(qadm.GetFirstPedestalSample()),
130 fLastPedestalSample(qadm.GetLastPedestalSample()),
131 fFirstPedestalSampleTRU(qadm.GetFirstPedestalSampleTRU()),
132 fLastPedestalSampleTRU(qadm.GetLastPedestalSampleTRU()),
133 fMinSignalLG(qadm.GetMinSignalLG()),
134 fMaxSignalLG(qadm.GetMaxSignalLG()),
135 fMinSignalHG(qadm.GetMinSignalHG()),
136 fMaxSignalHG(qadm.GetMaxSignalHG()),
137 fMinSignalTRU(qadm.GetMinSignalTRU()),
138 fMaxSignalTRU(qadm.GetMaxSignalTRU()),
139 fMinSignalLGLEDMon(qadm.GetMinSignalLGLEDMon()),
140 fMaxSignalLGLEDMon(qadm.GetMaxSignalLGLEDMon()),
141 fMinSignalHGLEDMon(qadm.GetMinSignalHGLEDMon()),
142 fMaxSignalHGLEDMon(qadm.GetMaxSignalHGLEDMon()),
143 fCalibRefHistoPro(NULL),
144 fCalibRefHistoH2F(NULL),
145 fLEDMonRefHistoPro(NULL),
146 fHighEmcHistoH2F(NULL)
147 // fTextSM(new TText*[fSuperModules]) ,
152 SetName((const char*)qadm.GetName()) ;
153 SetTitle((const char*)qadm.GetTitle());
154 SetFittingAlgorithm(qadm.GetFittingAlgorithm());
156 //fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
157 fRawAnalyzerTRU = (AliCaloRawAnalyzerLMS*)AliCaloRawAnalyzerFactory::CreateAnalyzer(kLMS);
158 fRawAnalyzerTRU->SetFixTau(kTRUE);
159 fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
160 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
161 // fTextSM[sm] = qadm.fTextSM[sm] ;
165 //__________________________________________________________________
166 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
169 this->~AliEMCALQADataMakerRec();
170 new(this) AliEMCALQADataMakerRec(qadm);
173 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
174 // fTextSM[sm] = qadm.fTextSM[sm] ;
179 //____________________________________________________________________________
180 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
182 //Detector specific actions at end of cycle
185 // GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
187 // do the QA checking
188 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
191 //____________________________________________________________________________
192 void AliEMCALQADataMakerRec::GetCalibRefFromOCDB()
194 //Get the reference histogram from OCDB
195 TString sName1("hHighEmcalRawMaxMinusMin") ;
196 TString sName2("hLowLEDMonEmcalRawMaxMinusMin") ;
197 sName1.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ;
198 sName2.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ;
200 TString refStorage(AliQAv1::GetQARefStorage()) ;
201 if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
202 AliFatal(Form("%s is not a valid location for reference data", refStorage.Data())) ;
204 AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::kRAWS) ;
205 AliQAv1::SetQARefDataDirName(AliRecoParam::kCalib) ;
206 if ( ! manQA->GetLock() ) {
207 manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
208 manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
209 manQA->SetRun(AliCDBManager::Instance()->GetRun()) ;
212 char * detOCDBDir = Form("%s/%s/%s", GetName(), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
213 AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
215 TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
216 if ( strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
217 AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), GetName())) ;
220 TObjArray * dirOCDB= NULL ;
222 dirOCDB = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib)))) ;
224 fCalibRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName1.Data())) ;
225 fLEDMonRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName2.Data())) ;
230 if(fCalibRefHistoPro && fLEDMonRefHistoPro){
232 //Defining histograms binning, each 2D histogram covers all SMs
233 Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
234 Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
235 Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
237 if(!fCalibRefHistoH2F)
238 fCalibRefHistoH2F = new TH2F("hCalibRefHisto", "hCalibRefHisto", nbinsZ, -0.5, nbinsZ - 0.5, nbinsPhi, -0.5, nbinsPhi -0.5);
239 ConvertProfile2H(fCalibRefHistoPro,fCalibRefHistoH2F) ;
241 AliFatal(Form("No reference object with name %s or %s found", sName1.Data(), sName2.Data())) ;
244 //____________________________________________________________________________
245 void AliEMCALQADataMakerRec::InitESDs()
247 //Create histograms to controll ESD
248 const Bool_t expert = kTRUE ;
249 const Bool_t image = kTRUE ;
251 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts", 200, 0., 100.) ;
253 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
255 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
257 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
259 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts", 500, 0., 50.) ;
261 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
263 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
265 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
269 //____________________________________________________________________________
270 void AliEMCALQADataMakerRec::InitDigits()
272 // create Digits histograms in Digits subdir
273 const Bool_t expert = kTRUE ;
274 const Bool_t image = kTRUE ;
276 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
278 Add2DigitsList(h0, 0, !expert, image) ;
279 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
281 Add2DigitsList(h1, 1, !expert, image) ;
284 //____________________________________________________________________________
285 void AliEMCALQADataMakerRec::InitRecPoints()
287 // create Reconstructed PoInt_ts histograms in RecPoints subdir
288 const Bool_t expert = kTRUE ;
289 const Bool_t image = kTRUE ;
291 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
293 Add2RecPointsList(h0,kRecPE, !expert, image);
295 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
297 Add2RecPointsList(h1,kRecPM, !expert, image);
299 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
301 Add2RecPointsList(h2,kRecPDigM, !expert, image);
305 //____________________________________________________________________________
306 void AliEMCALQADataMakerRec::InitRaws()
308 // create Raws histograms in Raws subdir
309 const Bool_t expert = kTRUE ;
310 const Bool_t saveCorr = kTRUE ;
311 const Bool_t image = kTRUE ;
312 const Option_t *profileOption = "s";
314 Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
315 Int_t nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
317 //Defining histograms binning, each 2D histogram covers all SMs
318 Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
319 Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
320 Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
322 // counter info: number of channels per event (bins are SM index)
323 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
324 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
325 Add2RawsList(h0, kNsmodLG, expert, !image, !saveCorr) ;
326 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
327 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
328 Add2RawsList(h1, kNsmodHG, expert, !image, !saveCorr) ;
330 // where did max sample occur? (bins are towers)
331 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
332 nTot, -0.5, nTot-0.5, profileOption) ;
333 Add2RawsList(h2, kTimeLG, expert, !image, !saveCorr) ;
334 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
335 nTot, -0.5, nTot-0.5, profileOption) ;
336 Add2RawsList(h3, kTimeHG, expert, !image, !saveCorr) ;
338 // how much above pedestal was the max sample? (bins are towers)
339 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
340 nTot, -0.5, nTot-0.5, profileOption) ;
341 Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
342 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
343 nTot, -0.5, nTot-0.5, profileOption) ;
344 Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
346 // total counter: channels per event
347 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
349 Add2RawsList(h6, kNtotLG, expert, !image, !saveCorr) ;
350 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
352 Add2RawsList(h7, kNtotHG, expert, !image, !saveCorr) ;
354 // pedestal (bins are towers)
355 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
356 nTot, -0.5, nTot-0.5, profileOption) ;
357 Add2RawsList(h8, kPedLG, expert, !image, !saveCorr) ;
358 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
359 nTot, -0.5, nTot-0.5, profileOption) ;
360 Add2RawsList(h9, kPedHG, expert, !image, !saveCorr) ;
362 //temp 2D amplitude histogram for the current run
363 fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
364 fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
365 //add ratio histograms: to comapre the current run with the reference data
366 TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5,
367 nbinsPhi, -0.5, nbinsPhi-0.5);
368 //settings for display in amore
369 h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}");
370 h15->SetMaximum(2.0);
371 h15->SetMinimum(0.1);
372 h15->SetOption("COLZ");
373 gStyle->SetOptStat(0);
374 Int_t color[] = {4,3,2} ;
375 gStyle->SetPalette(3,color);
376 h15->GetZaxis()->SetNdivisions(3);
377 h15->UseCurrentStyle();
378 h15->SetDirectory(0);
379 Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
381 TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
382 h16->SetMinimum(0.1);
383 h16->SetMaximum(100.);
384 gStyle->SetOptStat(0);
385 h16->UseCurrentStyle();
386 h16->SetDirectory(0);
387 Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
389 // now repeat the same for TRU and LEDMon data
390 Int_t nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
392 // counter info: number of channels per event (bins are SM index)
393 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
394 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
395 Add2RawsList(hT0, kNsmodTRU, expert, !image, !saveCorr) ;
397 // where did max sample occur? (bins are TRU channels)
398 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
399 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
400 Add2RawsList(hT1, kTimeTRU, expert, !image, !saveCorr) ;
402 // how much above pedestal was the max sample? (bins are TRU channels)
403 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
404 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
405 Add2RawsList(hT2, kSigTRU, expert, !image, !saveCorr) ;
407 // total counter: channels per event
408 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
410 Add2RawsList(hT3, kNtotTRU, expert, !image, !saveCorr) ;
412 // pedestal (bins are TRU channels)
413 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
414 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
415 Add2RawsList(hT4, kPedTRU, expert, !image, !saveCorr) ;
417 // L0 trigger hits: # of hits (bins are TRU channels)
418 TH1I * hT5 = new TH1I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated", nTot2x2, -0.5, nTot2x2);
420 Add2RawsList(hT5, kNL0TRU, expert, !image, !saveCorr);
422 // L0 trigger hits: average time (bins are TRU channels)
423 TProfile * hT6 = new TProfile("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTot2x2, -0.5, nTot2x2, profileOption);
424 Add2RawsList(hT6, kTimeL0TRU, expert, !image, !saveCorr);
426 // and also LED Mon..
427 // LEDMon has both high and low gain channels, just as regular FEE/towers
428 Int_t nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
430 // counter info: number of channels per event (bins are SM index)
431 TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
432 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
433 Add2RawsList(hL0, kNsmodLGLEDMon, expert, !image, !saveCorr) ;
434 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
435 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
436 Add2RawsList(hL1, kNsmodHGLEDMon, expert, !image, !saveCorr) ;
438 // where did max sample occur? (bins are strips)
439 TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
440 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
441 Add2RawsList(hL2, kTimeLGLEDMon, expert, !image, !saveCorr) ;
442 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
443 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
444 Add2RawsList(hL3, kTimeHGLEDMon, expert, !image, !saveCorr) ;
446 // how much above pedestal was the max sample? (bins are strips)
447 TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
448 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
449 Add2RawsList(hL4, kSigLGLEDMon, expert, !image, !saveCorr) ;
450 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
451 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
452 Add2RawsList(hL5, kSigHGLEDMon, expert, !image, !saveCorr) ;
454 // total counter: channels per event
455 TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
457 Add2RawsList(hL6, kNtotLGLEDMon, expert, !image, !saveCorr) ;
458 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
460 Add2RawsList(hL7, kNtotHGLEDMon, expert, !image, !saveCorr) ;
462 // pedestal (bins are strips)
463 TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
464 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
465 Add2RawsList(hL8, kPedLGLEDMon, expert, !image, !saveCorr) ;
466 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
467 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
468 Add2RawsList(hL9, kPedHGLEDMon, expert, !image, !saveCorr) ;
470 //add two histograms for shifter from the LED monitor system: comapre LED monitor with the reference run
471 //to be used for decision whether we need to change reference data
472 TH1F * hL10 = new TH1F("hMaxMinusMinLEDMonRatio", "LEDMon amplitude, Ratio to reference run", nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
473 //settings for display in amore
474 hL10->SetTitle("Amplitude_{LEDMon current}/Amplitude_{LEDMon reference}");
475 hL10->SetMaximum(2.0);
476 hL10->SetMinimum(0.1);
477 gStyle->SetOptStat(0);
478 hL10->UseCurrentStyle();
479 hL10->SetDirectory(0);
480 // hL10->SetOption("E");
481 Add2RawsList(hL10, kLEDMonRatio, !expert, image, !saveCorr) ;
483 TH1F * hL11 = new TH1F("hMaxMinusMinLEDMonRatioDist", "LEDMon amplitude, Ratio distribution", nTotLEDMon, 0, 2);
484 hL11->SetMinimum(0.1) ;
485 gStyle->SetOptStat(0);
486 hL11->UseCurrentStyle();
487 hL11->SetDirectory(0);
488 Add2RawsList(hL11, kLEDMonRatioDist, !expert, image, !saveCorr) ;
490 GetCalibRefFromOCDB();
493 //____________________________________________________________________________
494 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
496 // make QA data from ESDs
499 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
500 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
501 if( clu->IsEMCAL() ) {
502 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
506 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
509 AliESDCaloCells* cells = esd->GetEMCALCells();
510 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
512 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
513 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
518 //____________________________________________________________________________
519 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
521 // Check that all the reference histograms exist before we try to use them - otherwise call InitRaws
522 if (!fCalibRefHistoPro || !fCalibRefHistoH2F || !fLEDMonRefHistoPro || !fHighEmcHistoH2F) {
526 // make sure EMCal was readout during the event
527 Int_t emcID = AliDAQ::DetectorID("EMCAL"); // bit 18..
528 const UInt_t *detPattern = rawReader->GetDetectorPattern();
529 UInt_t emcInReadout = ( ((1 << emcID) & detPattern[0]) >> emcID);
530 if (! emcInReadout) return; // no poInt_t in looking at this event, if no EMCal data
534 AliCaloRawStreamV3 in(rawReader,"EMCAL");
535 rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's
537 AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
539 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
540 SetEventSpecie(AliRecoParam::kCalib) ;
543 const Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
544 const Int_t nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
545 const Int_t nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
546 const Int_t n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
547 const Int_t n2x2PerTRU = AliEMCALGeoParams::fgkEMCAL2x2PerTRU;
549 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
550 Int_t nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
551 Int_t nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
552 Int_t nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
553 Int_t nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
554 Int_t nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
556 const Int_t nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
557 Int_t iSM = 0; // SuperModule index
558 // start loop over input stream
559 while (in.NextDDL()) {
560 Int_t iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
561 fRawAnalyzer->SetIsZeroSuppressed( in.GetZeroSupp() );
563 while (in.NextChannel()) {
564 iSM = in.GetModule(); // SuperModule
565 //prInt_tf("iSM %d DDL %d", iSM, in.GetDDLNumber());
566 if (iSM>=0 && iSM<fSuperModules) { // valid module reading
569 vector<AliCaloBunchInfo> bunchlist;
570 while (in.NextBunch()) {
571 nsamples += in.GetBunchLength();
572 bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
575 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
578 // indices for pedestal calc.
579 Int_t firstPedSample = 0;
580 Int_t lastPedSample = 0;
581 bool isTRUL0IdData = false;
583 if (! in.IsTRUData() ) { // high gain, low gain, LED Mon data - all have the same shaper/sampling
584 AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
585 amp = fitResults.GetAmp();
586 time = fitResults.GetTof();
587 firstPedSample = fFirstPedestalSample;
588 lastPedSample = fLastPedestalSample;
590 else { // TRU data is special, needs its own analyzer
591 AliCaloFitResults fitResults = fRawAnalyzerTRU->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
592 amp = fitResults.GetAmp();
593 time = fitResults.GetTof();
594 firstPedSample = fFirstPedestalSampleTRU;
595 lastPedSample = fLastPedestalSampleTRU;
596 if (in.GetColumn() > n2x2PerTRU) {
597 isTRUL0IdData = true;
603 vector<Int_t> pedSamples;
605 // select earliest bunch
606 unsigned int bunchIndex = 0;
607 unsigned int startBin = bunchlist.at(0).GetStartBin();
608 if (bunchlist.size() > 0) {
609 for(unsigned int ui=1; ui < bunchlist.size(); ui++ ) {
610 if (startBin > bunchlist.at(ui).GetStartBin() ) {
611 startBin = bunchlist.at(ui).GetStartBin();
617 // check bunch for entries in the pedestal sample range
618 Int_t bunchLength = bunchlist.at(bunchIndex).GetLength();
619 const UShort_t *sig = bunchlist.at(bunchIndex).GetData();
622 if (! isTRUL0IdData) { // regular data, can look at pedestals
623 for (Int_t i = 0; i<bunchLength; i++) {
624 timebin = startBin--;
625 if ( firstPedSample<=timebin && timebin<=lastPedSample ) {
626 pedSamples.push_back( sig[i] );
631 else { // TRU L0 Id Data
632 // which TRU the channel belongs to?
633 Int_t iTRUId = in.GetModule()*3 + (iRCU*in.GetBranch() + iRCU);
635 for (Int_t i = 0; i< bunchLength; i++) {
636 for( Int_t j = 0; j < nTRUL0ChannelBits; j++ ){
637 // check if the bit j is 1
638 if( (sig[i] & ( 1 << j )) > 0 ){
639 Int_t iTRUIdInSM = (in.GetColumn() - n2x2PerTRU)*nTRUL0ChannelBits+j;
640 if(iTRUIdInSM < n2x2PerTRU) {
641 Int_t iTRUAbsId = iTRUIdInSM + n2x2PerTRU * iTRUId;
642 // Fill the histograms
643 GetRawsData(kNL0TRU)->Fill(iTRUAbsId);
644 GetRawsData(kTimeL0TRU)->Fill(iTRUAbsId, startBin);
653 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
654 Int_t towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
655 if ( in.IsLowGain() ) {
657 if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) {
658 GetRawsData(kSigLG)->Fill(towerId, amp);
659 GetRawsData(kTimeLG)->Fill(towerId, time);
662 for (Int_t i=0; i<nPed; i++) {
663 GetRawsData(kPedLG)->Fill(towerId, pedSamples[i]);
667 else if ( in.IsHighGain() ) {
669 if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) {
670 GetRawsData(kSigHG)->Fill(towerId, amp);
671 GetRawsData(kTimeHG)->Fill(towerId, time);
674 for (Int_t i=0; i<nPed; i++) {
675 GetRawsData(kPedHG)->Fill(towerId, pedSamples[i]);
679 } // low or high gain
681 else if ( in.IsTRUData() && in.GetColumn()<AliEMCALGeoParams::fgkEMCAL2x2PerTRU) {
682 // for TRU data, the mapping class holds the TRU Int_ternal 2x2 number (0..95) in the Column var..
683 Int_t iTRU = (iRCU*in.GetBranch() + iRCU); //TRU0 is from RCU0, TRU1 from RCU1, TRU2 is from branch B on RCU1
684 Int_t iTRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU
687 if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) {
688 GetRawsData(kSigTRU)->Fill(iTRU2x2Id, amp);
689 GetRawsData(kTimeTRU)->Fill(iTRU2x2Id, time);
692 for (Int_t i=0; i<nPed; i++) {
693 GetRawsData(kPedTRU)->Fill(iTRU2x2Id, pedSamples[i]);
698 else if ( in.IsLEDMonData() ) {
699 // for LED Mon data, the mapping class holds the gain info in the Row variable
700 // and the Strip number in the Column..
701 Int_t gain = in.GetRow();
702 Int_t stripId = iSM*nStripsPerSM + in.GetColumn();
705 nTotalSMLGLEDMon[iSM]++;
706 if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) {
707 GetRawsData(kSigLGLEDMon)->Fill(stripId, amp);
708 GetRawsData(kTimeLGLEDMon)->Fill(stripId, time);
711 for (Int_t i=0; i<nPed; i++) {
712 GetRawsData(kPedLGLEDMon)->Fill(stripId, pedSamples[i]);
716 else if ( gain == 1 ) {
717 nTotalSMHGLEDMon[iSM]++;
718 if ( (amp > fMinSignalHGLEDMon) && (amp < fMaxSignalHGLEDMon) ) {
719 GetRawsData(kSigHGLEDMon)->Fill(stripId, amp);
720 GetRawsData(kTimeHGLEDMon)->Fill(stripId, time);
723 for (Int_t i=0; i<nPed; i++) {
724 GetRawsData(kPedHGLEDMon)->Fill(stripId, pedSamples[i]);
727 } // low or high gain
732 } // nsamples>0 check, some data found for this channel; not only trailer/header
733 }// end while over channel
735 }//end while over DDL's, of input stream
737 // TProfile * p = dynamic_cast<TProfile *>(GetRawsData(kSigHG)) ;
738 ConvertProfile2H(dynamic_cast<TProfile *>(GetRawsData(kSigHG)), fHighEmcHistoH2F) ;
739 Double_t binContent = 0. ;
742 //calculate the ratio of the amplitude and fill the histograms, only if the events type is Calib
743 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
744 //reset ratio histograms
745 GetRawsData(k2DRatioAmp)->Reset("ICE");
746 GetRawsData(kRatioDist)->Reset("ICE");
747 GetRawsData(kLEDMonRatio)->Reset("ICE");
748 GetRawsData(kLEDMonRatioDist)->Reset("ICE");
749 GetRawsData(k2DRatioAmp)->ResetStats();
750 GetRawsData(kRatioDist)->ResetStats();
751 GetRawsData(kLEDMonRatio)->ResetStats();
752 GetRawsData(kLEDMonRatioDist)->ResetStats();
754 for(Int_t ix = 1; ix <= fHighEmcHistoH2F->GetNbinsX(); ix++) {
755 for(Int_t iy = 1; iy <= fHighEmcHistoH2F->GetNbinsY(); iy++) {
756 if(fCalibRefHistoH2F->GetBinContent(ix, iy))binContent = fHighEmcHistoH2F->GetBinContent(ix, iy)/fCalibRefHistoH2F->GetBinContent(ix, iy) ;
757 GetRawsData(k2DRatioAmp)->SetBinContent(ix, iy, binContent);
758 GetRawsData(kRatioDist)->Fill(GetRawsData(k2DRatioAmp)->GetBinContent(ix, iy));
762 //Now for LED monitor system, to calculate the ratio as well
763 Double_t binError = 0. ;
764 // for the binError, we add the relative errors, squared
765 Double_t relativeErrorSqr = 0. ;
767 for(int ib = 1; ib <= fLEDMonRefHistoPro->GetNbinsX(); ib++) {
769 if(fLEDMonRefHistoPro->GetBinContent(ib) != 0) {
770 binContent = GetRawsData(kSigLGLEDMon)->GetBinContent(ib) / fLEDMonRefHistoPro->GetBinContent(ib);
772 relativeErrorSqr = TMath::Power( (fLEDMonRefHistoPro->GetBinError(ib) / fLEDMonRefHistoPro->GetBinContent(ib)), 2);
773 if(GetRawsData(kSigLGLEDMon)->GetBinContent(ib) != 0) {
774 relativeErrorSqr += TMath::Power( (GetRawsData(kSigLGLEDMon)->GetBinError(ib)/GetRawsData(kSigLGLEDMon)->GetBinContent(ib)), 2);
779 relativeErrorSqr = 0;
781 GetRawsData(kLEDMonRatio)->SetBinContent(ib, binContent);
783 binError = sqrt(relativeErrorSqr) * binContent;
784 GetRawsData(kLEDMonRatio)->SetBinError(ib, binError);
785 GetRawsData(kLEDMonRatioDist)->Fill(GetRawsData(kLEDMonRatio)->GetBinContent(ib));
789 // let's also fill the SM and event counter histograms
793 Int_t nTotalHGLEDMon = 0;
794 Int_t nTotalLGLEDMon = 0;
795 for (iSM=0; iSM<fSuperModules; iSM++) {
796 nTotalLG += nTotalSMLG[iSM];
797 nTotalHG += nTotalSMHG[iSM];
798 nTotalTRU += nTotalSMTRU[iSM];
799 nTotalLGLEDMon += nTotalSMLGLEDMon[iSM];
800 nTotalHGLEDMon += nTotalSMHGLEDMon[iSM];
801 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
802 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
803 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
804 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
805 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
808 GetRawsData(kNtotLG)->Fill(nTotalLG);
809 GetRawsData(kNtotHG)->Fill(nTotalHG);
810 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
811 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
812 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
815 SetEventSpecie(saveSpecie) ;
816 // just in case the next rawreader consumer forgets to reset; let's do it here again..
822 //____________________________________________________________________________
823 void AliEMCALQADataMakerRec::MakeDigits()
825 // makes data from Digits
827 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
828 TIter next(fDigitsArray) ;
829 AliEMCALDigit * digit ;
830 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
831 GetDigitsData(0)->Fill( digit->GetAmplitude()) ;
836 //____________________________________________________________________________
837 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
839 // makes data from Digit Tree
841 fDigitsArray->Clear("C") ;
843 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
845 TBranch * branch = digitTree->GetBranch("EMCAL") ;
847 AliWarning("EMCAL branch in Digit Tree not found") ;
849 branch->SetAddress(&fDigitsArray) ;
850 branch->GetEntry(0) ;
856 //____________________________________________________________________________
857 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
859 // makes data from RecPoints
860 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
862 AliError("can't get the branch with the EMCAL clusters !");
866 TObjArray * emcRecPoints = new TObjArray(100) ;
867 emcbranch->SetAddress(&emcRecPoints);
868 emcbranch->GetEntry(0);
870 GetRecPointsData(kRecPM)->Fill(emcRecPoints->GetEntriesFast()) ;
871 TIter next(emcRecPoints) ;
872 AliEMCALRecPoint * rp ;
873 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
874 GetRecPointsData(kRecPE)->Fill(rp->GetEnergy()) ;
875 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
877 emcRecPoints->Delete();
882 //____________________________________________________________________________
883 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
885 //Detector specific actions at start of cycle
889 //____________________________________________________________________________
890 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)
892 //Set fitting algorithm and initialize it if this same algorithm was not set before.
893 //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
896 // fRawAnalyzer = AliCaloRawAnalyzerFactory::CreateAnalyzer(fitAlgo);
897 // fFittingAlgorithm = fitAlgo;
901 fRawAnalyzer = AliCaloRawAnalyzerFactory::CreateAnalyzer(kLMS);
902 fFittingAlgorithm = kLMS;
906 if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
907 //Do nothing, this same algorithm already set before.
908 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
911 //Initialize the requested algorithm
912 if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
913 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
915 fFittingAlgorithm = fitAlgo;
916 if (fRawAnalyzer) delete fRawAnalyzer; // delete prev. analyzer if existed.
918 if (fitAlgo == kFastFit) {
919 fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
921 else if (fitAlgo == kNeuralNet) {
922 fRawAnalyzer = new AliCaloRawAnalyzerNN();
924 else if (fitAlgo == kLMS) {
925 fRawAnalyzer = new AliCaloRawAnalyzerLMS();
927 else if (fitAlgo == kPeakFinder) {
928 fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
930 else if (fitAlgo == kCrude) {
931 fRawAnalyzer = new AliCaloRawAnalyzerCrude();
934 AliWarning("EMCAL QA invalid fit algorithm choice") ;
942 //_____________________________________________________________________________________
943 void AliEMCALQADataMakerRec::ConvertProfile2H(TProfile * p, TH2 * histo)
946 histo->Reset("ICE") ;
949 Int_t nbinsProf = p->GetNbinsX();
951 // loop through the TProfile p and fill the TH2F histo
954 Double_t binContent = 0;
955 Int_t towerNum = 0; // global tower Id
956 // i = 0; // tower Id within SuperModule
957 Int_t iSM = 0; // SuperModule index
958 Int_t iSMSide = 0; // 0=A, 1=C side
959 Int_t iSMSector = 0; // 2 SM's per sector
961 // indices for 2D plots
965 for (Int_t ibin = 1; ibin <= nbinsProf; ibin++) {
966 towerNum = (Int_t) p->GetBinCenter(ibin);
967 binContent = p->GetBinContent(ibin);
969 // figure out what the tower indices are: col, row within a SuperModule
970 iSM = towerNum/(AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols);
971 col = (towerNum/AliEMCALGeoParams::fgkEMCALRows) % (AliEMCALGeoParams::fgkEMCALCols);
972 row = towerNum % (AliEMCALGeoParams::fgkEMCALRows);
974 //DecodeTowerNum(towerNum, &SM, &col, &row);
975 // then we calculate what the global 2D coord are, based on which SM
980 if (iSMSide == 1) { // C side, shown to the right
981 col2d = col + AliEMCALGeoParams::fgkEMCALCols;
983 else { // A side, shown to the left
987 row2d = row + iSMSector * AliEMCALGeoParams::fgkEMCALRows;
989 histo->SetBinContent(col2d+1, row2d+1, binContent);