CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalDbService.cc
Go to the documentation of this file.
1 //
2 // F.Ratnikov (UMd), Aug. 9, 2005
3 //
4 
6 
10 
13 
15 
16 #include <cmath>
17 
19  mPedestals (0), mPedestalWidths (0),
20  mGains (0), mGainWidths (0),
21  mQieShapeCache (0), mQIEData(0),
22  mElectronicsMap(0),
23  mRespCorrs(0),
24  mL1TriggerObjects(0),
25  mTimeCorrs(0),
26  mLUTCorrs(0),
27  mPFCorrs(0),
28  mLutMetadata(0),
29  mUpdateCalibrations (true), mUpdateCalibWidths(true)
30  {}
31 
32 
34 {
36  {
38  mUpdateCalibrations = false;
39  }
40  return mCalibSet.getCalibrations(fId);
41 }
42 
44 {
45  if (mUpdateCalibWidths)
46  {
48  mUpdateCalibWidths = false;
49  }
51 }
52 
54  // we use the set of ids for pedestals as the master list
55  if ((!mPedestals) || (!mGains) || (!mQIEData) || (!mRespCorrs) || (!mTimeCorrs) || (!mLUTCorrs) ) return;
56 
57  std::vector<DetId> ids=mPedestals->getAllChannels();
58  bool pedsInADC = mPedestals->isADC();
59  // clear the calibrations set
60  mCalibSet.clear();
61  // loop!
62  HcalCalibrations tool;
63 
64  // std::cout << " length of id-vector: " << ids.size() << std::endl;
65  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
66  // make
67  bool ok=makeHcalCalibration(*id,&tool,pedsInADC);
68  // store
69  if (ok) mCalibSet.setCalibrations(*id,tool);
70  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
71  }
72  mCalibSet.sort();
73 }
74 
76  // we use the set of ids for pedestal widths as the master list
77  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
78 
79  std::vector<DetId> ids=mPedestalWidths->getAllChannels();
80  bool pedsInADC = mPedestalWidths->isADC();
81  // clear the calibrations set
83  // loop!
85 
86  // std::cout << " length of id-vector: " << ids.size() << std::endl;
87  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
88  // make
89  bool ok=makeHcalCalibrationWidth(*id,&tool,pedsInADC);
90  // store
91  if (ok) mCalibWidthSet.setCalibrationWidths(*id,tool);
92  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
93  }
95 }
96 
97 bool HcalDbService::makeHcalCalibration (const HcalGenericDetId& fId, HcalCalibrations* fObject, bool pedestalInADC) const {
98  if (fObject) {
99  const HcalPedestal* pedestal = getPedestal (fId);
100  const HcalGain* gain = getGain (fId);
101  const HcalRespCorr* respcorr = getHcalRespCorr (fId);
102  const HcalTimeCorr* timecorr = getHcalTimeCorr (fId);
103  const HcalLUTCorr* lutcorr = getHcalLUTCorr (fId);
104 
105  if (pedestalInADC) {
106  const HcalQIEShape* shape=getHcalShape();
107  const HcalQIECoder* coder=getHcalCoder(fId);
108  if (pedestal && gain && shape && coder && respcorr && timecorr && lutcorr) {
109  float pedTrue[4];
110  for (int i=0; i<4; i++) {
111  float x=pedestal->getValues()[i];
112  int x1=(int)std::floor(x);
113  int x2=(int)std::floor(x+1);
114  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
115  float y2=coder->charge(*shape,x2,i);
116  float y1=coder->charge(*shape,x1,i);
117  pedTrue[i]=(y2-y1)*(x-x1)+y1;
118  }
119  *fObject = HcalCalibrations (gain->getValues (), pedTrue, respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
120  return true;
121  }
122  } else {
123  if (pedestal && gain && respcorr && timecorr && lutcorr) {
124  *fObject = HcalCalibrations (gain->getValues (), pedestal->getValues (), respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
125  return true;
126  }
127  }
128  }
129  return false;
130 }
131 
133  HcalCalibrationWidths* fObject, bool pedestalInADC) const {
134  if (fObject) {
135  const HcalPedestalWidth* pedestalwidth = getPedestalWidth (fId);
136  const HcalGainWidth* gainwidth = getGainWidth (fId);
137  if (pedestalInADC) {
138  const HcalQIEShape* shape=getHcalShape();
139  const HcalQIECoder* coder=getHcalCoder(fId);
140  if (pedestalwidth && gainwidth && shape && coder) {
141  float pedTrueWidth[4];
142  for (int i=0; i<4; i++) {
143  float x=pedestalwidth->getWidth(i);
144  // assume QIE is linear in low range and use x1=0 and x2=1
145  // y = (y2-y1) * (x) [do not add any constant, only scale!]
146  float y2=coder->charge(*shape,1,i);
147  float y1=coder->charge(*shape,0,i);
148  pedTrueWidth[i]=(y2-y1)*x;
149  }
150  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
151  return true;
152  }
153  } else {
154  if (pedestalwidth && gainwidth) {
155  float pedestalWidth [4];
156  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
157  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedestalWidth);
158  return true;
159  }
160  }
161  }
162  return false;
163 }
164 
165 
167  if (mRespCorrs) {
168  return mRespCorrs->getValues (fId);
169  }
170  return 0;
171 }
172 
174  if (mPedestals) {
175  return mPedestals->getValues (fId);
176  }
177  return 0;
178 }
179 
181  if (mPedestalWidths) {
182  return mPedestalWidths->getValues (fId);
183  }
184  return 0;
185 }
186 
188  if (mGains) {
189  return mGains->getValues(fId);
190  }
191  return 0;
192 }
193 
195  if (mGainWidths) {
196  return mGainWidths->getValues (fId);
197  }
198  return 0;
199 }
200 
202  if (mQIEData) {
203  return mQIEData->getCoder (fId);
204  }
205  return 0;
206 }
207 
209  if (mQIEData) {
210  return &mQIEData->getShape ();
211  }
212  return 0;
213 }
214 
216  return mElectronicsMap;
217 }
218 
220 {
221  return mL1TriggerObjects->getValues (fId);
222 }
223 
225 {
226  return mChannelQuality->getValues (fId);
227 }
228 
230 {
231  return mZSThresholds->getValues (fId);
232 }
233 
235  if (mTimeCorrs) {
236  return mTimeCorrs->getValues (fId);
237  }
238  return 0;
239 }
240 
242  if (mLUTCorrs) {
243  return mLUTCorrs->getValues (fId);
244  }
245  return 0;
246 }
247 
249  if (mPFCorrs) {
250  return mPFCorrs->getValues (fId);
251  }
252  return 0;
253 }
254 
256  return mLutMetadata;
257 }
258 
const HcalLUTCorrs * mLUTCorrs
Definition: HcalDbService.h:84
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
int i
Definition: DBlmapReader.cc:9
bool mUpdateCalibWidths
Definition: HcalDbService.h:90
bool isADC() const
Definition: HcalPedestals.h:23
bool mUpdateCalibrations
Definition: HcalDbService.h:90
const HcalQIECoder * getCoder(DetId fId) const
get QIE parameters
Definition: HcalQIEData.h:38
const HcalChannelStatus * getHcalChannelStatus(const HcalGenericDetId &fId) const
const HcalL1TriggerObjects * mL1TriggerObjects
Definition: HcalDbService.h:82
const HcalPFCorrs * mPFCorrs
Definition: HcalDbService.h:85
const HcalGains * mGains
Definition: HcalDbService.h:74
const HcalQIEShape & getShape() const
get basic shape
Definition: HcalQIEData.h:36
const HcalChannelQuality * mChannelQuality
Definition: HcalDbService.h:78
const HcalLutMetadata * mLutMetadata
Definition: HcalDbService.h:86
const float * getValues() const
get value for all capId = 0..3
Definition: HcalGainWidth.h:17
const HcalPedestalWidth * getPedestalWidth(const HcalGenericDetId &fId) const
const HcalRespCorr * getHcalRespCorr(const HcalGenericDetId &fId) const
const HcalLUTCorr * getHcalLUTCorr(const HcalGenericDetId &fId) const
const HcalTimeCorr * getHcalTimeCorr(const HcalGenericDetId &fId) const
void buildCalibWidths() const
HcalCalibrationsSet mCalibSet
Definition: HcalDbService.h:88
std::vector< DetId > getAllChannels() const
const HcalL1TriggerObject * getHcalL1TriggerObject(const HcalGenericDetId &fId) const
bool makeHcalCalibration(const HcalGenericDetId &fId, HcalCalibrations *fObject, bool pedestalInADC) const
const HcalLutMetadata * getHcalLutMetadata() const
bool makeHcalCalibrationWidth(const HcalGenericDetId &fId, HcalCalibrationWidths *fObject, bool pedestalInADC) const
void setCalibrationWidths(const DetId id, const HcalCalibrationWidths &ca)
const float * getValues() const
get value for all capId = 0..3
Definition: HcalGain.h:18
const HcalPedestalWidths * mPedestalWidths
Definition: HcalDbService.h:73
const HcalRespCorrs * mRespCorrs
Definition: HcalDbService.h:80
const HcalCalibrationWidths & getHcalCalibrationWidths(const HcalGenericDetId &fId) const
HcalDbService(const edm::ParameterSet &)
const HcalZSThreshold * getHcalZSThreshold(const HcalGenericDetId &fId) const
const HcalElectronicsMap * mElectronicsMap
Definition: HcalDbService.h:79
const HcalGainWidths * mGainWidths
Definition: HcalDbService.h:75
const HcalTimeCorrs * mTimeCorrs
Definition: HcalDbService.h:83
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:97
const HcalCalibrationWidths & getCalibrationWidths(const DetId id) const
float getValue() const
Definition: HcalLUTCorr.h:18
HcalCalibrationWidthsSet mCalibWidthSet
Definition: HcalDbService.h:89
const HcalPFCorr * getHcalPFCorr(const HcalGenericDetId &fId) const
const HcalZSThresholds * mZSThresholds
Definition: HcalDbService.h:81
float getValue() const
Definition: HcalTimeCorr.h:18
const HcalQIEData * mQIEData
Definition: HcalDbService.h:77
const HcalCalibrations & getCalibrations(const DetId id) const
const HcalGain * getGain(const HcalGenericDetId &fId) const
float getWidth(int fCapId) const
get width (sqrt(sigma_i_i)) for capId = 0..3
const HcalQIECoder * getHcalCoder(const HcalGenericDetId &fId) const
void buildCalibrations() const
float getValue() const
Definition: HcalRespCorr.h:18
void setCalibrations(const DetId id, const HcalCalibrations &ca)
const HcalElectronicsMap * getHcalMapping() const
const float * getValues() const
get value for all capId = 0..3
Definition: HcalPedestal.h:17
const HcalPedestals * mPedestals
Definition: HcalDbService.h:72
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
const Item * getValues(DetId fId) const
const HcalPedestal * getPedestal(const HcalGenericDetId &fId) const
const HcalQIEShape * getHcalShape() const
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -&gt; fC conversion.
Definition: HcalQIECoder.cc:22