ROOT logo
#ifndef ALICALOCALIBPEDESTAL_H
#define ALICALOCALIBPEDESTAL_H

/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

/* $Id$ */


// \file AliCaloCalibPedestal.h
//   \brief Description:
//   A help class for monitoring and calibration tools: MOOD, AMORE etc.,
//   that can process events from a standard AliCaloRawStreamV3,
//   most usually from LED/pulser runs. It stores signal info as
//   typical (highest) amplitude, and pedestal info in geometrically-binned
//   2D profiles of the detectors (EMCAL and PHOS).
//   Comparisons (ratios and differences) can be done with references.

//   \author: Timo Alho (Jyvaskyla), original version. 
//   [Consultant: D. Silvermyr (ORNL)]
//   Partly based on AliTPCCalibPedestal.
//   
//   \version $Revision$
//   \date $Date$

#include "TProfile.h"
#include "TProfile2D.h"
#include "TH2.h"
#include "TObjArray.h"
#include "AliEMCALGeoParams.h"
class AliCaloRawStreamV3;
class AliCaloAltroMapping;
class AliRawReader;

class AliCaloCalibPedestal : public TObject {
  
 public:

  enum kDetType {kPhos, kEmCal, kNone};//The detector types
  enum kDeadMapEntry{kAlive = 0, kDead, kHot, kWarning, kResurrected, kRecentlyDeceased, kNumDeadMapStates};//The entries being put to the deadmap
  
  AliCaloCalibPedestal(kDetType detectorType = kPhos);
  virtual ~AliCaloCalibPedestal();

  // copy ctor, and '=' operator, are not fully tested/debugged yet
  // at least for now; the reference info is not copied from one to the other
  AliCaloCalibPedestal(AliCaloCalibPedestal &ped); 
  AliCaloCalibPedestal& operator = (AliCaloCalibPedestal &source);
  
  // Event processing methods:  
  Bool_t ProcessEvent(AliRawReader *rawReader);
  Bool_t ProcessEvent(AliCaloRawStreamV3    *in);
  
  // Mapping handling
  AliCaloAltroMapping **GetAltroMapping() const { return fMapping; };
  void  SetAltroMapping(AliCaloAltroMapping **mapp) { fMapping = mapp; };

  // Parameter/cut handling
  void SetParametersFromFile(const char *parameterFile);
  void WriteParametersToFile(const char *parameterFile);

  ////////////////////////////
  //Simple getters
  // Main profiles:
  TProfile2D * GetPedProfileLowGain(int i) {ValidateProfiles(); return (TProfile2D*)fPedestalLowGain[i];};	// Return a pointer to the low-gain pedestal profile
  TProfile2D * GetPedProfileHighGain(int i) {ValidateProfiles(); return (TProfile2D*)fPedestalHighGain[i];};	// Return a pointer to the high-gain pedestal profile
  TProfile * GetPedLEDRefProfileLowGain(int i) {ValidateProfiles(); return (TProfile*)fPedestalLEDRefLowGain[i];};	// Return a pointer to the low-gain LEDRef profile 
  TProfile * GetPedLEDRefProfileHighGain(int i) {ValidateProfiles(); return (TProfile*)fPedestalLEDRefHighGain[i];};	// Return a pointer to the high-gain LEDRef profile 
  TProfile2D * GetPeakProfileLowGain(int i) {ValidateProfiles(); return (TProfile2D*)fPeakMinusPedLowGain[i];};	// Return a pointer to the low-gain peak-pedestal profile
  TProfile2D * GetPeakProfileHighGain(int i) {ValidateProfiles(); return (TProfile2D*)fPeakMinusPedHighGain[i];};	// Return a pointer to the high-gain peak-pedestal profile
  
  // Differences to references:
  TProfile2D * GetPedProfileLowGainDiff(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPedestalLowGainDiff[i];};	// Return a pointer to the low-gain pedestal profile difference
  TProfile2D * GetPedProfileHighGainDiff(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPedestalHighGainDiff[i];};	// Return a pointer to the high-gain pedestal profile difference
  TProfile * GetPedLEDRefProfileLowGainDiff(int i) {ValidateComparisonProfiles(); return (TProfile*)fPedestalLEDRefLowGainDiff[i];};	// Return a pointer to the low-gain LEDRef profile difference
  TProfile * GetPedLEDRefProfileHighGainDiff(int i) {ValidateComparisonProfiles(); return (TProfile*)fPedestalLEDRefHighGainDiff[i];};	// Return a pointer to the high-gain LEDRef profile difference 
  TProfile2D * GetPeakProfileLowGainDiff(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPeakMinusPedLowGainDiff[i];};	// Return a pointer to the low-gain peak-pedestal profile difference
  TProfile2D * GetPeakProfileHighGainDiff(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPeakMinusPedHighGainDiff[i];};	// Return a pointer to the high-gain peak-pedestal profile difference
  
  // Ratio to references:
  TProfile2D * GetPedProfileLowGainRatio(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPedestalLowGainRatio[i];};	// Return a pointer to the low-gain pedestal profile ratio
  TProfile2D * GetPedProfileHighGainRatio(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPedestalHighGainRatio[i];};	// Return a pointer to the high-gain pedestal profile ratio
  TProfile * GetPedLEDRefProfileLowGainRatio(int i) {ValidateComparisonProfiles(); return (TProfile*)fPedestalLEDRefLowGainRatio[i];};	// Return a pointer to the low-gain LEDRef profile ratio
  TProfile * GetPedLEDRefProfileHighGainRatio(int i) {ValidateComparisonProfiles(); return (TProfile*)fPedestalLEDRefHighGainRatio[i];};	// Return a pointer to the high-gain LEDRef profile ratio 
  TProfile2D * GetPeakProfileLowGainRatio(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPeakMinusPedLowGainRatio[i];};	// Return a pointer to the low-gain peak-pedestal profile ratio
  TProfile2D * GetPeakProfileHighGainRatio(int i){ValidateComparisonProfiles(); return (TProfile2D*)fPeakMinusPedHighGainRatio[i];};	// Return a pointer to the high-gain peak-pedestal profile ratio
  
  TH2F * GetPeakHighGainHisto(int i) {ValidateProfiles(); return (TH2F*)fPeakMinusPedHighGainHisto[i];};	// Return a pointer to the high-gain peak-pedestal histo


  TH2D * GetDeadMap(int i) {ValidateProfiles(); return (TH2D*)fDeadMap[i];}
  //void SetDeadMap(int i, TH2D *h) const {((TH2D*)fDeadMap[i])=h;}
	
  Bool_t IsBadChannel(int imod, int icol, int irow) const; 
  void  SetChannelStatus(int imod, int icol, int irow, int status); 
  Int_t GetChannelStatus(int imod, int icol, int irow) const { return  (Int_t)((TH2D*)fDeadMap[imod])->GetBinContent(icol, irow);	}
	
  TObjArray GetDeadMap()  {ValidateProfiles(); return fDeadMap;}
  void SetDeadMap(TObjArray map) {fDeadMap = map;}
	
  // Basic info: getters  
  kDetType GetDetectorType() const {return fDetType;};//Returns if this is a PHOS or EMCAL object
  TString GetCaloString() const {return fCaloString;}; //Returns if this is a PHOS or EMCAL object
  
  int GetColumns() const {return fColumns;}; //The number of columns per module
  int GetRows() const {return fRows;}; //The number of rows per module
  int GetLEDRefs() const {return fLEDRefs;}; //The number of LED references/monitors per module
  int GetModules() const {return fModules;}; //The number of modules
  int GetRowMin() const {return fRowMin;}; //for histo def.
  int GetRowMax() const {return fRowMax;}; //for histo def.
  int GetRowMultiplier() const {return fRowMultiplier;}; //for histo filling

  // RunNumbers : setters and getters
  void SetRunNumber(int runNo) {fRunNumber = runNo;};
  int GetRunNumber() const {return fRunNumber;};
  int GetRefRunNumber() const {if (fReference) return fReference->GetRunNumber(); else return -1;};

  // Possibility to select only some samples for the pedestal calculation
  void SetSelectPedestalSamples(Bool_t flag = kFALSE) {fSelectPedestalSamples = flag;} // select to to use only some range of samples for pedestal calc.
  Bool_t GetSelectPedestalSamples() const {return fSelectPedestalSamples;} // select to to use only some range of samples for pedestal calc.
  void SetFirstPedestalSample(int i) {fFirstPedestalSample = i;} // first sample to use
  void SetLastPedestalSample(int i) {fLastPedestalSample = i;} // last sample to use
  int GetFirstPedestalSample() const {return fFirstPedestalSample;}; // first sample to use
  int GetLastPedestalSample() const {return fLastPedestalSample;}; // last sample to use

  //Set threshold/event fraction for tower warnings
  void SetDeadThreshold(int i) {fDeadThreshold = i;} // peak - pedestal dead threshold
  void SetWarningThreshold(int i) {fWarningThreshold = i;} // peak - pedestal warning threshold
  void SetWarningFraction(double d) {fWarningFraction = d;} // event fraction for warnings
  int GetDeadThreshold() const {return fDeadThreshold;}; // peak - pedestal dead threshold
  int GetWarningThreshold() const {return fWarningThreshold;}; // peak - pedestal warning threshold
  double GetWarningFraction() const {return fWarningFraction;}; // event fraction for warnings
  // hot towers
  void SetHotSigma(double d) {fHotSigma = d;} // rms away from normal
  double GetHotSigma() const {return fHotSigma;}; // rms away from normal

  // Basic counters
  int GetNEvents() const {return fNEvents;};
  int GetNChanFills() const {return fNChanFills;};
  
  /////////////////////////////
  //Analysis functions
  void SetDeadTowerCount(Int_t dead)  {fDeadTowers = dead;};//Returns the number of dead towers, by counting the bins in peak-pedestal smaller than threshold
  int GetDeadTowerCount() const {return fDeadTowers;};//Returns the number of dead towers, by counting the bins in peak-pedestal smaller than threshold
  double GetDeadTowerRatio() const {return fDeadTowers/(double)(fRows*fColumns);}; //returns the percentage of dead towers, relative to a full module
  int GetDeadTowerNew() const {return fNewDeadTowers;}; //return the new dead towers compared to the reference
  int GetDeadTowerResurrected() const {return fResurrectedTowers;}; //The the towers resurrected since the reference run

  void Reset();//Resets the whole class.
  Bool_t AddInfo(AliCaloCalibPedestal *ped);//picks up new info from supplied argument
  
  //////////////////////////////////////////////////////
  //Functions related to comparing this with another (reference) run.
  Bool_t LoadReferenceCalib(TString fileName, TString objectName); //Loads another AliCaloCalibPedestal by name "objectName" from the file "fileName", for reference
  void ComputeDiffAndRatio();//Actually computes the difference and ratio into the histo's in memory
  AliCaloCalibPedestal * GetReference() const {return fReference;}; //Get the reference object. Needed for debug, will probably be removed later
  Bool_t SetReference(AliCaloCalibPedestal *ref);
  void ComputeDeadTowers(const char * deadMapFile = 0);//Computes the dead tower values
  void ComputeHotAndWarningTowers(const char * hotMapFile = 0);//Computes the hot tower values

  //Saving functions
  Bool_t SaveHistograms(TString fileName, Bool_t saveEmptyHistos = kFALSE); //Saves the histograms to a .root file

  void Init() { ValidateProfiles(); } // do basic setup

 private:
  
  void ValidateProfiles(); //Makes sure that basic histos/profiles exist
  void CompressAndSetOwner(); //Makes sure that basic histos/profiles exist
  void ValidateComparisonProfiles(); //Makes sure that fPe..Diff and fPe..Ratio profiles exist
  
  //The histograms. We use a TObjArray instead of a simple array,because this gives automatic streaming properties for the
  //class. A TClonesArray would be more efficient, but it's a bit more difficult to use and it doesn't matter too much
  //since we have only one object per module in the array anyway.
  TObjArray fPedestalLowGain; // pedestal info for low gain
  TObjArray fPedestalHighGain; // pedestal info for high gain
  TObjArray fPedestalLEDRefLowGain; // pedestal LEDRef info for low gain
  TObjArray fPedestalLEDRefHighGain; // pedestal LEDRef info for high gain
  TObjArray fPeakMinusPedLowGain; // (peak-pedestal) info for low gain
  TObjArray fPeakMinusPedHighGain; // (peak-pedestal) info for high gain

  TObjArray fPeakMinusPedHighGainHisto; // (peak-pedestal TH2F) info for high gain, used for hot towers eveluation
  
  //The difference of profiles between this and the reference object
  TObjArray fPedestalLowGainDiff; //!
  TObjArray fPedestalHighGainDiff; //!
  TObjArray fPedestalLEDRefLowGainDiff; //!
  TObjArray fPedestalLEDRefHighGainDiff; //! 
  TObjArray fPeakMinusPedLowGainDiff; //!
  TObjArray fPeakMinusPedHighGainDiff; //!
  
  //The ratio of profiles between this and the reference object
  TObjArray fPedestalLowGainRatio; //!
  TObjArray fPedestalHighGainRatio; //!
  TObjArray fPedestalLEDRefLowGainRatio; //!
  TObjArray fPedestalLEDRefHighGainRatio; //! 
  TObjArray fPeakMinusPedLowGainRatio; //!
  TObjArray fPeakMinusPedHighGainRatio; //!
  
  TObjArray fDeadMap;//The deadmap

  // status counters
  int fNEvents; //# total events processed, 
  int fNChanFills; //# total channel fills (NChan * NEvents if not zero-suppressed)

  //The dead tower counts
  int fDeadTowers; // Number of towers found dead.
  int fNewDeadTowers; //! Towers that have died since the reference run
  int fResurrectedTowers; //! Towers that have been resurrected from the dead, compared to the reference
  
  AliCaloCalibPedestal * fReference; //! A reference object, for comparing the accumulated results to a previous run
  
  kDetType fDetType; //The detector type for this object
  int fColumns;	//The number of columns per module
  int fRows;	//The number of rows per module
  int fLEDRefs;	//The number of LED references/monitors per module
  int fModules;	//The number of modules
  int fRowMin; // Minimum Row number
  int fRowMax; // Maximum now number
  int fRowMultiplier; // Multiplication factor to get proper row range between PHOS and EMCAL
  TString fCaloString; // id for which detector type we have 
  AliCaloAltroMapping **fMapping;    //! Altro Mapping object
  int fRunNumber; //The run number. Needs to be set by the user.
  Bool_t fSelectPedestalSamples; // select to to use only some range of samples for pedestal calc.
  int fFirstPedestalSample; // first sample to use
  int fLastPedestalSample; // last sample to use

  int fDeadThreshold; // Peak - ped threshold used for dead towers evaluation
  int fWarningThreshold; // Peak - ped threshold used for warm/warning towers evaluation
  double fWarningFraction; //if(Peak - ped) > threshold in more than this fraction of event -> tower is assigned kWarning
  double fHotSigma; // if pedestal rms more than fHotSigma away from normal -> tower is assigned kHot

  //Constants needed by the class: EMCAL ones are kept in AliEMCALGeoParams.h
  static const int fgkPhosRows = 64; // number of rows per module for PHOS
  static const int fgkPhosCols = 56; // number of columns per module for PHOS
  static const int fgkPhosLEDRefs = 1; // no LED monitor channels for PHOS, set to 1 just to keep code simpler (also create LEDRef histos for PHOS)
  static const int fgkPhosModules = 5; // number of modules for PHOS

  ClassDef(AliCaloCalibPedestal, 8)

};
    
#endif
 AliCaloCalibPedestal.h:1
 AliCaloCalibPedestal.h:2
 AliCaloCalibPedestal.h:3
 AliCaloCalibPedestal.h:4
 AliCaloCalibPedestal.h:5
 AliCaloCalibPedestal.h:6
 AliCaloCalibPedestal.h:7
 AliCaloCalibPedestal.h:8
 AliCaloCalibPedestal.h:9
 AliCaloCalibPedestal.h:10
 AliCaloCalibPedestal.h:11
 AliCaloCalibPedestal.h:12
 AliCaloCalibPedestal.h:13
 AliCaloCalibPedestal.h:14
 AliCaloCalibPedestal.h:15
 AliCaloCalibPedestal.h:16
 AliCaloCalibPedestal.h:17
 AliCaloCalibPedestal.h:18
 AliCaloCalibPedestal.h:19
 AliCaloCalibPedestal.h:20
 AliCaloCalibPedestal.h:21
 AliCaloCalibPedestal.h:22
 AliCaloCalibPedestal.h:23
 AliCaloCalibPedestal.h:24
 AliCaloCalibPedestal.h:25
 AliCaloCalibPedestal.h:26
 AliCaloCalibPedestal.h:27
 AliCaloCalibPedestal.h:28
 AliCaloCalibPedestal.h:29
 AliCaloCalibPedestal.h:30
 AliCaloCalibPedestal.h:31
 AliCaloCalibPedestal.h:32
 AliCaloCalibPedestal.h:33
 AliCaloCalibPedestal.h:34
 AliCaloCalibPedestal.h:35
 AliCaloCalibPedestal.h:36
 AliCaloCalibPedestal.h:37
 AliCaloCalibPedestal.h:38
 AliCaloCalibPedestal.h:39
 AliCaloCalibPedestal.h:40
 AliCaloCalibPedestal.h:41
 AliCaloCalibPedestal.h:42
 AliCaloCalibPedestal.h:43
 AliCaloCalibPedestal.h:44
 AliCaloCalibPedestal.h:45
 AliCaloCalibPedestal.h:46
 AliCaloCalibPedestal.h:47
 AliCaloCalibPedestal.h:48
 AliCaloCalibPedestal.h:49
 AliCaloCalibPedestal.h:50
 AliCaloCalibPedestal.h:51
 AliCaloCalibPedestal.h:52
 AliCaloCalibPedestal.h:53
 AliCaloCalibPedestal.h:54
 AliCaloCalibPedestal.h:55
 AliCaloCalibPedestal.h:56
 AliCaloCalibPedestal.h:57
 AliCaloCalibPedestal.h:58
 AliCaloCalibPedestal.h:59
 AliCaloCalibPedestal.h:60
 AliCaloCalibPedestal.h:61
 AliCaloCalibPedestal.h:62
 AliCaloCalibPedestal.h:63
 AliCaloCalibPedestal.h:64
 AliCaloCalibPedestal.h:65
 AliCaloCalibPedestal.h:66
 AliCaloCalibPedestal.h:67
 AliCaloCalibPedestal.h:68
 AliCaloCalibPedestal.h:69
 AliCaloCalibPedestal.h:70
 AliCaloCalibPedestal.h:71
 AliCaloCalibPedestal.h:72
 AliCaloCalibPedestal.h:73
 AliCaloCalibPedestal.h:74
 AliCaloCalibPedestal.h:75
 AliCaloCalibPedestal.h:76
 AliCaloCalibPedestal.h:77
 AliCaloCalibPedestal.h:78
 AliCaloCalibPedestal.h:79
 AliCaloCalibPedestal.h:80
 AliCaloCalibPedestal.h:81
 AliCaloCalibPedestal.h:82
 AliCaloCalibPedestal.h:83
 AliCaloCalibPedestal.h:84
 AliCaloCalibPedestal.h:85
 AliCaloCalibPedestal.h:86
 AliCaloCalibPedestal.h:87
 AliCaloCalibPedestal.h:88
 AliCaloCalibPedestal.h:89
 AliCaloCalibPedestal.h:90
 AliCaloCalibPedestal.h:91
 AliCaloCalibPedestal.h:92
 AliCaloCalibPedestal.h:93
 AliCaloCalibPedestal.h:94
 AliCaloCalibPedestal.h:95
 AliCaloCalibPedestal.h:96
 AliCaloCalibPedestal.h:97
 AliCaloCalibPedestal.h:98
 AliCaloCalibPedestal.h:99
 AliCaloCalibPedestal.h:100
 AliCaloCalibPedestal.h:101
 AliCaloCalibPedestal.h:102
 AliCaloCalibPedestal.h:103
 AliCaloCalibPedestal.h:104
 AliCaloCalibPedestal.h:105
 AliCaloCalibPedestal.h:106
 AliCaloCalibPedestal.h:107
 AliCaloCalibPedestal.h:108
 AliCaloCalibPedestal.h:109
 AliCaloCalibPedestal.h:110
 AliCaloCalibPedestal.h:111
 AliCaloCalibPedestal.h:112
 AliCaloCalibPedestal.h:113
 AliCaloCalibPedestal.h:114
 AliCaloCalibPedestal.h:115
 AliCaloCalibPedestal.h:116
 AliCaloCalibPedestal.h:117
 AliCaloCalibPedestal.h:118
 AliCaloCalibPedestal.h:119
 AliCaloCalibPedestal.h:120
 AliCaloCalibPedestal.h:121
 AliCaloCalibPedestal.h:122
 AliCaloCalibPedestal.h:123
 AliCaloCalibPedestal.h:124
 AliCaloCalibPedestal.h:125
 AliCaloCalibPedestal.h:126
 AliCaloCalibPedestal.h:127
 AliCaloCalibPedestal.h:128
 AliCaloCalibPedestal.h:129
 AliCaloCalibPedestal.h:130
 AliCaloCalibPedestal.h:131
 AliCaloCalibPedestal.h:132
 AliCaloCalibPedestal.h:133
 AliCaloCalibPedestal.h:134
 AliCaloCalibPedestal.h:135
 AliCaloCalibPedestal.h:136
 AliCaloCalibPedestal.h:137
 AliCaloCalibPedestal.h:138
 AliCaloCalibPedestal.h:139
 AliCaloCalibPedestal.h:140
 AliCaloCalibPedestal.h:141
 AliCaloCalibPedestal.h:142
 AliCaloCalibPedestal.h:143
 AliCaloCalibPedestal.h:144
 AliCaloCalibPedestal.h:145
 AliCaloCalibPedestal.h:146
 AliCaloCalibPedestal.h:147
 AliCaloCalibPedestal.h:148
 AliCaloCalibPedestal.h:149
 AliCaloCalibPedestal.h:150
 AliCaloCalibPedestal.h:151
 AliCaloCalibPedestal.h:152
 AliCaloCalibPedestal.h:153
 AliCaloCalibPedestal.h:154
 AliCaloCalibPedestal.h:155
 AliCaloCalibPedestal.h:156
 AliCaloCalibPedestal.h:157
 AliCaloCalibPedestal.h:158
 AliCaloCalibPedestal.h:159
 AliCaloCalibPedestal.h:160
 AliCaloCalibPedestal.h:161
 AliCaloCalibPedestal.h:162
 AliCaloCalibPedestal.h:163
 AliCaloCalibPedestal.h:164
 AliCaloCalibPedestal.h:165
 AliCaloCalibPedestal.h:166
 AliCaloCalibPedestal.h:167
 AliCaloCalibPedestal.h:168
 AliCaloCalibPedestal.h:169
 AliCaloCalibPedestal.h:170
 AliCaloCalibPedestal.h:171
 AliCaloCalibPedestal.h:172
 AliCaloCalibPedestal.h:173
 AliCaloCalibPedestal.h:174
 AliCaloCalibPedestal.h:175
 AliCaloCalibPedestal.h:176
 AliCaloCalibPedestal.h:177
 AliCaloCalibPedestal.h:178
 AliCaloCalibPedestal.h:179
 AliCaloCalibPedestal.h:180
 AliCaloCalibPedestal.h:181
 AliCaloCalibPedestal.h:182
 AliCaloCalibPedestal.h:183
 AliCaloCalibPedestal.h:184
 AliCaloCalibPedestal.h:185
 AliCaloCalibPedestal.h:186
 AliCaloCalibPedestal.h:187
 AliCaloCalibPedestal.h:188
 AliCaloCalibPedestal.h:189
 AliCaloCalibPedestal.h:190
 AliCaloCalibPedestal.h:191
 AliCaloCalibPedestal.h:192
 AliCaloCalibPedestal.h:193
 AliCaloCalibPedestal.h:194
 AliCaloCalibPedestal.h:195
 AliCaloCalibPedestal.h:196
 AliCaloCalibPedestal.h:197
 AliCaloCalibPedestal.h:198
 AliCaloCalibPedestal.h:199
 AliCaloCalibPedestal.h:200
 AliCaloCalibPedestal.h:201
 AliCaloCalibPedestal.h:202
 AliCaloCalibPedestal.h:203
 AliCaloCalibPedestal.h:204
 AliCaloCalibPedestal.h:205
 AliCaloCalibPedestal.h:206
 AliCaloCalibPedestal.h:207
 AliCaloCalibPedestal.h:208
 AliCaloCalibPedestal.h:209
 AliCaloCalibPedestal.h:210
 AliCaloCalibPedestal.h:211
 AliCaloCalibPedestal.h:212
 AliCaloCalibPedestal.h:213
 AliCaloCalibPedestal.h:214
 AliCaloCalibPedestal.h:215
 AliCaloCalibPedestal.h:216
 AliCaloCalibPedestal.h:217
 AliCaloCalibPedestal.h:218
 AliCaloCalibPedestal.h:219
 AliCaloCalibPedestal.h:220
 AliCaloCalibPedestal.h:221
 AliCaloCalibPedestal.h:222
 AliCaloCalibPedestal.h:223
 AliCaloCalibPedestal.h:224
 AliCaloCalibPedestal.h:225
 AliCaloCalibPedestal.h:226
 AliCaloCalibPedestal.h:227
 AliCaloCalibPedestal.h:228
 AliCaloCalibPedestal.h:229
 AliCaloCalibPedestal.h:230
 AliCaloCalibPedestal.h:231
 AliCaloCalibPedestal.h:232
 AliCaloCalibPedestal.h:233
 AliCaloCalibPedestal.h:234
 AliCaloCalibPedestal.h:235
 AliCaloCalibPedestal.h:236
 AliCaloCalibPedestal.h:237
 AliCaloCalibPedestal.h:238
 AliCaloCalibPedestal.h:239
 AliCaloCalibPedestal.h:240
 AliCaloCalibPedestal.h:241
 AliCaloCalibPedestal.h:242
 AliCaloCalibPedestal.h:243