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
23 // --- ROOT system ---
24 #include <TClonesArray.h>
31 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
36 #include "AliESDCaloCells.h"
37 #include "AliESDEvent.h"
39 #include "AliEMCALQADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliEMCALDigit.h"
42 #include "AliEMCALRecPoint.h"
43 #include "AliEMCALRawUtils.h"
44 #include "AliEMCALReconstructor.h"
45 #include "AliEMCALRecParam.h"
46 #include "AliRawReader.h"
47 #include "AliCaloRawStreamV3.h"
49 ClassImp(AliEMCALQADataMakerRec)
51 //____________________________________________________________________________
52 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec() :
53 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
54 fSuperModules(4) // FIXME!!! number of SuperModules; 4 for 2009; update default to 12 for later runs..
59 //____________________________________________________________________________
60 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
61 AliQADataMakerRec(), fSuperModules()
64 SetName((const char*)qadm.GetName()) ;
65 SetTitle((const char*)qadm.GetTitle());
66 fSuperModules = qadm.GetSuperModules();
69 //__________________________________________________________________
70 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
73 this->~AliEMCALQADataMakerRec();
74 new(this) AliEMCALQADataMakerRec(qadm);
78 //____________________________________________________________________________
79 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
81 //Detector specific actions at end of cycle
83 AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;
86 //____________________________________________________________________________
87 void AliEMCALQADataMakerRec::InitESDs()
89 //Create histograms to controll ESD
90 const Bool_t expert = kTRUE ;
91 const Bool_t image = kTRUE ;
93 TH1F * h1 = new TH1F("hESDCaloClusterE", "ESDs CaloCluster energy in EMCAL;Energy [MeV];Counts", 200, 0., 20.) ;
95 Add2ESDsList(h1, kESDCaloClusE, !expert, image) ;
97 TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0, 100) ;
99 Add2ESDsList(h2, kESDCaloClusM, !expert, image) ;
101 TH1F * h3 = new TH1F("hESDCaloCellA", "ESDs CaloCell amplitude in EMCAL;Energy [MeV];Counts", 500, 0., 250.) ;
103 Add2ESDsList(h3, kESDCaloCellA, !expert, image) ;
105 TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0, 1000) ;
107 Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
111 //____________________________________________________________________________
112 void AliEMCALQADataMakerRec::InitDigits()
114 // create Digits histograms in Digits subdir
115 const Bool_t expert = kTRUE ;
116 const Bool_t image = kTRUE ;
118 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ;
120 Add2DigitsList(h0, 0, !expert, image) ;
121 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ;
123 Add2DigitsList(h1, 1, !expert, image) ;
126 //____________________________________________________________________________
127 void AliEMCALQADataMakerRec::InitRecPoints()
129 // create Reconstructed Points histograms in RecPoints subdir
130 const Bool_t expert = kTRUE ;
131 const Bool_t image = kTRUE ;
133 TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [MeV];Counts",200, 0.,20.); //GeV
135 Add2RecPointsList(h0,kRecPE, !expert, image);
137 TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
139 Add2RecPointsList(h1,kRecPM, !expert, image);
141 TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
143 Add2RecPointsList(h2,kRecPDigM, !expert, image);
147 //____________________________________________________________________________
148 void AliEMCALQADataMakerRec::InitRaws()
150 // create Raws histograms in Raws subdir
151 const Bool_t expert = kTRUE ;
152 const Bool_t saveCorr = kTRUE ;
153 const Bool_t image = kTRUE ;
155 int nTowersPerSM = 1152; // number of towers in a SuperModule; 24x48
156 int nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
158 // counter info: number of channels per event (bins are SM index)
159 TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
160 fSuperModules, -0.5, fSuperModules-0.5) ;
161 Add2RawsList(h0, kNsmodLG, !expert, image, !saveCorr) ;
162 TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
163 fSuperModules, -0.5, fSuperModules-0.5) ;
164 Add2RawsList(h1, kNsmodHG, !expert, image, !saveCorr) ;
166 // where did max sample occur? (bins are towers)
167 TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ns]",
168 nTot, -0.5, nTot-0.5) ;
169 Add2RawsList(h2, kTimeLG, !expert, image, !saveCorr) ;
170 TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ns]",
171 nTot, -0.5, nTot-0.5) ;
172 Add2RawsList(h3, kTimeHG, !expert, image, !saveCorr) ;
174 // how much above pedestal was the max sample? (bins are towers)
175 TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;??",
176 nTot, -0.5, nTot-0.5) ;
177 Add2RawsList(h4, kSigLG, !expert, image, !saveCorr) ;
178 TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;??",
179 nTot, -0.5, nTot-0.5) ;
180 Add2RawsList(h5, kSigHG, !expert, image, !saveCorr) ;
182 // total counter: channels per event
183 TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
185 Add2RawsList(h6, kNtotLG, !expert, image, !saveCorr) ;
186 TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
188 Add2RawsList(h7, kNtotHG, !expert, image, !saveCorr) ;
190 // pedestal (bins are towers)
191 TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
192 nTot, -0.5, nTot-0.5) ;
193 Add2RawsList(h8, kPedLG, !expert, image, !saveCorr) ;
194 TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
195 nTot, -0.5, nTot-0.5) ;
196 Add2RawsList(h9, kPedHG, !expert, image, !saveCorr) ;
198 // pedestal rms (standard dev = sqrt of variance estimator for pedestal) (bins are towers)
199 TProfile * h10 = new TProfile("hLowEmcalRawPedRMS", "Low Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
200 nTot, -0.5, nTot-0.5) ;
201 Add2RawsList(h10, kPedRMSLG, !expert, image, !saveCorr) ;
202 TProfile * h11 = new TProfile("hHighEmcalRawPedRMS", "High Gain EMC: Pedestal RMS vs towerId;Tower Id;Width [ADC counts]",
203 nTot, -0.5, nTot-0.5) ;
204 Add2RawsList(h11, kPedRMSHG, !expert, image, !saveCorr) ;
208 //____________________________________________________________________________
209 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
211 // make QA data from ESDs
214 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
215 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
216 if( clu->IsEMCAL() ) {
217 GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
221 GetESDsData(kESDCaloClusM)->Fill(nTot) ;
224 AliESDCaloCells* cells = esd->GetEMCALCells();
225 GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
227 for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
228 GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
233 //____________________________________________________________________________
234 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
236 //Fill prepared histograms with Raw digit properties
238 //Raw histogram filling not yet implemented
240 //Need to figure out how to get the info we want without having to
241 //actually run Raw2Digits twice.
242 //I suspect what we actually want is a raw digits method, not a true
243 //emcal raw data method, but this doesn't seem to be allowed in
244 //AliQADataMakerRec.h
246 // For now, to avoid redoing the expensive signal fits we just
247 // look at max vs min of the signal spextra, a la online usage in
248 // AliCaloCalibPedestal
251 AliCaloRawStreamV3 in(rawReader,"EMCAL");
254 int nTowersPerSM = 1152; // number of towers in a SuperModule; 24x48
255 int nRows = 24; // number of rows per SuperModule
257 int sampleMax = 0x3ff; // 1023 = 10-bit range
259 // for the pedestal calculation
260 Bool_t selectPedestalSamples = kTRUE;
261 int firstPedestalSample = 0;
262 int lastPedestalSample = 15;
264 // SM counters; decl. should be safe, assuming we don't get more than 12 SuperModules..
265 int nTotalSMLG[12] = {0};
266 int nTotalSMHG[12] = {0};
268 // indices for the reading
273 // counters, on sample level
274 int i = 0; // the sample number in current event.
279 double meanPed = 0, squaredMean = 0, rmsPed = 0;
281 // start loop over input stream
282 while (in.NextDDL()) {
283 while (in.NextChannel()) {
286 int max = sampleMin, min = sampleMax; // min and max sample values
288 // for the pedestal calculation
289 int sampleSum = 0; // sum of samples
290 int squaredSampleSum = 0; // sum of samples squared
291 int nSum = 0; // number of samples in sum
293 while (in.NextBunch()) {
294 const UShort_t *sig = in.GetSignals();
295 startBin = in.GetStartTimeBin();
296 for (i = 0; i < in.GetBunchLength(); i++) {
300 // check if it's a min or max value
301 if (sample < min) min = sample;
307 // should we add it for the pedestal calculation?
308 if ( (firstPedestalSample<=time && time<=lastPedestalSample) || // sample time in range
309 !selectPedestalSamples ) { // or we don't restrict the sample range.. - then we'll take all
311 squaredSampleSum += sample*sample;
315 } // loop over samples in bunch
316 } // loop over bunches
318 // calculate pedesstal estimate: mean of possibly selected samples
320 meanPed = sampleSum / (1.0 * nSum);
321 squaredMean = squaredSampleSum / (1.0 * nSum);
322 // The variance (rms squared) is equal to the mean of the squares minus the square of the mean..
323 rmsPed = sqrt(squaredMean - meanPed*meanPed);
331 //If we're here then we're done with this tower
332 gain = -1; // init to not valid value
333 if ( in.IsLowGain() ) {
336 else if ( in.IsHighGain() ) {
340 // it should be enough to check the SuperModule info for each DDL really, but let's keep it here for now
341 iSM = in.GetModule(); //The modules are numbered starting from 0
343 if (iSM>=0 && iSM<fSuperModules) { // valid module reading, can go on with filling
345 int towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
348 //fill the low gain histograms, and counters
349 nTotalSMLG[iSM]++; // one more channel found
350 GetRawsData(kSigLG)->Fill(towerId, max - min);
351 GetRawsData(kTimeLG)->Fill(towerId, maxTime);
352 if (nSum>0) { // only fill pedestal info in case it could be calculated
353 GetRawsData(kPedLG)->Fill(towerId, meanPed);
354 GetRawsData(kPedRMSLG)->Fill(towerId, rmsPed);
357 else if (gain == 1) {
358 //fill the high gain ones
359 nTotalSMHG[iSM]++; // one more channel found
360 GetRawsData(kSigHG)->Fill(towerId, max - min);
361 GetRawsData(kTimeHG)->Fill(towerId, maxTime);
362 if (nSum>0) { // only fill pedestal info in case it could be calculated
363 GetRawsData(kPedHG)->Fill(towerId, meanPed);
364 GetRawsData(kPedRMSHG)->Fill(towerId, rmsPed);
369 }// end while over channel
371 }//end while over DDL's, of input stream
373 // let's also fill the SM and event counter histograms
376 for (iSM=0; iSM<fSuperModules; iSM++) {
377 nTotalLG += nTotalSMLG[iSM];
378 nTotalHG += nTotalSMHG[iSM];
379 GetRawsData(kNsmodLG)->Fill(iSM, nTotalSMLG[iSM]);
380 GetRawsData(kNsmodHG)->Fill(iSM, nTotalSMHG[iSM]);
382 GetRawsData(kNtotLG)->Fill(nTotalLG);
383 GetRawsData(kNtotHG)->Fill(nTotalHG);
385 // just in case the next rawreader consumer forgets to reset; let's do it here again..
391 //____________________________________________________________________________
392 void AliEMCALQADataMakerRec::MakeDigits()
394 // makes data from Digits
396 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
397 TIter next(fDigitsArray) ;
398 AliEMCALDigit * digit ;
399 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
400 GetDigitsData(0)->Fill( digit->GetAmp()) ;
405 //____________________________________________________________________________
406 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
408 // makes data from Digit Tree
410 fDigitsArray->Clear() ;
412 fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ;
414 TBranch * branch = digitTree->GetBranch("EMCAL") ;
416 AliWarning("EMCAL branch in Digit Tree not found") ;
418 branch->SetAddress(&fDigitsArray) ;
419 branch->GetEntry(0) ;
425 //____________________________________________________________________________
426 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
428 // makes data from RecPoints
429 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
431 AliError("can't get the branch with the EMCAL clusters !");
435 TObjArray * emcrecpoints = new TObjArray(100) ;
436 emcbranch->SetAddress(&emcrecpoints);
437 emcbranch->GetEntry(0);
439 GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ;
440 TIter next(emcrecpoints) ;
441 AliEMCALRecPoint * rp ;
442 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
443 GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
444 GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
446 emcrecpoints->Delete();
451 //____________________________________________________________________________
452 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
454 //Detector specific actions at start of cycle