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 // $Id: CastorDbService.cc,v 1.4 2010/02/20 20:54:59 wmtan Exp $
6 
8 
12 
14 
15 #include <cmath>
16 
18  :
19  mQieShapeCache (0),
20  mPedestals (0),
21  mPedestalWidths (0),
22  mGains (0),
23  mGainWidths (0),
24  mQIEData(0),
25  mElectronicsMap(0)
26  {}
27 
28 bool CastorDbService::makeCastorCalibration (const HcalGenericDetId& fId, CastorCalibrations* fObject, bool pedestalInADC) const {
29  if (fObject) {
30  const CastorPedestal* pedestal = getPedestal (fId);
31  const CastorGain* gain = getGain (fId);
32 
33  if (pedestalInADC) {
34  const CastorQIEShape* shape=getCastorShape();
35  const CastorQIECoder* coder=getCastorCoder(fId);
36  if (pedestal && gain && shape && coder ) {
37  float pedTrue[4];
38  for (int i=0; i<4; i++) {
39  float x=pedestal->getValues()[i];
40  int x1=(int)std::floor(x);
41  int x2=(int)std::floor(x+1);
42  // y = (y2-y1)/(x2-x1) * (x - x1) + y1 [note: x2-x1=1]
43  float y2=coder->charge(*shape,x2,i);
44  float y1=coder->charge(*shape,x1,i);
45  pedTrue[i]=(y2-y1)*(x-x1)+y1;
46  }
47  *fObject = CastorCalibrations (gain->getValues (), pedTrue );
48  return true;
49  }
50  } else {
51  if (pedestal && gain ) {
52  *fObject = CastorCalibrations (gain->getValues (), pedestal->getValues () );
53  return true;
54  }
55  }
56  }
57  return false;
58 }
59 
61  // we use the set of ids for pedestals as the master list
62  if ((!mPedestals) || (!mGains) || (!mQIEData) ) return;
63  std::vector<DetId> ids=mPedestals->getAllChannels();
64  bool pedsInADC = mPedestals->isADC();
65  // clear the calibrations set
66  mCalibSet.clear();
67  // loop!
68  CastorCalibrations tool;
69 
70  // std::cout << " length of id-vector: " << ids.size() << std::endl;
71  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
72  // make
73  bool ok=makeCastorCalibration(*id,&tool,pedsInADC);
74  // store
75  if (ok) mCalibSet.setCalibrations(*id,tool);
76  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
77  }
78  mCalibSet.sort();
79 }
80 
82  // we use the set of ids for pedestal widths as the master list
83  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
84 
85  std::vector<DetId> ids=mPedestalWidths->getAllChannels();
86  bool pedsInADC = mPedestalWidths->isADC();
87  // clear the calibrations set
89  // loop!
91 
92  // std::cout << " length of id-vector: " << ids.size() << std::endl;
93  for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
94  // make
95  bool ok=makeCastorCalibrationWidth(*id,&tool,pedsInADC);
96  // store
97  if (ok) mCalibWidthSet.setCalibrationWidths(*id,tool);
98  // std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
99  }
101 }
102 
104  CastorCalibrationWidths* fObject, bool pedestalInADC) const {
105  if (fObject) {
106  const CastorPedestalWidth* pedestalwidth = getPedestalWidth (fId);
107  const CastorGainWidth* gainwidth = getGainWidth (fId);
108  if (pedestalInADC) {
109  const CastorQIEShape* shape=getCastorShape();
110  const CastorQIECoder* coder=getCastorCoder(fId);
111  if (pedestalwidth && gainwidth && shape && coder) {
112  float pedTrueWidth[4];
113  for (int i=0; i<4; i++) {
114  float x=pedestalwidth->getWidth(i);
115  // assume QIE is linear in low range and use x1=0 and x2=1
116  // y = (y2-y1) * (x) [do not add any constant, only scale!]
117  float y2=coder->charge(*shape,1,i);
118  float y1=coder->charge(*shape,0,i);
119  pedTrueWidth[i]=(y2-y1)*x;
120  }
121  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
122  return true;
123  }
124  } else {
125  if (pedestalwidth && gainwidth) {
126  float pedestalWidth [4];
127  for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
128  *fObject = CastorCalibrationWidths (gainwidth->getValues (), pedestalWidth);
129  return true;
130  }
131  }
132  }
133  return false;
134 }
135 
136 
138  if (mPedestals) {
139  return mPedestals->getValues (fId);
140  }
141  return 0;
142 }
143 
145  if (mPedestalWidths) {
146  return mPedestalWidths->getValues (fId);
147  }
148  return 0;
149 }
150 
152  if (mGains) {
153  return mGains->getValues(fId);
154  }
155  return 0;
156 }
157 
159  if (mGainWidths) {
160  return mGainWidths->getValues (fId);
161  }
162  return 0;
163 }
164 
166  if (mQIEData) {
167  return mQIEData->getCoder (fId);
168  }
169  return 0;
170 }
171 
173  if (mQIEData) {
174  return &mQIEData->getShape ();
175  }
176  return 0;
177 }
178 
180  return mElectronicsMap;
181 }
182 
184 {
185  return mChannelQuality->getValues (fId);
186 }
187 
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 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
const Item * getValues(DetId fId) const
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:97
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