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(4), // FIXME!!! number of SuperModules; 4 for 2009; update default to 12 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
310 //Defining histograms binning, each 2D histogram covers all SMs
311 Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
312 Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
313 Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
315 // counter info: number of channels per event (bins are SM index)
316 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
317 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
318 Add2RawsList(h0, kNsmodLG, expert, !image, !saveCorr) ;
319 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
320 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
321 Add2RawsList(h1, kNsmodHG, expert, !image, !saveCorr) ;
323 // where did max sample occur? (bins are towers)
324 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
325 nTot, -0.5, nTot-0.5, profileOption) ;
326 Add2RawsList(h2, kTimeLG, expert, !image, !saveCorr) ;
327 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]",
328 nTot, -0.5, nTot-0.5, profileOption) ;
329 Add2RawsList(h3, kTimeHG, expert, !image, !saveCorr) ;
331 // how much above pedestal was the max sample? (bins are towers)
332 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
333 nTot, -0.5, nTot-0.5, profileOption) ;
334 Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
335 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
336 nTot, -0.5, nTot-0.5, profileOption) ;
337 Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
339 // total counter: channels per event
340 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
342 Add2RawsList(h6, kNtotLG, expert, !image, !saveCorr) ;
343 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
345 Add2RawsList(h7, kNtotHG, expert, !image, !saveCorr) ;
347 // pedestal (bins are towers)
348 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
349 nTot, -0.5, nTot-0.5, profileOption) ;
350 Add2RawsList(h8, kPedLG, expert, !image, !saveCorr) ;
351 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
352 nTot, -0.5, nTot-0.5, profileOption) ;
353 Add2RawsList(h9, kPedHG, expert, !image, !saveCorr) ;
355 //temp 2D amplitude histogram for the current run
356 fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
357 fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
358 //add ratio histograms: to comapre the current run with the reference data
359 TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5,
360 nbinsPhi, -0.5, nbinsPhi-0.5);
361 //settings for display in amore
362 h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}");
363 h15->SetMaximum(2.0);
364 h15->SetMinimum(0.1);
365 h15->SetOption("COLZ");
366 gStyle->SetOptStat(0);
367 Int_t color[] = {4,3,2} ;
368 gStyle->SetPalette(3,color);
369 h15->GetZaxis()->SetNdivisions(3);
370 h15->UseCurrentStyle();
371 h15->SetDirectory(0);
372 Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
374 TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
375 h16->SetMinimum(0.1);
376 h16->SetMaximum(100.);
377 gStyle->SetOptStat(0);
378 h16->UseCurrentStyle();
379 h16->SetDirectory(0);
380 Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
382 // now repeat the same for TRU and LEDMon data
383 Int_t nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
385 // counter info: number of channels per event (bins are SM index)
386 TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
387 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
388 Add2RawsList(hT0, kNsmodTRU, expert, !image, !saveCorr) ;
390 // where did max sample occur? (bins are TRU channels)
391 TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]",
392 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
393 Add2RawsList(hT1, kTimeTRU, expert, !image, !saveCorr) ;
395 // how much above pedestal was the max sample? (bins are TRU channels)
396 TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]",
397 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
398 Add2RawsList(hT2, kSigTRU, expert, !image, !saveCorr) ;
400 // total counter: channels per event
401 TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
403 Add2RawsList(hT3, kNtotTRU, expert, !image, !saveCorr) ;
405 // pedestal (bins are TRU channels)
406 TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]",
407 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
408 Add2RawsList(hT4, kPedTRU, expert, !image, !saveCorr) ;
410 // L0 trigger hits: # of hits (bins are TRU channels)
411 TH1I * hT5 = new TH1I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated", nTot2x2, -0.5, nTot2x2);
413 Add2RawsList(hT5, kNL0TRU, expert, !image, !saveCorr);
415 // L0 trigger hits: average time (bins are TRU channels)
416 TProfile * hT6 = new TProfile("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTot2x2, -0.5, nTot2x2, profileOption);
417 Add2RawsList(hT6, kTimeL0TRU, expert, !image, !saveCorr);
419 // and also LED Mon..
420 // LEDMon has both high and low gain channels, just as regular FEE/towers
421 Int_t nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
423 // counter info: number of channels per event (bins are SM index)
424 TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
425 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
426 Add2RawsList(hL0, kNsmodLGLEDMon, expert, !image, !saveCorr) ;
427 TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
428 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
429 Add2RawsList(hL1, kNsmodHGLEDMon, expert, !image, !saveCorr) ;
431 // where did max sample occur? (bins are strips)
432 TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
433 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
434 Add2RawsList(hL2, kTimeLGLEDMon, expert, !image, !saveCorr) ;
435 TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]",
436 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
437 Add2RawsList(hL3, kTimeHGLEDMon, expert, !image, !saveCorr) ;
439 // how much above pedestal was the max sample? (bins are strips)
440 TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
441 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
442 Add2RawsList(hL4, kSigLGLEDMon, expert, !image, !saveCorr) ;
443 TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
444 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
445 Add2RawsList(hL5, kSigHGLEDMon, expert, !image, !saveCorr) ;
447 // total counter: channels per event
448 TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
450 Add2RawsList(hL6, kNtotLGLEDMon, expert, !image, !saveCorr) ;
451 TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
453 Add2RawsList(hL7, kNtotHGLEDMon, expert, !image, !saveCorr) ;
455 // pedestal (bins are strips)
456 TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
457 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
458 Add2RawsList(hL8, kPedLGLEDMon, expert, !image, !saveCorr) ;
459 TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
460 nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
461 Add2RawsList(hL9, kPedHGLEDMon, expert, !image, !saveCorr) ;
463 //add two histograms for shifter from the LED monitor system: comapre LED monitor with the reference run
464 //to be used for decision whether we need to change reference data
465 TH1F * hL10 = new TH1F("hMaxMinusMinLEDMonRatio", "LEDMon amplitude, Ratio to reference run", nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
466 //settings for display in amore
467 hL10->SetTitle("Amplitude_{LEDMon current}/Amplitude_{LEDMon reference}");
468 hL10->SetMaximum(2.0);
469 hL10->SetMinimum(0.1);
470 gStyle->SetOptStat(0);
471 hL10->UseCurrentStyle();
472 hL10->SetDirectory(0);
473 // hL10->SetOption("E");
474 Add2RawsList(hL10, kLEDMonRatio, !expert, image, !saveCorr) ;
476 TH1F * hL11 = new TH1F("hMaxMinusMinLEDMonRatioDist", "LEDMon amplitude, Ratio distribution", nTot, 0, 2);
477 hL11->SetMinimum(0.1) ;
478 hL11->SetMaximum(100.);
479 gStyle->SetOptStat(0);
480 hL11->UseCurrentStyle();
481 hL11->SetDirectory(0);
482 Add2RawsList(hL11, kLEDMonRatioDist, !expert, image, !saveCorr) ;
484 GetCalibRefFromOCDB() ;
487 //____________________________________________________________________________
488 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
490 // make QA data from ESDs
493 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
494 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
495 if( clu->IsEMCAL() ) {
496 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
500 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
503 AliESDCaloCells* cells = esd->GetEMCALCells();
504 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
506 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
507 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
512 //____________________________________________________________________________
513 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
515 // make sure EMCal was readout during the event
516 Int_t emcID = AliDAQ::DetectorID("EMCAL"); // bit 18..
517 const UInt_t *detPattern = rawReader->GetDetectorPattern();
518 UInt_t emcInReadout = ( ((1 << emcID) & detPattern[0]) >> emcID);
519 if (! emcInReadout) return; // no poInt_t in looking at this event, if no EMCal data
523 AliCaloRawStreamV3 in(rawReader,"EMCAL");
524 rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's
526 AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
528 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
529 SetEventSpecie(AliRecoParam::kCalib) ;
532 const Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
533 const Int_t nRows = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
534 const Int_t nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
535 const Int_t n2x2PerSM = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
536 const Int_t n2x2PerTRU = AliEMCALGeoParams::fgkEMCAL2x2PerTRU;
538 // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
539 Int_t nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules] = {0};
540 Int_t nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules] = {0};
541 Int_t nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules] = {0};
542 Int_t nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
543 Int_t nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
545 const Int_t nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
546 Int_t iSM = 0; // SuperModule index
547 // start loop over input stream
548 while (in.NextDDL()) {
549 Int_t iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
550 fRawAnalyzer->SetIsZeroSuppressed( in.GetZeroSupp() );
552 while (in.NextChannel()) {
553 iSM = in.GetModule(); // SuperModule
554 //prInt_tf("iSM %d DDL %d", iSM, in.GetDDLNumber());
555 if (iSM>=0 && iSM<fSuperModules) { // valid module reading
558 vector<AliCaloBunchInfo> bunchlist;
559 while (in.NextBunch()) {
560 nsamples += in.GetBunchLength();
561 bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
564 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
567 // indices for pedestal calc.
568 Int_t firstPedSample = 0;
569 Int_t lastPedSample = 0;
570 bool isTRUL0IdData = false;
572 if (! in.IsTRUData() ) { // high gain, low gain, LED Mon data - all have the same shaper/sampling
573 AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
574 amp = fitResults.GetAmp();
575 time = fitResults.GetTof();
576 firstPedSample = fFirstPedestalSample;
577 lastPedSample = fLastPedestalSample;
579 else { // TRU data is special, needs its own analyzer
580 AliCaloFitResults fitResults = fRawAnalyzerTRU->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
581 amp = fitResults.GetAmp();
582 time = fitResults.GetTof();
583 firstPedSample = fFirstPedestalSampleTRU;
584 lastPedSample = fLastPedestalSampleTRU;
585 if (in.GetColumn() > n2x2PerTRU) {
586 isTRUL0IdData = true;
592 vector<Int_t> pedSamples;
594 // select earliest bunch
595 unsigned int bunchIndex = 0;
596 unsigned int startBin = bunchlist.at(0).GetStartBin();
597 if (bunchlist.size() > 0) {
598 for(unsigned int ui=1; ui < bunchlist.size(); ui++ ) {
599 if (startBin > bunchlist.at(ui).GetStartBin() ) {
600 startBin = bunchlist.at(ui).GetStartBin();
606 // check bunch for entries in the pedestal sample range
607 Int_t bunchLength = bunchlist.at(bunchIndex).GetLength();
608 const UShort_t *sig = bunchlist.at(bunchIndex).GetData();
611 if (! isTRUL0IdData) { // regular data, can look at pedestals
612 for (Int_t i = 0; i<bunchLength; i++) {
613 timebin = startBin--;
614 if ( firstPedSample<=timebin && timebin<=lastPedSample ) {
615 pedSamples.push_back( sig[i] );
620 else { // TRU L0 Id Data
621 // which TRU the channel belongs to?
622 Int_t iTRUId = in.GetModule()*3 + (iRCU*in.GetBranch() + iRCU);
624 for (Int_t i = 0; i< bunchLength; i++) {
625 for( Int_t j = 0; j < nTRUL0ChannelBits; j++ ){
626 // check if the bit j is 1
627 if( (sig[i] & ( 1 << j )) > 0 ){
628 Int_t iTRUIdInSM = (in.GetColumn() - n2x2PerTRU)*nTRUL0ChannelBits+j;
629 if(iTRUIdInSM < n2x2PerTRU) {
630 Int_t iTRUAbsId = iTRUIdInSM + n2x2PerTRU * iTRUId;
631 // Fill the histograms
632 GetRawsData(kNL0TRU)->Fill(iTRUAbsId);
633 GetRawsData(kTimeL0TRU)->Fill(iTRUAbsId, startBin);
642 if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
643 Int_t towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
644 if ( in.IsLowGain() ) {
646 if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) {
647 GetRawsData(kSigLG)->Fill(towerId, amp);
648 GetRawsData(kTimeLG)->Fill(towerId, time);
651 for (Int_t i=0; i<nPed; i++) {
652 GetRawsData(kPedLG)->Fill(towerId, pedSamples[i]);
656 else if ( in.IsHighGain() ) {
658 if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) {
659 GetRawsData(kSigHG)->Fill(towerId, amp);
660 GetRawsData(kTimeHG)->Fill(towerId, time);
663 for (Int_t i=0; i<nPed; i++) {
664 GetRawsData(kPedHG)->Fill(towerId, pedSamples[i]);
668 } // low or high gain
670 else if ( in.IsTRUData() && in.GetColumn()<AliEMCALGeoParams::fgkEMCAL2x2PerTRU) {
671 // for TRU data, the mapping class holds the TRU Int_ternal 2x2 number (0..95) in the Column var..
672 Int_t iTRU = (iRCU*in.GetBranch() + iRCU); //TRU0 is from RCU0, TRU1 from RCU1, TRU2 is from branch B on RCU1
673 Int_t iTRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU
676 if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) {
677 GetRawsData(kSigTRU)->Fill(iTRU2x2Id, amp);
678 GetRawsData(kTimeTRU)->Fill(iTRU2x2Id, time);
681 for (Int_t i=0; i<nPed; i++) {
682 GetRawsData(kPedTRU)->Fill(iTRU2x2Id, pedSamples[i]);
687 else if ( in.IsLEDMonData() ) {
688 // for LED Mon data, the mapping class holds the gain info in the Row variable
689 // and the Strip number in the Column..
690 Int_t gain = in.GetRow();
691 Int_t stripId = iSM*nStripsPerSM + in.GetColumn();
694 nTotalSMLGLEDMon[iSM]++;
695 if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) {
696 GetRawsData(kSigLGLEDMon)->Fill(stripId, amp);
697 GetRawsData(kTimeLGLEDMon)->Fill(stripId, time);
700 for (Int_t i=0; i<nPed; i++) {
701 GetRawsData(kPedLGLEDMon)->Fill(stripId, pedSamples[i]);
705 else if ( gain == 1 ) {
706 nTotalSMHGLEDMon[iSM]++;
707 if ( (amp > fMinSignalHGLEDMon) && (amp < fMaxSignalHGLEDMon) ) {
708 GetRawsData(kSigHGLEDMon)->Fill(stripId, amp);
709 GetRawsData(kTimeHGLEDMon)->Fill(stripId, time);
712 for (Int_t i=0; i<nPed; i++) {
713 GetRawsData(kPedHGLEDMon)->Fill(stripId, pedSamples[i]);
716 } // low or high gain
721 } // nsamples>0 check, some data found for this channel; not only trailer/header
722 }// end while over channel
724 }//end while over DDL's, of input stream
726 // TProfile * p = dynamic_cast<TProfile *>(GetRawsData(kSigHG)) ;
727 ConvertProfile2H(dynamic_cast<TProfile *>(GetRawsData(kSigHG)), fHighEmcHistoH2F) ;
728 Double_t binContent = 0. ;
730 //reset ratio histograms
731 GetRawsData(k2DRatioAmp)->Reset("ICE");
732 GetRawsData(kRatioDist)->Reset("ICE");
733 GetRawsData(kLEDMonRatio)->Reset("ICE");
734 GetRawsData(kLEDMonRatioDist)->Reset("ICE");
735 GetRawsData(k2DRatioAmp)->ResetStats();
736 GetRawsData(kRatioDist)->ResetStats();
737 GetRawsData(kLEDMonRatio)->ResetStats();
738 GetRawsData(kLEDMonRatioDist)->ResetStats();
740 //calculate the ratio of the amplitude and fill the histograms, only if the events type is Calib
741 if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) {
742 for(Int_t ix = 1; ix <= fHighEmcHistoH2F->GetNbinsX(); ix++) {
743 for(Int_t iy = 1; iy <= fHighEmcHistoH2F->GetNbinsY(); iy++) {
744 if(fCalibRefHistoH2F->GetBinContent(ix, iy))binContent = fHighEmcHistoH2F->GetBinContent(ix, iy)/fCalibRefHistoH2F->GetBinContent(ix, iy) ;
745 GetRawsData(k2DRatioAmp)->SetBinContent(ix, iy, binContent);
746 GetRawsData(kRatioDist)->Fill(GetRawsData(k2DRatioAmp)->GetBinContent(ix, iy));
751 //Now for LED monitor system, to calculate the ratio as well
752 Double_t binError = 0. ;
753 // for the binError, we add the relative errors, squared
754 Double_t relativeErrorSqr = 0. ;
756 for(int ib = 1; ib <= fLEDMonRefHistoPro->GetNbinsX(); ib++) {
758 if(fLEDMonRefHistoPro->GetBinContent(ib) != 0) {
759 binContent = GetRawsData(kSigLGLEDMon)->GetBinContent(ib) / fLEDMonRefHistoPro->GetBinContent(ib);
760 relativeErrorSqr = TMath::Power(fLEDMonRefHistoPro->GetBinError(ib) / fLEDMonRefHistoPro->GetBinContent(ib), 2);
764 relativeErrorSqr = 0;
766 GetRawsData(kLEDMonRatio)->SetBinContent(ib, binContent);
768 if(GetRawsData(kLEDMonRatio)->GetBinContent(ib) != 0) {
769 relativeErrorSqr += TMath::Power(GetRawsData(kSigLGLEDMon)->GetBinError(ib)/GetRawsData(kLEDMonRatio)->GetBinContent(ib), 2);
772 binError = sqrt(relativeErrorSqr) * binContent;
773 GetRawsData(kLEDMonRatio)->SetBinError(ib, binError);
774 GetRawsData(kLEDMonRatioDist)->Fill(GetRawsData(kLEDMonRatio)->GetBinContent(ib));
778 // let's also fill the SM and event counter histograms
782 Int_t nTotalHGLEDMon = 0;
783 Int_t nTotalLGLEDMon = 0;
784 for (iSM=0; iSM<fSuperModules; iSM++) {
785 nTotalLG += nTotalSMLG[iSM];
786 nTotalHG += nTotalSMHG[iSM];
787 nTotalTRU += nTotalSMTRU[iSM];
788 nTotalLGLEDMon += nTotalSMLGLEDMon[iSM];
789 nTotalHGLEDMon += nTotalSMHGLEDMon[iSM];
790 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
791 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
792 GetRawsData(kNsmodTRU)->Fill(iSM, nTotalSMTRU[iSM]);
793 GetRawsData(kNsmodLGLEDMon)->Fill(iSM, nTotalSMLGLEDMon[iSM]);
794 GetRawsData(kNsmodHGLEDMon)->Fill(iSM, nTotalSMHGLEDMon[iSM]);
797 GetRawsData(kNtotLG)->Fill(nTotalLG);
798 GetRawsData(kNtotHG)->Fill(nTotalHG);
799 GetRawsData(kNtotTRU)->Fill(nTotalTRU);
800 GetRawsData(kNtotLGLEDMon)->Fill(nTotalLGLEDMon);
801 GetRawsData(kNtotHGLEDMon)->Fill(nTotalHGLEDMon);
804 SetEventSpecie(saveSpecie) ;
805 // just in case the next rawreader consumer forgets to reset; let's do it here again..
811 //____________________________________________________________________________
812 void AliEMCALQADataMakerRec::MakeDigits()
814 // makes data from Digits
816 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
817 TIter next(fDigitsArray) ;
818 AliEMCALDigit * digit ;
819 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
820 GetDigitsData(0)->Fill( digit->GetAmplitude()) ;
825 //____________________________________________________________________________
826 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
828 // makes data from Digit Tree
830 fDigitsArray->Clear() ;
832 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
834 TBranch * branch = digitTree->GetBranch("EMCAL") ;
836 AliWarning("EMCAL branch in Digit Tree not found") ;
838 branch->SetAddress(&fDigitsArray) ;
839 branch->GetEntry(0) ;
845 //____________________________________________________________________________
846 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
848 // makes data from RecPoints
849 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
851 AliError("can't get the branch with the EMCAL clusters !");
855 TObjArray * emcRecPoints = new TObjArray(100) ;
856 emcbranch->SetAddress(&emcRecPoints);
857 emcbranch->GetEntry(0);
859 GetRecPointsData(kRecPM)->Fill(emcRecPoints->GetEntriesFast()) ;
860 TIter next(emcRecPoints) ;
861 AliEMCALRecPoint * rp ;
862 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
863 GetRecPointsData(kRecPE)->Fill(rp->GetEnergy()) ;
864 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
866 emcRecPoints->Delete();
871 //____________________________________________________________________________
872 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
874 //Detector specific actions at start of cycle
878 //____________________________________________________________________________
879 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)
881 //Set fitting algorithm and initialize it if this same algorithm was not set before.
882 //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
884 if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
885 //Do nothing, this same algorithm already set before.
886 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
889 //Initialize the requested algorithm
890 if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
891 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
893 fFittingAlgorithm = fitAlgo;
894 if (fRawAnalyzer) delete fRawAnalyzer; // delete prev. analyzer if existed.
896 if (fitAlgo == kFastFit) {
897 fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
899 else if (fitAlgo == kNeuralNet) {
900 fRawAnalyzer = new AliCaloRawAnalyzerNN();
902 else if (fitAlgo == kLMS) {
903 fRawAnalyzer = new AliCaloRawAnalyzerLMS();
905 else if (fitAlgo == kPeakFinder) {
906 fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
908 else if (fitAlgo == kCrude) {
909 fRawAnalyzer = new AliCaloRawAnalyzerCrude();
912 AliWarning("EMCAL QA invalid fit algorithm choice") ;
919 //_____________________________________________________________________________________
920 void AliEMCALQADataMakerRec::ConvertProfile2H(TProfile * p, TH2 * histo)
922 // set some histogram defaults
924 //histo->SetStats(kFALSE); // no statistics box shown
925 Int_t nbinsProf = p->GetNbinsX();
927 // loop through the TProfile p and fill the TH2F histo
930 Double_t binContent = 0;
931 Int_t towerNum = 0; // global tower Id
932 // i = 0; // tower Id within SuperModule
933 Int_t iSM = 0; // SuperModule index
934 Int_t iSMSide = 0; // 0=A, 1=C side
935 Int_t iSMSector = 0; // 2 SM's per sector
937 // indices for 2D plots
941 for (Int_t ibin = 1; ibin <= nbinsProf; ibin++) {
942 towerNum = (Int_t) p->GetBinCenter(ibin);
943 binContent = p->GetBinContent(ibin);
945 // figure out what the tower indices are: col, row within a SuperModule
946 iSM = towerNum/(AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols);
947 col = (towerNum/AliEMCALGeoParams::fgkEMCALRows) % (AliEMCALGeoParams::fgkEMCALCols);
948 row = towerNum % (AliEMCALGeoParams::fgkEMCALRows);
950 //DecodeTowerNum(towerNum, &SM, &col, &row);
951 // then we calculate what the global 2D coord are, based on which SM
956 if (iSMSide == 1) { // C side, shown to the right
957 col2d = col + AliEMCALGeoParams::fgkEMCALCols;
959 else { // A side, shown to the left
963 row2d = row + iSMSector * AliEMCALGeoParams::fgkEMCALRows;
965 histo->SetBinContent(col2d+1, row2d+1, binContent);