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"
74 ClassImp(AliEMCALQADataMakerRec)
76 //____________________________________________________________________________
77 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) :
78 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
82 fSuperModules(10), // FIXME!!! number of SuperModules; 10 for 2011; update default for later runs
83 fFirstPedestalSample(0),
84 fLastPedestalSample(3),
85 fFirstPedestalSampleTRU(0),
86 fLastPedestalSampleTRU(3),
88 fMaxSignalLG(AliEMCALGeoParams::fgkSampleMax),
90 fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax),
92 fMaxSignalTRU(AliEMCALGeoParams::fgkSampleMax),
93 fMinSignalLGLEDMon(0),
94 fMaxSignalLGLEDMon(AliEMCALGeoParams::fgkSampleMax),
95 fMinSignalHGLEDMon(0),
96 fMaxSignalHGLEDMon(AliEMCALGeoParams::fgkSampleMax),
97 fCalibRefHistoPro(NULL),
98 fCalibRefHistoH2F(NULL),
99 fLEDMonRefHistoPro(NULL),
100 fHighEmcHistoH2F(NULL)
101 // fTextSM(new TText*[fSuperModules]) ,
107 SetFittingAlgorithm(fitAlgo);
108 fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
109 fRawAnalyzerTRU->SetFixTau(kTRUE);
110 fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
111 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
112 // fTextSM[sm] = NULL ;
116 //____________________________________________________________________________
117 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
119 fFittingAlgorithm(0),
122 fSuperModules(qadm.GetSuperModules()),
123 fFirstPedestalSample(qadm.GetFirstPedestalSample()),
124 fLastPedestalSample(qadm.GetLastPedestalSample()),
125 fFirstPedestalSampleTRU(qadm.GetFirstPedestalSampleTRU()),
126 fLastPedestalSampleTRU(qadm.GetLastPedestalSampleTRU()),
127 fMinSignalLG(qadm.GetMinSignalLG()),
128 fMaxSignalLG(qadm.GetMaxSignalLG()),
129 fMinSignalHG(qadm.GetMinSignalHG()),
130 fMaxSignalHG(qadm.GetMaxSignalHG()),
131 fMinSignalTRU(qadm.GetMinSignalTRU()),
132 fMaxSignalTRU(qadm.GetMaxSignalTRU()),
133 fMinSignalLGLEDMon(qadm.GetMinSignalLGLEDMon()),
134 fMaxSignalLGLEDMon(qadm.GetMaxSignalLGLEDMon()),
135 fMinSignalHGLEDMon(qadm.GetMinSignalHGLEDMon()),
136 fMaxSignalHGLEDMon(qadm.GetMaxSignalHGLEDMon()),
137 fCalibRefHistoPro(NULL),
138 fCalibRefHistoH2F(NULL),
139 fLEDMonRefHistoPro(NULL),
140 fHighEmcHistoH2F(NULL)
141 // fTextSM(new TText*[fSuperModules]) ,
146 SetName((const char*)qadm.GetName()) ;
147 SetTitle((const char*)qadm.GetTitle());
148 SetFittingAlgorithm(qadm.GetFittingAlgorithm());
149 fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
150 fRawAnalyzerTRU->SetFixTau(kTRUE);
151 fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
152 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
153 // fTextSM[sm] = qadm.fTextSM[sm] ;
157 //__________________________________________________________________
158 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
161 this->~AliEMCALQADataMakerRec();
162 new(this) AliEMCALQADataMakerRec(qadm);
165 // for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
166 // fTextSM[sm] = qadm.fTextSM[sm] ;
171 //____________________________________________________________________________
172 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
174 //Detector specific actions at end of cycle
177 // GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
179 // do the QA checking
180 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
183 //____________________________________________________________________________
184 void AliEMCALQADataMakerRec::GetCalibRefFromOCDB()
186 //Get the reference histogram from OCDB
187 TString sName1("hHighEmcalRawMaxMinusMin") ;
188 TString sName2("hLowLEDMonEmcalRawMaxMinusMin") ;
189 sName1.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ;
190 sName2.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ;
192 TString refStorage(AliQAv1::GetQARefStorage()) ;
193 if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
194 AliFatal(Form("%s is not a valid location for reference data", refStorage.Data())) ;
196 AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::kRAWS) ;
197 AliQAv1::SetQARefDataDirName(AliRecoParam::kCalib) ;
198 if ( ! manQA->GetLock() ) {
199 manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
200 manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
201 manQA->SetRun(AliCDBManager::Instance()->GetRun()) ;
204 char * detOCDBDir = Form("%s/%s/%s", GetName(), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
205 AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
207 TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
208 if ( strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
209 AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), GetName())) ;
212 TObjArray * dirOCDB= NULL ;
214 dirOCDB = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib)))) ;
216 fCalibRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName1.Data())) ;
217 fLEDMonRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName2.Data())) ;
222 if(fCalibRefHistoPro && fLEDMonRefHistoPro){
224 //Defining histograms binning, each 2D histogram covers all SMs
225 Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
226 Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
227 Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
229 if(!fCalibRefHistoH2F)
230 fCalibRefHistoH2F = new TH2F("hCalibRefHisto", "hCalibRefHisto", nbinsZ, -0.5, nbinsZ - 0.5, nbinsPhi, -0.5, nbinsPhi -0.5);
231 ConvertProfile2H(fCalibRefHistoPro,fCalibRefHistoH2F) ;
233 AliFatal(Form("No reference object with name %s or %s found", sName1.Data(), sName2.Data())) ;
236 //____________________________________________________________________________
237 void AliEMCALQADataMakerRec::InitESDs()
239 //Create histograms to controll ESD
240 const Bool_t expert = kTRUE ;
241 const Bool_t image = kTRUE ;
243 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts", 200, 0., 100.) ;
245 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
247 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
249 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
251 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts", 500, 0., 50.) ;
253 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
255 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
257 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
261 //____________________________________________________________________________
262 void AliEMCALQADataMakerRec::InitDigits()
264 // create Digits histograms in Digits subdir
265 const Bool_t expert = kTRUE ;
266 const Bool_t image = kTRUE ;
268 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
270 Add2DigitsList(h0, 0, !expert, image) ;
271 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
273 Add2DigitsList(h1, 1, !expert, image) ;
276 //____________________________________________________________________________
277 void AliEMCALQADataMakerRec::InitRecPoints()
279 // create Reconstructed PoInt_ts histograms in RecPoints subdir
280 const Bool_t expert = kTRUE ;
281 const Bool_t image = kTRUE ;
283 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
285 Add2RecPointsList(h0,kRecPE, !expert, image);
287 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
289 Add2RecPointsList(h1,kRecPM, !expert, image);
291 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
293 Add2RecPointsList(h2,kRecPDigM, !expert, image);
297 //____________________________________________________________________________
298 void AliEMCALQADataMakerRec::InitRaws()
300 // create Raws histograms in Raws subdir
301 const Bool_t expert = kTRUE ;
302 const Bool_t saveCorr = kTRUE ;
303 const Bool_t image = kTRUE ;
304 const Option_t *profileOption = "s";
306 Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
307 Int_t nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
309 //Defining histograms binning, each 2D histogram covers all SMs
310 Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
311 Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
312 Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
314 // counter info: number of channels per event (bins are SM index)
315 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
316 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
317 Add2RawsList(h0, kNsmodLG, expert, !image, !saveCorr) ;
318 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
319 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
320 Add2RawsList(h1, kNsmodHG, expert, !image, !saveCorr) ;
322 // where did max sample occur? (bins are towers)
323 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
324 nTot, -0.5, nTot-0.5, profileOption) ;
325 Add2RawsList(h2, kTimeLG, expert, !image, !saveCorr) ;
326 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
327 nTot, -0.5, nTot-0.5, profileOption) ;
328 Add2RawsList(h3, kTimeHG, expert, !image, !saveCorr) ;
330 // how much above pedestal was the max sample? (bins are towers)
331 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
332 nTot, -0.5, nTot-0.5, profileOption) ;
333 Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
334 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
335 nTot, -0.5, nTot-0.5, profileOption) ;
336 Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
338 // total counter: channels per event
339 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
341 Add2RawsList(h6, kNtotLG, expert, !image, !saveCorr) ;
342 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
344 Add2RawsList(h7, kNtotHG, expert, !image, !saveCorr) ;
346 // pedestal (bins are towers)
347 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
348 nTot, -0.5, nTot-0.5, profileOption) ;
349 Add2RawsList(h8, kPedLG, expert, !image, !saveCorr) ;
350 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
351 nTot, -0.5, nTot-0.5, profileOption) ;
352 Add2RawsList(h9, kPedHG, expert, !image, !saveCorr) ;
354 //temp 2D amplitude histogram for the current run
355 fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
356 fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
357 //add ratio histograms: to comapre the current run with the reference data
358 TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5,
359 nbinsPhi, -0.5, nbinsPhi-0.5);
360 //settings for display in amore
361 h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}");
362 h15->SetMaximum(2.0);
363 h15->SetMinimum(0.1);
364 h15->SetOption("COLZ");
365 gStyle->SetOptStat(0);
366 Int_t color[] = {4,3,2} ;
367 gStyle->SetPalette(3,color);
368 h15->GetZaxis()->SetNdivisions(3);
369 h15->UseCurrentStyle();
370 h15->SetDirectory(0);
371 Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
373 TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
374 h16->SetMinimum(0.1);
375 h16->SetMaximum(100.);
376 gStyle->SetOptStat(0);
377 h16->UseCurrentStyle();
378 h16->SetDirectory(0);
379 Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
381 // now repeat the same for TRU and LEDMon data
382 Int_t nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
384 // counter info: number of channels per event (bins are SM index)
385 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
386 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
387 Add2RawsList(hT0, kNsmodTRU, expert, !image, !saveCorr) ;
389 // where did max sample occur? (bins are TRU channels)
390 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
391 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
392 Add2RawsList(hT1, kTimeTRU, expert, !image, !saveCorr) ;
394 // how much above pedestal was the max sample? (bins are TRU channels)
395 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
396 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
397 Add2RawsList(hT2, kSigTRU, expert, !image, !saveCorr) ;
399 // total counter: channels per event
400 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
402 Add2RawsList(hT3, kNtotTRU, expert, !image, !saveCorr) ;
404 // pedestal (bins are TRU channels)
405 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
406 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
407 Add2RawsList(hT4, kPedTRU, expert, !image, !saveCorr) ;
409 // L0 trigger hits: # of hits (bins are TRU channels)
410 TH1I * hT5 = new TH1I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated", nTot2x2, -0.5, nTot2x2);
412 Add2RawsList(hT5, kNL0TRU, expert, !image, !saveCorr);
414 // L0 trigger hits: average time (bins are TRU channels)
415 TProfile * hT6 = new TProfile("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTot2x2, -0.5, nTot2x2, profileOption);
416 Add2RawsList(hT6, kTimeL0TRU, expert, !image, !saveCorr);
418 // and also LED Mon..
419 // LEDMon has both high and low gain channels, just as regular FEE/towers
420 Int_t nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
422 // counter info: number of channels per event (bins are SM index)
423 TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
424 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
425 Add2RawsList(hL0, kNsmodLGLEDMon, expert, !image, !saveCorr) ;
426 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
427 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
428 Add2RawsList(hL1, kNsmodHGLEDMon, expert, !image, !saveCorr) ;
430 // where did max sample occur? (bins are strips)
431 TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
432 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
433 Add2RawsList(hL2, kTimeLGLEDMon, expert, !image, !saveCorr) ;
434 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
435 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
436 Add2RawsList(hL3, kTimeHGLEDMon, expert, !image, !saveCorr) ;
438 // how much above pedestal was the max sample? (bins are strips)
439 TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
440 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
441 Add2RawsList(hL4, kSigLGLEDMon, expert, !image, !saveCorr) ;
442 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
443 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
444 Add2RawsList(hL5, kSigHGLEDMon, expert, !image, !saveCorr) ;
446 // total counter: channels per event
447 TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
449 Add2RawsList(hL6, kNtotLGLEDMon, expert, !image, !saveCorr) ;
450 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
452 Add2RawsList(hL7, kNtotHGLEDMon, expert, !image, !saveCorr) ;
454 // pedestal (bins are strips)
455 TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
456 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
457 Add2RawsList(hL8, kPedLGLEDMon, expert, !image, !saveCorr) ;
458 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
459 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
460 Add2RawsList(hL9, kPedHGLEDMon, expert, !image, !saveCorr) ;
462 //add two histograms for shifter from the LED monitor system: comapre LED monitor with the reference run
463 //to be used for decision whether we need to change reference data
464 TH1F * hL10 = new TH1F("hMaxMinusMinLEDMonRatio", "LEDMon amplitude, Ratio to reference run", nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
465 //settings for display in amore
466 hL10->SetTitle("Amplitude_{LEDMon current}/Amplitude_{LEDMon reference}");
467 hL10->SetMaximum(2.0);
468 hL10->SetMinimum(0.1);
469 gStyle->SetOptStat(0);
470 hL10->UseCurrentStyle();
471 hL10->SetDirectory(0);
472 // hL10->SetOption("E");
473 Add2RawsList(hL10, kLEDMonRatio, !expert, image, !saveCorr) ;
475 TH1F * hL11 = new TH1F("hMaxMinusMinLEDMonRatioDist", "LEDMon amplitude, Ratio distribution", nTotLEDMon, 0, 2);
476 hL11->SetMinimum(0.1) ;
477 gStyle->SetOptStat(0);
478 hL11->UseCurrentStyle();
479 hL11->SetDirectory(0);
480 Add2RawsList(hL11, kLEDMonRatioDist, !expert, image, !saveCorr) ;
482 GetCalibRefFromOCDB();
485 //____________________________________________________________________________
486 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
488 // make QA data from ESDs
491 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
492 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
493 if( clu->IsEMCAL() ) {
494 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
498 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
501 AliESDCaloCells* cells = esd->GetEMCALCells();
502 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
504 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
505 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
510 //____________________________________________________________________________
511 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
513 // Check that all the reference histograms exist before we try to use them - otherwise call InitRaws
514 if (!fCalibRefHistoPro || !fCalibRefHistoH2F || !fLEDMonRefHistoPro || !fHighEmcHistoH2F) {
518 // make sure EMCal was readout during the event
519 Int_t emcID = AliDAQ::DetectorID("EMCAL"); // bit 18..
520 const UInt_t *detPattern = rawReader->GetDetectorPattern();
521 UInt_t emcInReadout = ( ((1 << emcID) & detPattern[0]) >> emcID);
522 if (! emcInReadout) return; // no poInt_t in looking at this event, if no EMCal data
526 AliCaloRawStreamV3 in(rawReader,"EMCAL");
527 rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's
529 AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
531 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
532 SetEventSpecie(AliRecoParam::kCalib) ;
535 const Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
536 const Int_t nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
537 const Int_t nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
538 const Int_t n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
539 const Int_t n2x2PerTRU = AliEMCALGeoParams::fgkEMCAL2x2PerTRU;
541 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
542 Int_t nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
543 Int_t nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
544 Int_t nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
545 Int_t nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
546 Int_t nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
548 const Int_t nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
549 Int_t iSM = 0; // SuperModule index
550 // start loop over input stream
551 while (in.NextDDL()) {
552 Int_t iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
553 fRawAnalyzer->SetIsZeroSuppressed( in.GetZeroSupp() );
555 while (in.NextChannel()) {
556 iSM = in.GetModule(); // SuperModule
557 //prInt_tf("iSM %d DDL %d", iSM, in.GetDDLNumber());
558 if (iSM>=0 && iSM<fSuperModules) { // valid module reading
561 vector<AliCaloBunchInfo> bunchlist;
562 while (in.NextBunch()) {
563 nsamples += in.GetBunchLength();
564 bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
567 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
570 // indices for pedestal calc.
571 Int_t firstPedSample = 0;
572 Int_t lastPedSample = 0;
573 bool isTRUL0IdData = false;
575 if (! in.IsTRUData() ) { // high gain, low gain, LED Mon data - all have the same shaper/sampling
576 AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
577 amp = fitResults.GetAmp();
578 time = fitResults.GetTof();
579 firstPedSample = fFirstPedestalSample;
580 lastPedSample = fLastPedestalSample;
582 else { // TRU data is special, needs its own analyzer
583 AliCaloFitResults fitResults = fRawAnalyzerTRU->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
584 amp = fitResults.GetAmp();
585 time = fitResults.GetTof();
586 firstPedSample = fFirstPedestalSampleTRU;
587 lastPedSample = fLastPedestalSampleTRU;
588 if (in.GetColumn() > n2x2PerTRU) {
589 isTRUL0IdData = true;
595 vector<Int_t> pedSamples;
597 // select earliest bunch
598 unsigned int bunchIndex = 0;
599 unsigned int startBin = bunchlist.at(0).GetStartBin();
600 if (bunchlist.size() > 0) {
601 for(unsigned int ui=1; ui < bunchlist.size(); ui++ ) {
602 if (startBin > bunchlist.at(ui).GetStartBin() ) {
603 startBin = bunchlist.at(ui).GetStartBin();
609 // check bunch for entries in the pedestal sample range
610 Int_t bunchLength = bunchlist.at(bunchIndex).GetLength();
611 const UShort_t *sig = bunchlist.at(bunchIndex).GetData();
614 if (! isTRUL0IdData) { // regular data, can look at pedestals
615 for (Int_t i = 0; i<bunchLength; i++) {
616 timebin = startBin--;
617 if ( firstPedSample<=timebin && timebin<=lastPedSample ) {
618 pedSamples.push_back( sig[i] );
623 else { // TRU L0 Id Data
624 // which TRU the channel belongs to?
625 Int_t iTRUId = in.GetModule()*3 + (iRCU*in.GetBranch() + iRCU);
627 for (Int_t i = 0; i< bunchLength; i++) {
628 for( Int_t j = 0; j < nTRUL0ChannelBits; j++ ){
629 // check if the bit j is 1
630 if( (sig[i] & ( 1 << j )) > 0 ){
631 Int_t iTRUIdInSM = (in.GetColumn() - n2x2PerTRU)*nTRUL0ChannelBits+j;
632 if(iTRUIdInSM < n2x2PerTRU) {
633 Int_t iTRUAbsId = iTRUIdInSM + n2x2PerTRU * iTRUId;
634 // Fill the histograms
635 GetRawsData(kNL0TRU)->Fill(iTRUAbsId);
636 GetRawsData(kTimeL0TRU)->Fill(iTRUAbsId, startBin);
645 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
646 Int_t towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
647 if ( in.IsLowGain() ) {
649 if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) {
650 GetRawsData(kSigLG)->Fill(towerId, amp);
651 GetRawsData(kTimeLG)->Fill(towerId, time);
654 for (Int_t i=0; i<nPed; i++) {
655 GetRawsData(kPedLG)->Fill(towerId, pedSamples[i]);
659 else if ( in.IsHighGain() ) {
661 if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) {
662 GetRawsData(kSigHG)->Fill(towerId, amp);
663 GetRawsData(kTimeHG)->Fill(towerId, time);
666 for (Int_t i=0; i<nPed; i++) {
667 GetRawsData(kPedHG)->Fill(towerId, pedSamples[i]);
671 } // low or high gain
673 else if ( in.IsTRUData() && in.GetColumn()<AliEMCALGeoParams::fgkEMCAL2x2PerTRU) {
674 // for TRU data, the mapping class holds the TRU Int_ternal 2x2 number (0..95) in the Column var..
675 Int_t iTRU = (iRCU*in.GetBranch() + iRCU); //TRU0 is from RCU0, TRU1 from RCU1, TRU2 is from branch B on RCU1
676 Int_t iTRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU
679 if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) {
680 GetRawsData(kSigTRU)->Fill(iTRU2x2Id, amp);
681 GetRawsData(kTimeTRU)->Fill(iTRU2x2Id, time);
684 for (Int_t i=0; i<nPed; i++) {
685 GetRawsData(kPedTRU)->Fill(iTRU2x2Id, pedSamples[i]);
690 else if ( in.IsLEDMonData() ) {
691 // for LED Mon data, the mapping class holds the gain info in the Row variable
692 // and the Strip number in the Column..
693 Int_t gain = in.GetRow();
694 Int_t stripId = iSM*nStripsPerSM + in.GetColumn();
697 nTotalSMLGLEDMon[iSM]++;
698 if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) {
699 GetRawsData(kSigLGLEDMon)->Fill(stripId, amp);
700 GetRawsData(kTimeLGLEDMon)->Fill(stripId, time);
703 for (Int_t i=0; i<nPed; i++) {
704 GetRawsData(kPedLGLEDMon)->Fill(stripId, pedSamples[i]);
708 else if ( gain == 1 ) {
709 nTotalSMHGLEDMon[iSM]++;
710 if ( (amp > fMinSignalHGLEDMon) && (amp < fMaxSignalHGLEDMon) ) {
711 GetRawsData(kSigHGLEDMon)->Fill(stripId, amp);
712 GetRawsData(kTimeHGLEDMon)->Fill(stripId, time);
715 for (Int_t i=0; i<nPed; i++) {
716 GetRawsData(kPedHGLEDMon)->Fill(stripId, pedSamples[i]);
719 } // low or high gain
724 } // nsamples>0 check, some data found for this channel; not only trailer/header
725 }// end while over channel
727 }//end while over DDL's, of input stream
729 // TProfile * p = dynamic_cast<TProfile *>(GetRawsData(kSigHG)) ;
730 ConvertProfile2H(dynamic_cast<TProfile *>(GetRawsData(kSigHG)), fHighEmcHistoH2F) ;
731 Double_t binContent = 0. ;
733 //reset ratio histograms
734 GetRawsData(k2DRatioAmp)->Reset("ICE");
735 GetRawsData(kRatioDist)->Reset("ICE");
736 GetRawsData(kLEDMonRatio)->Reset("ICE");
737 GetRawsData(kLEDMonRatioDist)->Reset("ICE");
738 GetRawsData(k2DRatioAmp)->ResetStats();
739 GetRawsData(kRatioDist)->ResetStats();
740 GetRawsData(kLEDMonRatio)->ResetStats();
741 GetRawsData(kLEDMonRatioDist)->ResetStats();
743 //calculate the ratio of the amplitude and fill the histograms, only if the events type is Calib
744 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
745 for(Int_t ix = 1; ix <= fHighEmcHistoH2F->GetNbinsX(); ix++) {
746 for(Int_t iy = 1; iy <= fHighEmcHistoH2F->GetNbinsY(); iy++) {
747 if(fCalibRefHistoH2F->GetBinContent(ix, iy))binContent = fHighEmcHistoH2F->GetBinContent(ix, iy)/fCalibRefHistoH2F->GetBinContent(ix, iy) ;
748 GetRawsData(k2DRatioAmp)->SetBinContent(ix, iy, binContent);
749 GetRawsData(kRatioDist)->Fill(GetRawsData(k2DRatioAmp)->GetBinContent(ix, iy));
754 //Now for LED monitor system, to calculate the ratio as well
755 Double_t binError = 0. ;
756 // for the binError, we add the relative errors, squared
757 Double_t relativeErrorSqr = 0. ;
759 for(int ib = 1; ib <= fLEDMonRefHistoPro->GetNbinsX(); ib++) {
761 if(fLEDMonRefHistoPro->GetBinContent(ib) != 0) {
762 binContent = GetRawsData(kSigLGLEDMon)->GetBinContent(ib) / fLEDMonRefHistoPro->GetBinContent(ib);
764 relativeErrorSqr = TMath::Power( (fLEDMonRefHistoPro->GetBinError(ib) / fLEDMonRefHistoPro->GetBinContent(ib)), 2);
765 if(GetRawsData(kSigLGLEDMon)->GetBinContent(ib) != 0) {
766 relativeErrorSqr += TMath::Power( (GetRawsData(kSigLGLEDMon)->GetBinError(ib)/GetRawsData(kSigLGLEDMon)->GetBinContent(ib)), 2);
771 relativeErrorSqr = 0;
773 GetRawsData(kLEDMonRatio)->SetBinContent(ib, binContent);
775 binError = sqrt(relativeErrorSqr) * binContent;
776 GetRawsData(kLEDMonRatio)->SetBinError(ib, binError);
777 GetRawsData(kLEDMonRatioDist)->Fill(GetRawsData(kLEDMonRatio)->GetBinContent(ib));
781 // let's also fill the SM and event counter histograms
785 Int_t nTotalHGLEDMon = 0;
786 Int_t nTotalLGLEDMon = 0;
787 for (iSM=0; iSM<fSuperModules; iSM++) {
788 nTotalLG += nTotalSMLG[iSM];
789 nTotalHG += nTotalSMHG[iSM];
790 nTotalTRU += nTotalSMTRU[iSM];
791 nTotalLGLEDMon += nTotalSMLGLEDMon[iSM];
792 nTotalHGLEDMon += nTotalSMHGLEDMon[iSM];
793 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
794 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
795 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
796 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
797 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
800 GetRawsData(kNtotLG)->Fill(nTotalLG);
801 GetRawsData(kNtotHG)->Fill(nTotalHG);
802 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
803 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
804 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
807 SetEventSpecie(saveSpecie) ;
808 // just in case the next rawreader consumer forgets to reset; let's do it here again..
814 //____________________________________________________________________________
815 void AliEMCALQADataMakerRec::MakeDigits()
817 // makes data from Digits
819 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
820 TIter next(fDigitsArray) ;
821 AliEMCALDigit * digit ;
822 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
823 GetDigitsData(0)->Fill( digit->GetAmplitude()) ;
828 //____________________________________________________________________________
829 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
831 // makes data from Digit Tree
833 fDigitsArray->Clear("C") ;
835 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
837 TBranch * branch = digitTree->GetBranch("EMCAL") ;
839 AliWarning("EMCAL branch in Digit Tree not found") ;
841 branch->SetAddress(&fDigitsArray) ;
842 branch->GetEntry(0) ;
848 //____________________________________________________________________________
849 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
851 // makes data from RecPoints
852 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
854 AliError("can't get the branch with the EMCAL clusters !");
858 TObjArray * emcRecPoints = new TObjArray(100) ;
859 emcbranch->SetAddress(&emcRecPoints);
860 emcbranch->GetEntry(0);
862 GetRecPointsData(kRecPM)->Fill(emcRecPoints->GetEntriesFast()) ;
863 TIter next(emcRecPoints) ;
864 AliEMCALRecPoint * rp ;
865 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
866 GetRecPointsData(kRecPE)->Fill(rp->GetEnergy()) ;
867 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
869 emcRecPoints->Delete();
874 //____________________________________________________________________________
875 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
877 //Detector specific actions at start of cycle
881 //____________________________________________________________________________
882 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)
884 //Set fitting algorithm and initialize it if this same algorithm was not set before.
885 //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
887 if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
888 //Do nothing, this same algorithm already set before.
889 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
892 //Initialize the requested algorithm
893 if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
894 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
896 fFittingAlgorithm = fitAlgo;
897 if (fRawAnalyzer) delete fRawAnalyzer; // delete prev. analyzer if existed.
899 if (fitAlgo == kFastFit) {
900 fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
902 else if (fitAlgo == kNeuralNet) {
903 fRawAnalyzer = new AliCaloRawAnalyzerNN();
905 else if (fitAlgo == kLMS) {
906 fRawAnalyzer = new AliCaloRawAnalyzerLMS();
908 else if (fitAlgo == kPeakFinder) {
909 fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
911 else if (fitAlgo == kCrude) {
912 fRawAnalyzer = new AliCaloRawAnalyzerCrude();
915 AliWarning("EMCAL QA invalid fit algorithm choice") ;
922 //_____________________________________________________________________________________
923 void AliEMCALQADataMakerRec::ConvertProfile2H(TProfile * p, TH2 * histo)
926 histo->Reset("ICE") ;
929 Int_t nbinsProf = p->GetNbinsX();
931 // loop through the TProfile p and fill the TH2F histo
934 Double_t binContent = 0;
935 Int_t towerNum = 0; // global tower Id
936 // i = 0; // tower Id within SuperModule
937 Int_t iSM = 0; // SuperModule index
938 Int_t iSMSide = 0; // 0=A, 1=C side
939 Int_t iSMSector = 0; // 2 SM's per sector
941 // indices for 2D plots
945 for (Int_t ibin = 1; ibin <= nbinsProf; ibin++) {
946 towerNum = (Int_t) p->GetBinCenter(ibin);
947 binContent = p->GetBinContent(ibin);
949 // figure out what the tower indices are: col, row within a SuperModule
950 iSM = towerNum/(AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols);
951 col = (towerNum/AliEMCALGeoParams::fgkEMCALRows) % (AliEMCALGeoParams::fgkEMCALCols);
952 row = towerNum % (AliEMCALGeoParams::fgkEMCALRows);
954 //DecodeTowerNum(towerNum, &SM, &col, &row);
955 // then we calculate what the global 2D coord are, based on which SM
960 if (iSMSide == 1) { // C side, shown to the right
961 col2d = col + AliEMCALGeoParams::fgkEMCALCols;
963 else { // A side, shown to the left
967 row2d = row + iSMSector * AliEMCALGeoParams::fgkEMCALRows;
969 histo->SetBinContent(col2d+1, row2d+1, binContent);