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