CMS 3D CMS Logo

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 (nullptr), mPedestalWidths (nullptr),
18  mGains (nullptr), mGainWidths (nullptr),
19  mQIEData(nullptr),
20  mQIETypes(nullptr),
21  mElectronicsMap(nullptr), mFrontEndMap(nullptr),
22  mRespCorrs(nullptr),
23  mL1TriggerObjects(nullptr),
24  mTimeCorrs(nullptr),
25  mLUTCorrs(nullptr),
26  mPFCorrs(nullptr),
27  mLutMetadata(nullptr),
28  mSiPMParameters(nullptr), mSiPMCharacteristics(nullptr),
29  mTPChannelParameters(nullptr), mTPParameters(nullptr),
30  mMCParams(nullptr),
31  mCalibSet(nullptr), mCalibWidthSet(nullptr)
32  {}
33 
35  delete mCalibSet.load();
36  delete mCalibWidthSet.load();
37 }
38 
40  if (mPedestals && mPedestals->topo()) return mPedestals->topo();
41  if (mGains && mGains->topo()) return mGains->topo();
42  if (mRespCorrs && mRespCorrs->topo()) return mRespCorrs->topo();
43  if (mQIETypes && mQIETypes->topo()) return mQIETypes->topo();
45  if (mLutMetadata && mLutMetadata->topo()) return mLutMetadata->topo();
46  return nullptr;
47 }
48 
49 
51 {
53  return (*mCalibSet.load(std::memory_order_acquire)).getCalibrations(fId);
54 }
55 
57 {
59  return (*mCalibWidthSet.load(std::memory_order_acquire)).getCalibrationWidths(fId);
60 }
61 
63 {
65  return mCalibSet.load(std::memory_order_acquire);
66 }
67 
69 {
71  return mCalibWidthSet.load(std::memory_order_acquire);
72 }
73 
75  // we use the set of ids for pedestals as the master list
76  if ((!mPedestals) || (!mGains) || (!mQIEData) || (!mQIETypes) || (!mRespCorrs) || (!mTimeCorrs) || (!mLUTCorrs) ) return;
77 
78  if (!mCalibSet.load(std::memory_order_acquire)) {
79 
80  auto ptr = new HcalCalibrationsSet();
81 
82  std::vector<DetId> ids=mPedestals->getAllChannels();
83  bool pedsInADC = mPedestals->isADC();
84  // loop!
85  HcalCalibrations tool;
86 
87  // std::cout << " length of id-vector: " << ids.size() << std::endl;
88  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
89  // make
90  bool ok=makeHcalCalibration(*id,&tool,pedsInADC);
91  // store
92  if (ok) ptr->setCalibrations(*id,tool);
93  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
94  }
95  HcalCalibrationsSet const * cptr = ptr;
96  HcalCalibrationsSet const * expect = nullptr;
97  bool exchanged = mCalibSet.compare_exchange_strong(expect, cptr, std::memory_order_acq_rel);
98  if(!exchanged) {
99  delete ptr;
100  }
101  }
102 }
103 
105  // we use the set of ids for pedestal widths as the master list
106  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
107 
108  if (!mCalibWidthSet.load(std::memory_order_acquire)) {
109 
110  auto ptr = new HcalCalibrationWidthsSet();
111 
112  const std::vector<DetId>& ids=mPedestalWidths->getAllChannels();
113  bool pedsInADC = mPedestalWidths->isADC();
114  // loop!
116 
117  // std::cout << " length of id-vector: " << ids.size() << std::endl;
118  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
119  // make
120  bool ok=makeHcalCalibrationWidth(*id,&tool,pedsInADC);
121  // store
122  if (ok) ptr->setCalibrationWidths(*id,tool);
123  // std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
124  }
125  HcalCalibrationWidthsSet const * cptr = ptr;
126  HcalCalibrationWidthsSet const * expect = nullptr;
127  bool exchanged = mCalibWidthSet.compare_exchange_strong(expect, cptr, std::memory_order_acq_rel);
128  if(!exchanged) {
129  delete ptr;
130  }
131  }
132 }
133 
134 bool HcalDbService::makeHcalCalibration (const HcalGenericDetId& fId, HcalCalibrations* fObject, bool pedestalInADC) const {
135  if (fObject) {
136  const HcalPedestal* pedestal = getPedestal (fId);
137  const HcalGain* gain = getGain (fId);
138  const HcalRespCorr* respcorr = getHcalRespCorr (fId);
139  const HcalTimeCorr* timecorr = getHcalTimeCorr (fId);
140  const HcalLUTCorr* lutcorr = getHcalLUTCorr (fId);
141 
142  if (pedestalInADC) {
143  const HcalQIECoder* coder=getHcalCoder(fId);
144  const HcalQIEShape* shape=getHcalShape(coder);
145  if (pedestal && gain && shape && coder && respcorr && timecorr && lutcorr) {
146  float pedTrue[4];
147  for (int i=0; i<4; i++) {
148  float x=pedestal->getValues()[i];
149  int x1=(int)std::floor(x);
150  int x2=(int)std::floor(x+1);
151  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
152  float y2=coder->charge(*shape,x2,i);
153  float y1=coder->charge(*shape,x1,i);
154  pedTrue[i]=(y2-y1)*(x-x1)+y1;
155  }
156  *fObject = HcalCalibrations (gain->getValues (), pedTrue, respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
157  return true;
158  }
159  } else {
160  if (pedestal && gain && respcorr && timecorr && lutcorr) {
161  *fObject = HcalCalibrations (gain->getValues (), pedestal->getValues (), respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
162  return true;
163  }
164  }
165  }
166  return false;
167 }
168 
170  HcalCalibrationWidths* fObject, bool pedestalInADC) const {
171  if (fObject) {
172  const HcalPedestalWidth* pedestalwidth = getPedestalWidth (fId);
173  const HcalGainWidth* gainwidth = getGainWidth (fId);
174  if (pedestalInADC) {
175  const HcalQIECoder* coder=getHcalCoder(fId);
176  const HcalQIEShape* shape=getHcalShape(coder);
177  if (pedestalwidth && gainwidth && shape && coder) {
178  float pedTrueWidth[4];
179  for (int i=0; i<4; i++) {
180  float x=pedestalwidth->getWidth(i);
181  // assume QIE is linear in low range and use x1=0 and x2=1
182  // y = (y2-y1) * (x) [do not add any constant, only scale!]
183  float y2=coder->charge(*shape,1,i);
184  float y1=coder->charge(*shape,0,i);
185  pedTrueWidth[i]=(y2-y1)*x;
186  }
187  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
188  return true;
189  }
190  } else {
191  if (pedestalwidth && gainwidth) {
192  float pedestalWidth [4];
193  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
194  *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedestalWidth);
195  return true;
196  }
197  }
198  }
199  return false;
200 }
201 
203  if (mQIETypes) {
204  return mQIETypes->getValues (fId);
205  }
206  return nullptr;
207 }
208 
210  if (mRespCorrs) {
211  return mRespCorrs->getValues (fId);
212  }
213  return nullptr;
214 }
215 
217  if (mPedestals) {
218  return mPedestals->getValues (fId);
219  }
220  return nullptr;
221 }
222 
224  if (mPedestalWidths) {
225  return mPedestalWidths->getValues (fId);
226  }
227  return nullptr;
228 }
229 
231  if (mGains) {
232  return mGains->getValues(fId);
233  }
234  return nullptr;
235 }
236 
238  if (mGainWidths) {
239  return mGainWidths->getValues (fId);
240  }
241  return nullptr;
242 }
243 
245  if (mQIEData) {
246  return mQIEData->getCoder (fId);
247  }
248  return nullptr;
249 }
250 
252  if (mQIEData && mQIETypes) {
253  //currently 3 types of QIEs exist: QIE8, QIE10, QIE11
254  int qieType = mQIETypes->getValues(fId)->getValue();
255  //QIE10 and QIE11 have same shape (ADC ladder)
256  if(qieType>0) qieType = 1;
257  return &mQIEData->getShape(qieType);
258  }
259  return nullptr;
260 }
261 
263  HcalGenericDetId fId(coder->rawId());
264  return getHcalShape(fId);
265 }
266 
268  return mElectronicsMap;
269 }
270 
272  return mFrontEndMap;
273 }
274 
276 {
277  return mL1TriggerObjects->getValues (fId);
278 }
279 
281 {
282  return mChannelQuality->getValues (fId);
283 }
284 
286 {
287  return mZSThresholds->getValues (fId);
288 }
289 
291  if (mTimeCorrs) {
292  return mTimeCorrs->getValues (fId);
293  }
294  return nullptr;
295 }
296 
298  if (mLUTCorrs) {
299  return mLUTCorrs->getValues (fId);
300  }
301  return nullptr;
302 }
303 
305  if (mPFCorrs) {
306  return mPFCorrs->getValues (fId);
307  }
308  return nullptr;
309 }
310 
312  return mLutMetadata;
313 }
314 
316  if (mSiPMParameters) {
317  return mSiPMParameters->getValues (fId);
318  }
319  return nullptr;
320 }
321 
323  return mSiPMCharacteristics;
324 }
325 
327  if (mTPChannelParameters) {
328  return mTPChannelParameters->getValues (fId);
329  }
330  return nullptr;
331 }
332 
334  if (mMCParams) {
335  return mMCParams->getValues (fId);
336  }
337  return nullptr;
338 }
339 
341  return mTPParameters;
342 }
343 
const HcalLUTCorrs * mLUTCorrs
const HcalQIETypes * mQIETypes
Definition: HcalDbService.h:97
const HcalQIEShape & getShape(int qieType) const
get basic shape
Definition: HcalQIEData.h:37
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
const HcalFrontEndMap * getHcalFrontEndMapping() const
const HcalSiPMParameters * mSiPMParameters
int getValue() const
Definition: HcalQIEType.h:20
bool isADC() const
Definition: HcalPedestals.h:28
const HcalMCParams * mMCParams
const HcalTPParameters * mTPParameters
std::atomic< HcalCalibrationWidthsSet const * > mCalibWidthSet
const HcalTPChannelParameter * getHcalTPChannelParameter(const HcalGenericDetId &fId) const
const HcalQIECoder * getCoder(DetId fId) const
get QIE parameters
Definition: HcalQIEData.h:39
const HcalChannelStatus * getHcalChannelStatus(const HcalGenericDetId &fId) const
const HcalL1TriggerObjects * mL1TriggerObjects
const HcalPFCorrs * mPFCorrs
const Item * getValues(DetId fId, bool throwOnFail=true) const
const HcalGains * mGains
Definition: HcalDbService.h:94
#define nullptr
const HcalChannelQuality * mChannelQuality
Definition: HcalDbService.h:98
const HcalLutMetadata * mLutMetadata
const HcalTPChannelParameters * mTPChannelParameters
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
const HcalMCParam * getHcalMCParam(const HcalGenericDetId &fId) 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:93
const HcalRespCorrs * mRespCorrs
const HcalCalibrationWidths & getHcalCalibrationWidths(const HcalGenericDetId &fId) const
HcalDbService(const edm::ParameterSet &)
const HcalZSThreshold * getHcalZSThreshold(const HcalGenericDetId &fId) const
const HcalElectronicsMap * mElectronicsMap
Definition: HcalDbService.h:99
const HcalGainWidths * mGainWidths
Definition: HcalDbService.h:95
const HcalCalibrationsSet * getHcalCalibrationsSet() const
const HcalTimeCorrs * mTimeCorrs
#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
float getValue() const
Definition: HcalTimeCorr.h:20
const HcalQIEData * mQIEData
Definition: HcalDbService.h:96
const HcalGain * getGain(const HcalGenericDetId &fId) const
const HcalTPParameters * getHcalTPParameters() const
float getWidth(int fCapId) const
get width (sqrt(sigma_i_i)) for capId = 0..3
const HcalSiPMCharacteristics * getHcalSiPMCharacteristics() const
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
const float * getValues() const
get value for all capId = 0..3
Definition: HcalPedestal.h:19
const HcalFrontEndMap * mFrontEndMap
const HcalPedestals * mPedestals
Definition: HcalDbService.h:92
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
const HcalPedestal * getPedestal(const HcalGenericDetId &fId) const
const HcalSiPMParameter * getHcalSiPMParameter(const HcalGenericDetId &fId) const
const HcalSiPMCharacteristics * mSiPMCharacteristics
const HcalCalibrationWidthsSet * getHcalCalibrationWidthsSet() const
const HcalTopology * topo() const
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -> fC conversion.
Definition: HcalQIECoder.cc:22