CMS 3D CMS Logo

CastorDbService.cc
Go to the documentation of this file.
1 //
2 // F.Ratnikov (UMd), Aug. 9, 2005
3 // Adapted for CASTOR by L. Mundim
4 //
5 
7 
11 
13 
14 #include <cmath>
15 
17  :
18  mPedestals (nullptr),
19  mPedestalWidths (nullptr),
20  mGains (nullptr),
21  mGainWidths (nullptr),
22  mQIEData(nullptr),
23  mElectronicsMap(nullptr)
24  {}
25 
26 bool CastorDbService::makeCastorCalibration (const HcalGenericDetId& fId, CastorCalibrations* fObject, bool pedestalInADC) const {
27  if (fObject) {
28  const CastorPedestal* pedestal = getPedestal (fId);
29  const CastorGain* gain = getGain (fId);
30 
31  if (pedestalInADC) {
32  const CastorQIEShape* shape=getCastorShape();
33  const CastorQIECoder* coder=getCastorCoder(fId);
34  if (pedestal && gain && shape && coder ) {
35  float pedTrue[4];
36  for (int i=0; i<4; i++) {
37  float x=pedestal->getValues()[i];
38  int x1=(int)std::floor(x);
39  int x2=(int)std::floor(x+1);
40  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
41  float y2=coder->charge(*shape,x2,i);
42  float y1=coder->charge(*shape,x1,i);
43  pedTrue[i]=(y2-y1)*(x-x1)+y1;
44  }
45  *fObject = CastorCalibrations (gain->getValues (), pedTrue );
46  return true;
47  }
48  } else {
49  if (pedestal && gain ) {
50  *fObject = CastorCalibrations (gain->getValues (), pedestal->getValues () );
51  return true;
52  }
53  }
54  }
55  return false;
56 }
57 
59  // we use the set of ids for pedestals as the master list
60  if ((!mPedestals) || (!mGains) || (!mQIEData) ) return;
61  std::vector<DetId> ids=mPedestals->getAllChannels();
62  bool pedsInADC = mPedestals->isADC();
63  // clear the calibrations set
64  mCalibSet.clear();
65  // loop!
66  CastorCalibrations tool;
67 
68  // std::cout << " length of id-vector: " << ids.size() << std::endl;
69  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
70  // make
71  bool ok=makeCastorCalibration(*id,&tool,pedsInADC);
72  // store
73  if (ok) mCalibSet.setCalibrations(*id,tool);
74  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
75  }
76  mCalibSet.sort();
77 }
78 
80  // we use the set of ids for pedestal widths as the master list
81  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
82 
83  std::vector<DetId> ids=mPedestalWidths->getAllChannels();
84  bool pedsInADC = mPedestalWidths->isADC();
85  // clear the calibrations set
87  // loop!
89 
90  // std::cout << " length of id-vector: " << ids.size() << std::endl;
91  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
92  // make
93  bool ok=makeCastorCalibrationWidth(*id,&tool,pedsInADC);
94  // store
95  if (ok) mCalibWidthSet.setCalibrationWidths(*id,tool);
96  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
97  }
99 }
100 
102  CastorCalibrationWidths* fObject, bool pedestalInADC) const {
103  if (fObject) {
104  const CastorPedestalWidth* pedestalwidth = getPedestalWidth (fId);
105  const CastorGainWidth* gainwidth = getGainWidth (fId);
106  if (pedestalInADC) {
107  const CastorQIEShape* shape=getCastorShape();
108  const CastorQIECoder* coder=getCastorCoder(fId);
109  if (pedestalwidth && gainwidth && shape && coder) {
110  float pedTrueWidth[4];
111  for (int i=0; i<4; i++) {
112  float x=pedestalwidth->getWidth(i);
113  // assume QIE is linear in low range and use x1=0 and x2=1
114  // y = (y2-y1) * (x) [do not add any constant, only scale!]
115  float y2=coder->charge(*shape,1,i);
116  float y1=coder->charge(*shape,0,i);
117  pedTrueWidth[i]=(y2-y1)*x;
118  }
119  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
120  return true;
121  }
122  } else {
123  if (pedestalwidth && gainwidth) {
124  float pedestalWidth [4];
125  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
126  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedestalWidth);
127  return true;
128  }
129  }
130  }
131  return false;
132 }
133 
134 
136  if (mPedestals) {
137  return mPedestals->getValues (fId);
138  }
139  return nullptr;
140 }
141 
143  if (mPedestalWidths) {
144  return mPedestalWidths->getValues (fId);
145  }
146  return nullptr;
147 }
148 
150  if (mGains) {
151  return mGains->getValues(fId);
152  }
153  return nullptr;
154 }
155 
157  if (mGainWidths) {
158  return mGainWidths->getValues (fId);
159  }
160  return nullptr;
161 }
162 
164  if (mQIEData) {
165  return mQIEData->getCoder (fId);
166  }
167  return nullptr;
168 }
169 
171  if (mQIEData) {
172  return &mQIEData->getShape ();
173  }
174  return nullptr;
175 }
176 
178  return mElectronicsMap;
179 }
180 
182 {
183  return mChannelQuality->getValues (fId);
184 }
185 
bool makeCastorCalibration(const HcalGenericDetId &fId, CastorCalibrations *fObject, bool pedestalInADC) const
const CastorPedestal * getPedestal(const HcalGenericDetId &fId) const
void setCalibrationWidths(const DetId id, const CastorCalibrationWidths &ca)
const CastorQIEShape & getShape() const
get basic shape
Definition: CastorQIEData.h:36
const CastorChannelStatus * getCastorChannelStatus(const HcalGenericDetId &fId) const
const CastorQIEData * mQIEData
std::vector< DetId > getAllChannels() const
const CastorPedestals * mPedestals
#define nullptr
const Item * getValues(DetId fId, bool throwOnFail=true) const
const CastorGainWidths * mGainWidths
float getWidth(int fCapId) const
get width (sqrt(sigma_i_i)) for capId = 0..3
CastorDbService(const edm::ParameterSet &)
bool makeCastorCalibrationWidth(const HcalGenericDetId &fId, CastorCalibrationWidths *fObject, bool pedestalInADC) const
const CastorElectronicsMap * mElectronicsMap
CastorCalibrationsSet mCalibSet
const CastorQIEShape * getCastorShape() const
void setCalibrations(const DetId id, const CastorCalibrations &ca)
float charge(const CastorQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -> fC conversion.
const CastorGain * getGain(const HcalGenericDetId &fId) const
const float * getValues() const
get value for all capId = 0..3
const CastorElectronicsMap * getCastorMapping() const
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:96
const CastorQIECoder * getCoder(DetId fId) const
get QIE parameters
Definition: CastorQIEData.h:38
bool isADC() const
const CastorPedestalWidths * mPedestalWidths
const CastorGainWidth * getGainWidth(const HcalGenericDetId &fId) const
const CastorPedestalWidth * getPedestalWidth(const HcalGenericDetId &fId) const
const CastorChannelQuality * mChannelQuality
const CastorGains * mGains
const float * getValues() const
get value for all capId = 0..3
CastorCalibrationWidthsSet mCalibWidthSet
const float * getValues() const
get value for all capId = 0..3
Definition: CastorGain.h:17
const CastorQIECoder * getCastorCoder(const HcalGenericDetId &fId) const