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 
11 
13 
14 #include <cmath>
15 
17  mPedestals (0), mPedestalWidths (0),
18  mGains (0), mGainWidths (0),
19  mQIEData(0),
20  mElectronicsMap(0),
21  mRespCorrs(0),
22  mL1TriggerObjects(0),
23  mTimeCorrs(0),
24  mLUTCorrs(0),
25  mPFCorrs(0),
26  mLutMetadata(0),
27  mCalibSet(nullptr), mCalibWidthSet(nullptr)
28  {}
29 
31  if (mPedestals && mPedestals->topo()) return mPedestals->topo();
32  if (mGains && mGains->topo()) return mGains->topo();
33  if (mRespCorrs && mRespCorrs->topo()) return mRespCorrs->topo();
35  if (mLutMetadata && mLutMetadata->topo()) return mLutMetadata->topo();
36  return 0;
37 }
38 
39 
41 {
43  return (*mCalibSet.load(std::memory_order_acquire)).getCalibrations(fId);
44 }
45 
47 {
49  return (*mCalibWidthSet.load(std::memory_order_acquire)).getCalibrationWidths(fId);
50 }
51 
53  // we use the set of ids for pedestals as the master list
54  if ((!mPedestals) || (!mGains) || (!mQIEData) || (!mRespCorrs) || (!mTimeCorrs) || (!mLUTCorrs) ) return;
55 
56  if (!mCalibSet.load(std::memory_order_acquire)) {
57 
58  auto ptr = new HcalCalibrationsSet();
59 
60  std::vector<DetId> ids=mPedestals->getAllChannels();
61  bool pedsInADC = mPedestals->isADC();
62  // loop!
63  HcalCalibrations tool;
64 
65  // std::cout << " length of id-vector: " << ids.size() << std::endl;
66  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
67  // make
68  bool ok=makeHcalCalibration(*id,&tool,pedsInADC);
69  // store
70  if (ok) ptr->setCalibrations(*id,tool);
71  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
72  }
73  ptr->sort();
74 
75  HcalCalibrationsSet* expect = nullptr;
76  bool exchanged = mCalibSet.compare_exchange_strong(expect, ptr, std::memory_order_acq_rel);
77  if(!exchanged) {
78  delete ptr;
79  }
80  }
81 }
82 
84  // we use the set of ids for pedestal widths as the master list
85  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
86 
87  if (!mCalibWidthSet.load(std::memory_order_acquire)) {
88 
89  auto ptr = new HcalCalibrationWidthsSet();
90 
91  const std::vector<DetId>& ids=mPedestalWidths->getAllChannels();
92  bool pedsInADC = mPedestalWidths->isADC();
93  // loop!
95 
96  // std::cout << " length of id-vector: " << ids.size() << std::endl;
97  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
98  // make
99  bool ok=makeHcalCalibrationWidth(*id,&tool,pedsInADC);
100  // store
101  if (ok) ptr->setCalibrationWidths(*id,tool);
102  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
103  }
104  ptr->sort();
105 
106  HcalCalibrationWidthsSet* expect = nullptr;
107  bool exchanged = mCalibWidthSet.compare_exchange_strong(expect, ptr, std::memory_order_acq_rel);
108  if(!exchanged) {
109  delete ptr;
110  }
111  }
112 }
113 
114 bool HcalDbService::makeHcalCalibration (const HcalGenericDetId& fId, HcalCalibrations* fObject, bool pedestalInADC) const {
115  if (fObject) {
116  const HcalPedestal* pedestal = getPedestal (fId);
117  const HcalGain* gain = getGain (fId);
118  const HcalRespCorr* respcorr = getHcalRespCorr (fId);
119  const HcalTimeCorr* timecorr = getHcalTimeCorr (fId);
120  const HcalLUTCorr* lutcorr = getHcalLUTCorr (fId);
121 
122  if (pedestalInADC) {
123  const HcalQIECoder* coder=getHcalCoder(fId);
124  const HcalQIEShape* shape=getHcalShape(coder);
125  if (pedestal && gain && shape && coder && respcorr && timecorr && lutcorr) {
126  float pedTrue[4];
127  for (int i=0; i<4; i++) {
128  float x=pedestal->getValues()[i];
129  int x1=(int)std::floor(x);
130  int x2=(int)std::floor(x+1);
131  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
132  float y2=coder->charge(*shape,x2,i);
133  float y1=coder->charge(*shape,x1,i);
134  pedTrue[i]=(y2-y1)*(x-x1)+y1;
135  }
136  *fObject = HcalCalibrations (gain->getValues (), pedTrue, respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
137  return true;
138  }
139  } else {
140  if (pedestal && gain && respcorr && timecorr && lutcorr) {
141  *fObject = HcalCalibrations (gain->getValues (), pedestal->getValues (), respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
142  return true;
143  }
144  }
145  }
146  return false;
147 }
148 
150  HcalCalibrationWidths* fObject, bool pedestalInADC) const {
151  if (fObject) {
152  const HcalPedestalWidth* pedestalwidth = getPedestalWidth (fId);
153  const HcalGainWidth* gainwidth = getGainWidth (fId);
154  if (pedestalInADC) {
155  const HcalQIECoder* coder=getHcalCoder(fId);
156  const HcalQIEShape* shape=getHcalShape(coder);
157  if (pedestalwidth && gainwidth && shape && coder) {
158  float pedTrueWidth[4];
159  for (int i=0; i<4; i++) {
160  float x=pedestalwidth->getWidth(i);
161  // assume QIE is linear in low range and use x1=0 and x2=1
162  // y = (y2-y1) * (x) [do not add any constant, only scale!]
163  float y2=coder->charge(*shape,1,i);
164  float y1=coder->charge(*shape,0,i);
165  pedTrueWidth[i]=(y2-y1)*x;
166  }
167  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
168  return true;
169  }
170  } else {
171  if (pedestalwidth && gainwidth) {
172  float pedestalWidth [4];
173  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
174  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedestalWidth);
175  return true;
176  }
177  }
178  }
179  return false;
180 }
181 
182 
184  if (mRespCorrs) {
185  return mRespCorrs->getValues (fId);
186  }
187  return 0;
188 }
189 
191  if (mPedestals) {
192  return mPedestals->getValues (fId);
193  }
194  return 0;
195 }
196 
198  if (mPedestalWidths) {
199  return mPedestalWidths->getValues (fId);
200  }
201  return 0;
202 }
203 
205  if (mGains) {
206  return mGains->getValues(fId);
207  }
208  return 0;
209 }
210 
212  if (mGainWidths) {
213  return mGainWidths->getValues (fId);
214  }
215  return 0;
216 }
217 
219  if (mQIEData) {
220  return mQIEData->getCoder (fId);
221  }
222  return 0;
223 }
224 
226  if (mQIEData) {
227  return &mQIEData->getShape (fId);
228  }
229  return 0;
230 }
231 
233  if (mQIEData) {
234  return &mQIEData->getShape(coder);
235  }
236  return 0;
237 }
238 
239 
241  return mElectronicsMap;
242 }
243 
245 {
246  return mL1TriggerObjects->getValues (fId);
247 }
248 
250 {
251  return mChannelQuality->getValues (fId);
252 }
253 
255 {
256  return mZSThresholds->getValues (fId);
257 }
258 
260  if (mTimeCorrs) {
261  return mTimeCorrs->getValues (fId);
262  }
263  return 0;
264 }
265 
267  if (mLUTCorrs) {
268  return mLUTCorrs->getValues (fId);
269  }
270  return 0;
271 }
272 
274  if (mPFCorrs) {
275  return mPFCorrs->getValues (fId);
276  }
277  return 0;
278 }
279 
281  return mLutMetadata;
282 }
283 
const HcalLUTCorrs * mLUTCorrs
Definition: HcalDbService.h:88
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
int i
Definition: DBlmapReader.cc:9
bool isADC() const
Definition: HcalPedestals.h:28
const HcalQIECoder * getCoder(DetId fId) const
get QIE parameters
Definition: HcalQIEData.h:40
const HcalChannelStatus * getHcalChannelStatus(const HcalGenericDetId &fId) const
#define nullptr
const HcalL1TriggerObjects * mL1TriggerObjects
Definition: HcalDbService.h:86
const HcalPFCorrs * mPFCorrs
Definition: HcalDbService.h:89
const Item * getValues(DetId fId, bool throwOnFail=true) const
const HcalGains * mGains
Definition: HcalDbService.h:79
const HcalChannelQuality * mChannelQuality
Definition: HcalDbService.h:82
const HcalLutMetadata * mLutMetadata
Definition: HcalDbService.h:90
const float * getValues() const
get value for all capId = 0..3
Definition: HcalGainWidth.h:19
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
std::vector< DetId > getAllChannels() const
const HcalL1TriggerObject * getHcalL1TriggerObject(const HcalGenericDetId &fId) const
bool makeHcalCalibration(const HcalGenericDetId &fId, HcalCalibrations *fObject, bool pedestalInADC) const
const HcalQIEShape & getShape(DetId fId) const
get basic shape
Definition: HcalQIEData.h:37
const HcalLutMetadata * getHcalLutMetadata() const
bool makeHcalCalibrationWidth(const HcalGenericDetId &fId, HcalCalibrationWidths *fObject, bool pedestalInADC) const
const float * getValues() const
get value for all capId = 0..3
Definition: HcalGain.h:20
const HcalPedestalWidths * mPedestalWidths
Definition: HcalDbService.h:78
const HcalRespCorrs * mRespCorrs
Definition: HcalDbService.h:84
const HcalCalibrationWidths & getHcalCalibrationWidths(const HcalGenericDetId &fId) const
HcalDbService(const edm::ParameterSet &)
const HcalZSThreshold * getHcalZSThreshold(const HcalGenericDetId &fId) const
const HcalElectronicsMap * mElectronicsMap
Definition: HcalDbService.h:83
const HcalGainWidths * mGainWidths
Definition: HcalDbService.h:80
const HcalTimeCorrs * mTimeCorrs
Definition: HcalDbService.h:87
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:96
float getValue() const
Definition: HcalLUTCorr.h:20
const HcalPFCorr * getHcalPFCorr(const HcalGenericDetId &fId) const
const HcalTopology * getTopologyUsed() const
const HcalZSThresholds * mZSThresholds
Definition: HcalDbService.h:85
float getValue() const
Definition: HcalTimeCorr.h:20
const HcalQIEData * mQIEData
Definition: HcalDbService.h:81
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
const HcalQIEShape * getHcalShape(const HcalGenericDetId &fId) const
void buildCalibrations() const
float getValue() const
Definition: HcalRespCorr.h:20
const HcalElectronicsMap * getHcalMapping() const
const float * getValues() const
get value for all capId = 0..3
Definition: HcalPedestal.h:19
const HcalPedestals * mPedestals
Definition: HcalDbService.h:77
Definition: DDAxes.h:10
std::atomic< HcalCalibrationsSet * > mCalibSet
Definition: HcalDbService.h:92
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
std::atomic< HcalCalibrationWidthsSet * > mCalibWidthSet
Definition: HcalDbService.h:93
const HcalPedestal * getPedestal(const HcalGenericDetId &fId) const
const HcalTopology * topo() const
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -&gt; fC conversion.
Definition: HcalQIECoder.cc:22