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 **************************************************************************/
17 //________________________________________________________________________
19 // A help class for monitoring and calibration tools: MOOD, AMORE etc.,
20 // It can be created and used a la (ctor):
22 //Create the object for making the histograms
23 fPedestals = new AliCaloCalibPedestal( fDetType );
24 // AliCaloCalibPedestal knows how many modules we have for PHOS or EMCAL
25 fNumModules = fPedestals->GetModules();
28 // fPedestals->ProcessEvent(fCaloRawStream);
29 // asked to draw histograms:
30 // fPedestals->GetDeadMap(i)->Draw("col");
32 // fPedestals->GetPeakProfileHighGainRatio((i < fNumModules) ? i : fVisibleModule)->Draw("colz");
34 // The pseudo-code examples above were from the first implementation in MOOD (summer 2007).
35 //________________________________________________________________________
43 #include "AliCaloRawStream.h"
46 #include "AliCaloCalibPedestal.h"
48 ClassImp(AliCaloCalibPedestal)
52 // ctor; initialize everything in order to avoid compiler warnings
53 AliCaloCalibPedestal::AliCaloCalibPedestal(kDetType detectorType) :
57 fPeakMinusPedLowGain(),
58 fPeakMinusPedHighGain(),
59 fPedestalLowGainDiff(),
60 fPedestalHighGainDiff(),
61 fPeakMinusPedLowGainDiff(),
62 fPeakMinusPedHighGainDiff(),
63 fPedestalLowGainRatio(),
64 fPedestalHighGainRatio(),
65 fPeakMinusPedLowGainRatio(),
66 fPeakMinusPedHighGainRatio(),
72 fResurrectedTowers(0),
80 //Default constructor. First we set the detector-type related constants.
81 if (detectorType == kPhos) {
82 fColumns = fgkPhosCols;
84 fModules = fgkPhosModules;
87 //We'll just trust the enum to keep everything in line, so that if detectorType
88 //isn't kPhos then it is kEmCal. Note, however, that this is not necessarily the
89 //case, if someone intentionally gives another number
90 fColumns = fgkEmCalCols;
92 fModules = fgkEmCalModules;
94 fDetType = detectorType;
96 //Then, loop for the requested number of modules
98 for (int i = 0; i < fModules; i++) {
100 name = "hPedlowgain";
102 title = "Pedestals, low gain, module ";
104 fPedestalLowGain.Add(new TProfile2D(name, title,
105 fColumns, 0.0, fColumns,
106 fRows, -fRows, 0.0));
108 //Pedestals, high gain
109 name = "hPedhighgain";
111 title = "Pedestals, high gain, module ";
113 fPedestalHighGain.Add(new TProfile2D(name, title,
114 fColumns, 0.0, fColumns,
115 fRows, -fRows, 0.0));
117 //Peak-Pedestals, low gain
118 name = "hPeakMinusPedlowgain";
120 title = "Peak-Pedestal, low gain, module ";
122 fPeakMinusPedLowGain.Add(new TProfile2D(name, title,
123 fColumns, 0.0, fColumns,
124 fRows, -fRows, 0.0));
126 //Peak-Pedestals, high gain
127 name = "hPeakMinusPedhighgain";
129 title = "Peak-Pedestal, high gain, module ";
131 fPeakMinusPedHighGain.Add(new TProfile2D(name, title,
132 fColumns, 0.0, fColumns,
133 fRows, -fRows, 0.0));
137 title = "Dead map, module ";
139 fDeadMap.Add(new TH2D(name, title, fColumns, 0.0, fColumns,
140 fRows, -fRows, 0.0));
142 }//end for nModules create the histograms
144 //Compress the arrays, in order to remove the empty objects (a 16 slot array is created by default)
145 fPedestalLowGain.Compress();
146 fPedestalHighGain.Compress();
147 fPeakMinusPedLowGain.Compress();
148 fPeakMinusPedHighGain.Compress();
150 //Make them the owners of the profiles, so we don't need to care about deleting them
151 //fPedestalLowGain.SetOwner();
152 //fPedestalHighGain.SetOwner();
153 //fPeakMinusPedLowGain.SetOwner();
154 //fPeakMinusPedHighGain.SetOwner();
159 //_____________________________________________________________________
160 AliCaloCalibPedestal::~AliCaloCalibPedestal()
162 if (fReference) delete fReference;//Delete the reference object, if it has been loaded
163 //TObjArray will delete the histos/profiles when it is deleted.
167 //_____________________________________________________________________
168 AliCaloCalibPedestal::AliCaloCalibPedestal(const AliCaloCalibPedestal &ped) :
172 fPeakMinusPedLowGain(),
173 fPeakMinusPedHighGain(),
174 fPedestalLowGainDiff(),
175 fPedestalHighGainDiff(),
176 fPeakMinusPedLowGainDiff(),
177 fPeakMinusPedHighGainDiff(),
178 fPedestalLowGainRatio(),
179 fPedestalHighGainRatio(),
180 fPeakMinusPedLowGainRatio(),
181 fPeakMinusPedHighGainRatio(),
183 fNEvents(ped.GetNEvents()),
184 fNChanFills(ped.GetNChanFills()),
185 fDeadTowers(ped.GetDeadTowerCount()),
186 fNewDeadTowers(ped.GetDeadTowerNew()),
187 fResurrectedTowers(ped.GetDeadTowerResurrected()),
188 fReference( 0 ), //! note that we do not try to copy the reference info here
189 fDetType(ped.GetDetectorType()),
190 fColumns(ped.GetColumns()),
191 fRows(ped.GetRows()),
192 fModules(ped.GetModules()),
193 fRunNumber(ped.GetRunNumber())
195 // Then the ObjArray ones; we add the histograms rather than trying TObjArray = assignment
196 //DS: this has not really been tested yet..
197 for (int i = 0; i < fModules; i++) {
198 fPedestalLowGain.Add( ped.GetPedProfileLowGain(i) );
199 fPedestalHighGain.Add( ped.GetPedProfileHighGain(i) );
200 fPeakMinusPedLowGain.Add( ped.GetPeakProfileLowGain(i) );
201 fPeakMinusPedHighGain.Add( ped.GetPeakProfileHighGain(i) );
203 fDeadMap.Add( ped.GetDeadMap(i) );
206 //Compress the arrays, in order to remove the empty objects (a 16 slot array is created by default)
207 fPedestalLowGain.Compress();
208 fPedestalHighGain.Compress();
209 fPeakMinusPedLowGain.Compress();
210 fPeakMinusPedHighGain.Compress();
214 // assignment operator; use copy ctor to make life easy..
215 //_____________________________________________________________________
216 AliCaloCalibPedestal& AliCaloCalibPedestal::operator = (const AliCaloCalibPedestal &source)
218 // assignment operator; use copy ctor
219 if (&source == this) return *this;
221 new (this) AliCaloCalibPedestal(source);
225 //_____________________________________________________________________
226 void AliCaloCalibPedestal::Reset()
228 // Reset all arrays/histograms
229 for (int i = 0; i < fModules; i++) {
230 GetPedProfileLowGain(i)->Reset();
231 GetPedProfileHighGain(i)->Reset();
232 GetPeakProfileLowGain(i)->Reset();
233 GetPeakProfileHighGain(i)->Reset();
234 GetDeadMap(i)->Reset();
236 if (!fPedestalLowGainDiff.IsEmpty()) {
237 //This means that the comparison profiles have been created.
239 GetPedProfileLowGainDiff(i)->Reset();
240 GetPedProfileHighGainDiff(i)->Reset();
241 GetPeakProfileLowGainDiff(i)->Reset();
242 GetPeakProfileHighGainDiff(i)->Reset();
244 GetPedProfileLowGainRatio(i)->Reset();
245 GetPedProfileHighGainRatio(i)->Reset();
246 GetPeakProfileLowGainRatio(i)->Reset();
247 GetPeakProfileHighGainRatio(i)->Reset();
254 fResurrectedTowers = 0;
256 //To think about: should fReference be deleted too?... let's not do it this time, at least...
259 //_____________________________________________________________________
260 Bool_t AliCaloCalibPedestal::AddInfo(const AliCaloCalibPedestal *ped)
262 // just do this for the basic histograms/profiles that get filled in ProcessEvent
263 // may not have data for all modules, but let's just Add everything..
264 for (int i = 0; i < fModules; i++) {
265 GetPedProfileLowGain(i)->Add( ped->GetPedProfileLowGain(i) );
266 GetPedProfileHighGain(i)->Add( ped->GetPedProfileHighGain(i) );
267 GetPeakProfileLowGain(i)->Add( ped->GetPeakProfileLowGain(i) );
268 GetPeakProfileHighGain(i)->Add( ped->GetPeakProfileHighGain(i) );
271 // DeadMap; Diff profiles etc would need to be redone after this operation
273 return kTRUE;//We succesfully added info from the supplied object
276 //_____________________________________________________________________
277 Bool_t AliCaloCalibPedestal::ProcessEvent(AliCaloRawStream *in)
279 // Method to process=analyze one event in the data stream
280 if (!in) return kFALSE; //Return right away if there's a null pointer
282 fNEvents++; // one more event
283 int sample, i = 0; //The sample temp, and the sample number in current event.
284 int max = fgkSampleMin, min = fgkSampleMax;//Use these for picking the pedestal
288 sample = in->GetSignal(); //Get the adc signal
289 if (sample < min) min = sample;
290 if (sample > max) max = sample;
292 if ( i >= in->GetTimeLength()) {
293 //If we're here then we're done with this tower
294 gain = 1 - in->IsLowGain();
296 int arrayPos = in->GetModule(); //The modules are numbered starting from 0
297 if (arrayPos >= fModules) {
298 //TODO: return an error message, if appopriate (perhaps if debug>0?)
303 if (arrayPos < 0 || arrayPos >= fModules) {
304 printf("Oh no: arrayPos = %i.\n", arrayPos);
307 fNChanFills++; // one more channel found, and profile to be filled
308 //NOTE: coordinates are (column, row) for the profiles
310 //fill the low gain histograms
311 ((TProfile2D*)fPedestalLowGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, min);
312 ((TProfile2D*)fPeakMinusPedLowGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, max - min);
314 else {//fill the high gain ones
315 ((TProfile2D*)fPedestalHighGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, min);
316 ((TProfile2D*)fPeakMinusPedHighGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, max - min);
319 max = fgkSampleMin; min = fgkSampleMax;
322 }//End if end of tower
324 }//end while, of stream
329 //_____________________________________________________________________
330 Bool_t AliCaloCalibPedestal::SaveHistograms(TString fileName, Bool_t saveEmptyHistos)
332 //Saves all the histograms (or profiles, to be accurate) to the designated file
334 TFile destFile(fileName, "recreate");
336 if (destFile.IsZombie()) {
342 for (int i = 0; i < fModules; i++) {
343 if( ((TProfile2D *)fPeakMinusPedLowGain[i])->GetEntries() || saveEmptyHistos) {
344 fPeakMinusPedLowGain[i]->Write();
346 if( ((TProfile2D *)fPeakMinusPedHighGain[i])->GetEntries() || saveEmptyHistos) {
347 fPeakMinusPedHighGain[i]->Write();
349 if( ((TProfile2D *)fPedestalLowGain[i])->GetEntries() || saveEmptyHistos) {
350 fPedestalLowGain[i]->Write();
352 if( ((TProfile2D *)fPedestalHighGain[i])->GetEntries() || saveEmptyHistos) {
353 fPedestalHighGain[i]->Write();
362 //_____________________________________________________________________
363 Bool_t AliCaloCalibPedestal::LoadReferenceCalib(TString fileName, TString objectName)
366 //Make sure that the histograms created when loading the object are not destroyed as the file object is destroyed
367 TH1::AddDirectory(kFALSE);
369 TFile *sourceFile = new TFile(fileName);
370 if (sourceFile->IsZombie()) {
371 return kFALSE;//We couldn't load the reference
374 if (fReference) delete fReference;//Delete the reference object, if it already exists
377 fReference = (AliCaloCalibPedestal*)sourceFile->Get(objectName);
379 if (!fReference || !(fReference->InheritsFrom(AliCaloCalibPedestal::Class())) || (fReference->GetDetectorType() != fDetType)) {
380 if (fReference) delete fReference;//Delete the object, in case we had an object of the wrong type
387 //Reset the histogram ownership behaviour. NOTE: a better workaround would be good, since this may accidentally set AddDirectory to true, even
388 //if we are called by someone who has set it to false...
389 TH1::AddDirectory(kTRUE);
391 return kTRUE;//We succesfully loaded the object
394 //_____________________________________________________________________
395 void AliCaloCalibPedestal::ValidateComparisonProfiles()
397 //Make sure the comparison histos exist
398 if (!fPedestalLowGainDiff.IsEmpty()) return; //The profiles already exist. We just check one, because they're all created at
402 //Then, loop for the requested number of modules
404 for (int i = 0; i < fModules; i++) {
405 //Pedestals, low gain
406 name = "hPedlowgainDiff";
408 title = "Pedestals difference, low gain, module ";
410 fPedestalLowGainDiff.Add(new TProfile2D(name, title,
411 fColumns, 0.0, fColumns,
412 fRows, -fRows, 0.0));
414 //Pedestals, high gain
415 name = "hPedhighgainDiff";
417 title = "Pedestals difference, high gain, module ";
419 fPedestalHighGainDiff.Add(new TProfile2D(name, title,
420 fColumns, 0.0, fColumns,
421 fRows, -fRows, 0.0));
423 //Peak-Pedestals, low gain
424 name = "hPeakMinusPedlowgainDiff";
426 title = "Peak-Pedestal difference, low gain, module ";
428 fPeakMinusPedLowGainDiff.Add(new TProfile2D(name, title,
429 fColumns, 0.0, fColumns,
430 fRows, -fRows, 0.0));
432 //Peak-Pedestals, high gain
433 name = "hPeakMinusPedhighgainDiff";
435 title = "Peak-Pedestal difference, high gain, module ";
437 fPeakMinusPedHighGainDiff.Add(new TProfile2D(name, title,
438 fColumns, 0.0, fColumns,
439 fRows, -fRows, 0.0));
441 //Pedestals, low gain
442 name = "hPedlowgainRatio";
444 title = "Pedestals ratio, low gain, module ";
446 fPedestalLowGainRatio.Add(new TProfile2D(name, title,
447 fColumns, 0.0, fColumns,
448 fRows, -fRows, 0.0));
450 //Pedestals, high gain
451 name = "hPedhighgainRatio";
453 title = "Pedestals ratio, high gain, module ";
455 fPedestalHighGainRatio.Add(new TProfile2D(name, title,
456 fColumns, 0.0, fColumns,
457 fRows, -fRows, 0.0));
459 //Peak-Pedestals, low gain
460 name = "hPeakMinusPedlowgainRatio";
462 title = "Peak-Pedestal ratio, low gain, module ";
464 fPeakMinusPedLowGainRatio.Add(new TProfile2D(name, title,
465 fColumns, 0.0, fColumns,
466 fRows, -fRows, 0.0));
468 //Peak-Pedestals, high gain
469 name = "hPeakMinusPedhighgainRatio";
471 title = "Peak-Pedestal ratio, high gain, module ";
473 fPeakMinusPedHighGainRatio.Add(new TProfile2D(name, title,
474 fColumns, 0.0, fColumns,
475 fRows, -fRows, 0.0));
477 }//end for nModules create the histograms
480 //_____________________________________________________________________
481 void AliCaloCalibPedestal::ComputeDiffAndRatio()
483 // calculate differences and ratios relative to a reference
484 ValidateComparisonProfiles();//Make sure the comparison histos exist
487 return;//Return if the reference object isn't loaded
490 for (int i = 0; i < fModules; i++) {
491 //Compute the ratio of the histograms
493 ((TProfile2D*)fPedestalLowGainRatio[i])->Divide(GetPedProfileLowGain(i), fReference->GetPedProfileLowGain(i), 1.0, 1.0);
494 ((TProfile2D*)fPedestalHighGainRatio[i])->Divide(GetPedProfileHighGain(i), fReference->GetPedProfileHighGain(i), 1.0, 1.0);
495 ((TProfile2D*)fPeakMinusPedLowGainRatio[i])->Divide(GetPeakProfileLowGain(i), fReference->GetPeakProfileLowGain(i), 1.0, 1.0);
496 ((TProfile2D*)fPeakMinusPedHighGainRatio[i])->Divide(GetPeakProfileHighGain(i), fReference->GetPeakProfileHighGain(i), 1.0, 1.0);
498 //For computing the difference, we cannot simply do TProfile2D->Add(), because that subtracts the sum of all entries,
499 //which means that the mean of the new profile will not be the difference of the means. So do it by hand:
500 for (int j = 0; j <= fColumns; j++) {
501 for (int k = 0; k <= fRows; k++) {
502 int bin = ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->GetBin(j+1, k+1);//Note that we assume here that all histos have the same structure...
503 double diff = fReference->GetPeakProfileHighGain(i)->GetBinContent(bin) - GetPeakProfileHighGain(i)->GetBinContent(bin);
504 ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->SetBinContent(j+1, k+1, diff);
505 ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->SetBinEntries(bin, 1);
507 diff = fReference->GetPeakProfileLowGain(i)->GetBinContent(bin) - GetPeakProfileLowGain(i)->GetBinContent(bin);
508 ((TProfile2D*)fPeakMinusPedLowGainDiff[i])->SetBinContent(j+1, k+1, diff);
509 ((TProfile2D*)fPeakMinusPedLowGainDiff[i])->SetBinEntries(bin, 1);
511 diff = fReference->GetPedProfileHighGain(i)->GetBinContent(bin) - GetPedProfileHighGain(i)->GetBinContent(bin);
512 ((TProfile2D*)fPedestalHighGainDiff[i])->SetBinContent(j+1, k+1, diff);
513 ((TProfile2D*)fPedestalHighGainDiff[i])->SetBinEntries(bin, 1);
515 diff = fReference->GetPedProfileLowGain(i)->GetBinContent(bin) - GetPedProfileLowGain(i)->GetBinContent(bin);
516 ((TProfile2D*)fPedestalLowGainDiff[i])->SetBinContent(j+1, k+1, diff);
517 ((TProfile2D*)fPedestalLowGainDiff[i])->SetBinEntries(bin, 1);
526 //_____________________________________________________________________
527 void AliCaloCalibPedestal::ComputeDeadTowers(int threshold, const char * deadMapFile)
529 //Computes the number of dead towers etc etc into memory, after this you can call the GetDead... -functions
535 char name[512];//Quite a long temp buffer, just in case the filename includes a path
538 snprintf(name, 512, "%s.txt", deadMapFile);
539 fout = new ofstream(name);
540 snprintf(name, 512, "%sdiff.txt", deadMapFile);
541 diff = new ofstream(name);
542 if (!fout->is_open()) {
544 fout = 0;//Set the pointer to empty if the file was not opened
546 if (!diff->is_open()) {
548 fout = 0;//Set the pointer to empty if the file was not opened
552 for (int i = 0; i < fModules; i++) {
553 if (GetPeakProfileHighGain(i)->GetEntries() > 0) { //don't care about empty histos
554 for (int j = 1; j <= fColumns; j++) {
555 for (int k = 1; k <= fRows; k++) {
557 if (GetPeakProfileHighGain(i)->GetBinContent(j, k) < threshold) {//It's dead
558 countTot++;//One more dead total
561 << (fRows - k) << " "
564 << "0" << endl;//Write the status to the deadmap file, if the file is open.
567 if (fReference && fReference->GetPeakProfileHighGain(i)->GetBinContent(j, k) >= threshold) {
568 ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kRecentlyDeceased);
569 countNew++;//This tower wasn't dead before!
572 << (fRows - k) << " "
575 << "0" << endl;//Write the status to the deadmap difference file, if the file is open.
579 ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kDead);//This has been dead before. Nothing new
582 else { //It's ALIVE!!
583 //Don't bother with writing the live ones.
585 // (*fout) << i << " "
586 // << (fRows - k) << " "
589 // << "1" << endl;//Write the status to the deadmap file, if the file is open.
590 if (fReference && fReference->GetPeakProfileHighGain(i)->GetBinContent(j, k) < threshold) {
591 ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kResurrected);
592 countRes++; //This tower was dead before => it's a miracle! :P
595 << (fRows - k) << " "
598 << "1" << endl;//Write the status to the deadmap difference file, if the file is open.
602 ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kAlive);
608 }//end if GetEntries >= 0
617 fDeadTowers = countTot;
618 fNewDeadTowers = countNew;
619 fResurrectedTowers = countRes;