CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  mQieShapeCache (0),
19  mPedestals (0),
20  mPedestalWidths (0),
21  mGains (0),
22  mGainWidths (0),
23  mQIEData(0),
24  mElectronicsMap(0)
25  {}
26 
27 bool CastorDbService::makeCastorCalibration (const HcalGenericDetId& fId, CastorCalibrations* fObject, bool pedestalInADC) const {
28  if (fObject) {
29  const CastorPedestal* pedestal = getPedestal (fId);
30  const CastorGain* gain = getGain (fId);
31 
32  if (pedestalInADC) {
33  const CastorQIEShape* shape=getCastorShape();
34  const CastorQIECoder* coder=getCastorCoder(fId);
35  if (pedestal && gain && shape && coder ) {
36  float pedTrue[4];
37  for (int i=0; i<4; i++) {
38  float x=pedestal->getValues()[i];
39  int x1=(int)std::floor(x);
40  int x2=(int)std::floor(x+1);
41  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
42  float y2=coder->charge(*shape,x2,i);
43  float y1=coder->charge(*shape,x1,i);
44  pedTrue[i]=(y2-y1)*(x-x1)+y1;
45  }
46  *fObject = CastorCalibrations (gain->getValues (), pedTrue );
47  return true;
48  }
49  } else {
50  if (pedestal && gain ) {
51  *fObject = CastorCalibrations (gain->getValues (), pedestal->getValues () );
52  return true;
53  }
54  }
55  }
56  return false;
57 }
58 
60  // we use the set of ids for pedestals as the master list
61  if ((!mPedestals) || (!mGains) || (!mQIEData) ) return;
62  std::vector<DetId> ids=mPedestals->getAllChannels();
63  bool pedsInADC = mPedestals->isADC();
64  // clear the calibrations set
65  mCalibSet.clear();
66  // loop!
67  CastorCalibrations tool;
68 
69  // std::cout << " length of id-vector: " << ids.size() << std::endl;
70  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
71  // make
72  bool ok=makeCastorCalibration(*id,&tool,pedsInADC);
73  // store
74  if (ok) mCalibSet.setCalibrations(*id,tool);
75  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
76  }
77  mCalibSet.sort();
78 }
79 
81  // we use the set of ids for pedestal widths as the master list
82  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
83 
84  std::vector<DetId> ids=mPedestalWidths->getAllChannels();
85  bool pedsInADC = mPedestalWidths->isADC();
86  // clear the calibrations set
88  // loop!
90 
91  // std::cout << " length of id-vector: " << ids.size() << std::endl;
92  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
93  // make
94  bool ok=makeCastorCalibrationWidth(*id,&tool,pedsInADC);
95  // store
96  if (ok) mCalibWidthSet.setCalibrationWidths(*id,tool);
97  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
98  }
100 }
101 
103  CastorCalibrationWidths* fObject, bool pedestalInADC) const {
104  if (fObject) {
105  const CastorPedestalWidth* pedestalwidth = getPedestalWidth (fId);
106  const CastorGainWidth* gainwidth = getGainWidth (fId);
107  if (pedestalInADC) {
108  const CastorQIEShape* shape=getCastorShape();
109  const CastorQIECoder* coder=getCastorCoder(fId);
110  if (pedestalwidth && gainwidth && shape && coder) {
111  float pedTrueWidth[4];
112  for (int i=0; i<4; i++) {
113  float x=pedestalwidth->getWidth(i);
114  // assume QIE is linear in low range and use x1=0 and x2=1
115  // y = (y2-y1) * (x) [do not add any constant, only scale!]
116  float y2=coder->charge(*shape,1,i);
117  float y1=coder->charge(*shape,0,i);
118  pedTrueWidth[i]=(y2-y1)*x;
119  }
120  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
121  return true;
122  }
123  } else {
124  if (pedestalwidth && gainwidth) {
125  float pedestalWidth [4];
126  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
127  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedestalWidth);
128  return true;
129  }
130  }
131  }
132  return false;
133 }
134 
135 
137  if (mPedestals) {
138  return mPedestals->getValues (fId);
139  }
140  return 0;
141 }
142 
144  if (mPedestalWidths) {
145  return mPedestalWidths->getValues (fId);
146  }
147  return 0;
148 }
149 
151  if (mGains) {
152  return mGains->getValues(fId);
153  }
154  return 0;
155 }
156 
158  if (mGainWidths) {
159  return mGainWidths->getValues (fId);
160  }
161  return 0;
162 }
163 
165  if (mQIEData) {
166  return mQIEData->getCoder (fId);
167  }
168  return 0;
169 }
170 
172  if (mQIEData) {
173  return &mQIEData->getShape ();
174  }
175  return 0;
176 }
177 
179  return mElectronicsMap;
180 }
181 
183 {
184  return mChannelQuality->getValues (fId);
185 }
186 
int i
Definition: DBlmapReader.cc:9
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:37
const CastorChannelStatus * getCastorChannelStatus(const HcalGenericDetId &fId) const
const CastorQIEData * mQIEData
std::vector< DetId > getAllChannels() const
const CastorPedestals * mPedestals
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] -&gt; 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:39
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
Definition: DDAxes.h:10
const float * getValues() const
get value for all capId = 0..3
Definition: CastorGain.h:15
const CastorQIECoder * getCastorCoder(const HcalGenericDetId &fId) const