CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCConditions.cc
Go to the documentation of this file.
20 
22 : theGains(),
23  theCrosstalk(),
24  thePedestals(),
25  theNoiseMatrix(),
26  theBadStrips(),
27  theBadWires(),
28  theBadChambers(),
29  theChipCorrections(),
30  theChamberTimingCorrections(),
31  readBadChannels_(false), readBadChambers_(false),useTimingCorrections_(false),
32  theAverageGain( -1.0 )
33 {
34  readBadChannels_ = ps.getParameter<bool>("readBadChannels");
35  readBadChambers_ = ps.getParameter<bool>("readBadChambers");
36  useTimingCorrections_ = ps.getParameter<bool>("CSCUseTimingCorrections");
37  // set size to hold all layers, using enum defined in .h
38  badStripWords.resize( MAX_LAYERS, 0 );
39  badWireWords.resize( MAX_LAYERS, 0 );
40 }
41 
42 
44 {
45 }
46 
48 {
49  // Strip gains
50  es.get<CSCDBGainsRcd>().get( theGains );
51  // Strip X-talk
52  es.get<CSCDBCrosstalkRcd>().get( theCrosstalk );
53  // Strip pedestals
54  es.get<CSCDBPedestalsRcd>().get( thePedestals );
55  // Strip autocorrelation noise matrix
57 
58  if ( useTimingCorrections()){
59  // Buckeye chip speeds
61  // Cable lengths from chambers to peripheral crate and additional chamber level timing correction
63  }
64 
65  if ( readBadChannels() ) {
66  // Bad strip channels
67  es.get<CSCBadStripsRcd>().get( theBadStrips );
68  // Bad wiregroup channels
69  es.get<CSCBadWiresRcd>().get( theBadWires );
70 
71  //@@ if( badStripsWatcher_.check( es ) ) {
73  //@@ }
74  //@@ if( badWiresWatcher_.check( es ) ) {
76  //@ }
77 
78  }
79 
80  // Has GainsRcd changed?
81  if( gainsWatcher_.check( es ) ) { // Yes...
82  theAverageGain = -1.0; // ...reset, so next access will recalculate it
83  }
84 
85  if ( readBadChambers() ) {
86  // Entire bad chambers
88  }
89 
90 // print();
91 }
92 
94  // reset existing values
95  badStripWords.assign( MAX_LAYERS, 0 );
96  if ( readBadChannels() ) {
97  // unpack what we've read from theBadStrips
98 
99  // chambers is a vector<BadChamber>
100  // channels is a vector<BadChannel>
101  // Each BadChamber contains its index (1-468 or 540 w. ME42), the no. of bad channels,
102  // and the index within vector<BadChannel> where this chamber's bad channels start.
103 
104  CSCIndexer indexer;
105 
106  for ( size_t i=0; i<theBadStrips->chambers.size(); ++i ) { // loop over bad chambers
107  int indexc = theBadStrips->chambers[i].chamber_index;
108  int start = theBadStrips->chambers[i].pointer; // where this chamber's bad channels start in vector<BadChannel>
109  int nbad = theBadStrips->chambers[i].bad_channels;
110 
111  CSCDetId id = indexer.detIdFromChamberIndex( indexc ); // We need this to build layer index (1-2808)
112 
113  for ( int j=start-1; j<start-1+nbad; ++j ) { // bad channels in this chamber
114  short lay = theBadStrips->channels[j].layer; // value 1-6
115  short chan = theBadStrips->channels[j].channel; // value 1-80
116  // short f1 = theBadStrips->channels[j].flag1;
117  // short f2 = theBadStrips->channels[j].flag2;
118  // short f3 = theBadStrips->channels[j].flag3;
119  int indexl = indexer.layerIndex( id.endcap(), id.station(), id.ring(), id.chamber(), lay );
120  badStripWords[indexl-1].set( chan-1, 1 ); // set bit 0-79 in 80-bit bitset representing this layer
121  } // j
122  } // i
123 
124  }
125 }
126 
128  // reset existing values
129  badWireWords.assign( MAX_LAYERS, 0 );
130  if ( readBadChannels() ) {
131  // unpack what we've read from theBadWires
132  CSCIndexer indexer;
133 
134  for ( size_t i=0; i<theBadWires->chambers.size(); ++i ) { // loop over bad chambers
135  int indexc = theBadWires->chambers[i].chamber_index;
136  int start = theBadWires->chambers[i].pointer; // where this chamber's bad channels start in vector<BadChannel>
137  int nbad = theBadWires->chambers[i].bad_channels;
138 
139  CSCDetId id = indexer.detIdFromChamberIndex( indexc ); // We need this to build layer index (1-2808)
140 
141  for ( int j=start-1; j<start-1+nbad; ++j ) { // bad channels in this chamber
142  short lay = theBadWires->channels[j].layer; // value 1-6
143  short chan = theBadWires->channels[j].channel; // value 1-80
144  // short f1 = theBadWires->channels[j].flag1;
145  // short f2 = theBadWires->channels[j].flag2;
146  // short f3 = theBadWires->channels[j].flag3;
147  int indexl = indexer.layerIndex( id.endcap(), id.station(), id.ring(), id.chamber(), lay );
148  badWireWords[indexl-1].set( chan-1, 1 ); // set bit 0-111 in 112-bit bitset representing this layer
149  } // j
150  } // i
151 
152  }
153 }
154 
155 bool CSCConditions::isInBadChamber( const CSCDetId& id ) const {
156  if ( readBadChambers() ) return theBadChambers->isInBadChamber( id );
157  else return false;
158 }
159 
161  //@@ NEEDS THOROUGH UPDATING
162 {
163 /*
164  std::cout << "SIZES: GAINS: " << theGains->gains.size()
165  << " PEDESTALS: " << thePedestals->pedestals.size()
166  << " NOISES " << theNoiseMatrix->matrix.size() << std::endl;;
167 
168  std::map< int,std::vector<CSCDBGains::Item> >::const_iterator layerGainsItr = theGains->gains.begin(),
169  lastGain = theGains->gains.end();
170  for( ; layerGainsItr != lastGain; ++layerGainsItr)
171  {
172  std::cout << "GAIN " << layerGainsItr->first
173  << " STRIPS " << layerGainsItr->second.size() << " "
174  << layerGainsItr->second[0].gain_slope
175  << " " << layerGainsItr->second[0].gain_intercept << std::endl;
176  }
177 
178  std::map< int,std::vector<CSCDBPedestals::Item> >::const_iterator pedestalItr = thePedestals->pedestals.begin(),
179  lastPedestal = thePedestals->pedestals.end();
180  for( ; pedestalItr != lastPedestal; ++pedestalItr)
181  {
182  std::cout << "PEDS " << pedestalItr->first << " "
183  << " STRIPS " << pedestalItr->second.size() << " ";
184  for(int i = 1; i < 80; ++i)
185  {
186  std::cout << pedestalItr->second[i-1].rms << " " ;
187  }
188  std::cout << std::endl;
189  }
190 
191  std::map< int,std::vector<CSCDBCrosstalk::Item> >::const_iterator crosstalkItr = theCrosstalk->crosstalk.begin(),
192  lastCrosstalk = theCrosstalk->crosstalk.end();
193  for( ; crosstalkItr != lastCrosstalk; ++crosstalkItr)
194  {
195  std::cout << "XTALKS " << crosstalkItr->first
196  << " STRIPS " << crosstalkItr->second.size() << " "
197  << crosstalkItr->second[5].xtalk_slope_left << " "
198  << crosstalkItr->second[5].xtalk_slope_right << " "
199  << crosstalkItr->second[5].xtalk_intercept_left << " "
200  << crosstalkItr->second[5].xtalk_intercept_right << std::endl;
201  }
202 */
203 }
204 
205 float CSCConditions::gain(const CSCDetId & detId, int channel) const
206 {
207  assert(theGains.isValid());
208  return float( theGains->item(detId, channel).gain_slope )/theGains->factor_gain;
209 }
210 
211 
212 float CSCConditions::pedestal(const CSCDetId & detId, int channel) const
213 {
214  assert(thePedestals.isValid());
215  return float ( thePedestals->item(detId, channel).ped )/thePedestals->factor_ped;
216 }
217 
218 
219 float CSCConditions::pedestalSigma(const CSCDetId&detId, int channel) const
220 {
221  assert(thePedestals.isValid());
222  return float ( thePedestals->item(detId, channel).rms )/thePedestals->factor_rms;
223 }
224 
225 
226 float CSCConditions::crosstalkIntercept(const CSCDetId&detId, int channel, bool leftRight) const
227 {
228  assert(theCrosstalk.isValid());
229  const CSCDBCrosstalk::Item & item = theCrosstalk->item(detId, channel);
230 
231  // resistive fraction is at the peak, where t=0
232  return leftRight ? float ( item.xtalk_intercept_right )/theCrosstalk->factor_intercept
233  : float ( item.xtalk_intercept_left )/theCrosstalk->factor_intercept ;
234 }
235 
236 
237 float CSCConditions::crosstalkSlope(const CSCDetId&detId, int channel, bool leftRight) const
238 {
239  assert(theCrosstalk.isValid());
240  const CSCDBCrosstalk::Item & item = theCrosstalk->item(detId, channel);
241 
242  // resistive fraction is at the peak, where t=0
243  return leftRight ? float ( item.xtalk_slope_right )/theCrosstalk->factor_slope
244  : float ( item.xtalk_slope_left )/theCrosstalk->factor_slope ;
245 }
246 
247 const CSCDBNoiseMatrix::Item & CSCConditions::noiseMatrix(const CSCDetId& detId, int channel) const
248 {
249  assert(theNoiseMatrix.isValid());
250  return theNoiseMatrix->item(detId, channel);
251 }
252 
253 void CSCConditions::noiseMatrixElements( const CSCDetId& id, int channel, std::vector<float>& me ) const {
254  assert(me.size()>11);
255  const CSCDBNoiseMatrix::Item& item = noiseMatrix(id, channel);
256  me[0] = float ( item.elem33 )/theNoiseMatrix->factor_noise;
257  me[1] = float ( item.elem34 )/theNoiseMatrix->factor_noise;
258  me[2] = float ( item.elem35 )/theNoiseMatrix->factor_noise;
259  me[3] = float ( item.elem44 )/theNoiseMatrix->factor_noise;
260  me[4] = float ( item.elem45 )/theNoiseMatrix->factor_noise;
261  me[5] = float ( item.elem46 )/theNoiseMatrix->factor_noise;
262  me[6] = float ( item.elem55 )/theNoiseMatrix->factor_noise;
263  me[7] = float ( item.elem56 )/theNoiseMatrix->factor_noise;
264  me[8] = float ( item.elem57 )/theNoiseMatrix->factor_noise;
265  me[9] = float ( item.elem66 )/theNoiseMatrix->factor_noise;
266  me[10] = float ( item.elem67 )/theNoiseMatrix->factor_noise;
267  me[11] = float ( item.elem77 )/theNoiseMatrix->factor_noise;
268 }
269 
270 void CSCConditions::crossTalk( const CSCDetId& id, int channel, std::vector<float>& ct ) const {
271  assert(theCrosstalk.isValid());
272  const CSCDBCrosstalk::Item & item = theCrosstalk->item(id, channel);
273  ct[0] = float ( item.xtalk_slope_left )/theCrosstalk->factor_slope;
274  ct[1] = float ( item.xtalk_intercept_left )/theCrosstalk->factor_intercept;
275  ct[2] = float ( item.xtalk_slope_right )/theCrosstalk->factor_slope;
276  ct[3] = float ( item.xtalk_intercept_right )/theCrosstalk->factor_intercept;
277 }
278 
279 float CSCConditions::chipCorrection(const CSCDetId & detId, int stripChannel) const
280 {
281  if ( useTimingCorrections() ){
282  assert(theChipCorrections.isValid());
283  CSCIndexer indexer;
284  int chip = indexer.chipIndex(stripChannel); //Converts 1-80(64) in a layer to 1-5(4), expects ME1/1a to be channel 65-80
285  //printf("CSCCondition e:%d s:%d r:%d c:%d l:%d strip:%d chip: %d\n",detId.endcap(),detId.station(), detId.ring(),detId.chamber(),detId.layer(),stripChannel,chip);
286  return float ( theChipCorrections->item(detId,chip).speedCorr )/theChipCorrections->factor_speedCorr;
287  }
288  else
289  return 0;
290 }
292 {
293  if ( useTimingCorrections() ){
295  return float ( theChamberTimingCorrections->item(detId).cfeb_tmb_skew_delay*1./theChamberTimingCorrections->factor_precision
296  + theChamberTimingCorrections->item(detId).cfeb_timing_corr*1./theChamberTimingCorrections->factor_precision
297  + (theChamberTimingCorrections->item(detId).cfeb_cable_delay*25.)
298 );
299  }
300  else
301  return 0;
302 }
303 float CSCConditions::anodeBXoffset(const CSCDetId & detId) const
304 {
305  if ( useTimingCorrections() ){
307  return float ( theChamberTimingCorrections->item(detId).anode_bx_offset*1./theChamberTimingCorrections->factor_precision);
308  }
309  else
310  return 0;
311 }
312 
313 const std::bitset<80>& CSCConditions::badStripWord( const CSCDetId& id ) const {
314  CSCIndexer indexer;
315  return badStripWords[indexer.layerIndex(id) - 1];
316 }
317 
318 const std::bitset<112>& CSCConditions::badWireWord( const CSCDetId& id ) const {
319  CSCIndexer indexer;
320  return badWireWords[indexer.layerIndex(id) - 1];
321 }
322 
328 
329  const float loEdge = 5.0; // consider gains above this
330  const float hiEdge = 10.0; // consider gains below this
331  const float loLimit = 6.0; // lowest acceptable average gain
332  const float hiLimit = 9.0; // highest acceptable average gain
333  const float expectedAverage = 7.5; // default average gain
334 
335  if ( theAverageGain > 0. ) return theAverageGain; // only recalculate if necessary
336 
337  int n_strip = 0;
338  float gain_tot = 0.;
339 
340  CSCDBGains::GainContainer::const_iterator it;
341  for ( it=theGains->gains.begin(); it!=theGains->gains.end(); ++it ) {
342  float the_gain = float( it->gain_slope )/theGains->factor_gain;
343  if (the_gain > loEdge && the_gain < hiEdge ) {
344  gain_tot += the_gain;
345  ++n_strip;
346  }
347  }
348 
349  // Average gain
350  if ( n_strip > 0 ) {
351  theAverageGain = gain_tot / n_strip;
352  }
353 
354  // Average gain has been around 7.5 in real data
355  if ( theAverageGain < loLimit || theAverageGain > hiLimit ) {
356  // LogTrace("CSC") << "Average CSC strip gain = "
357  // << theAverageGain << " is reset to expected value " << expectedAverage;
358  theAverageGain = expectedAverage;
359  }
360 
361  return theAverageGain;
362 }
T getParameter(std::string const &) const
const std::bitset< 80 > & badStripWord(const CSCDetId &id) const
return bad channel words per CSCLayer - 1 bit per channel
int i
Definition: DBlmapReader.cc:9
edm::ESHandle< CSCDBNoiseMatrix > theNoiseMatrix
edm::ESHandle< CSCBadWires > theBadWires
float pedestalSigma(const CSCDetId &detId, int channel) const
void noiseMatrixElements(const CSCDetId &id, int channel, std::vector< float > &me) const
fill vector (dim 12, must be allocated by caller) with noise matrix elements (scaled to float) ...
void fillBadStripWords()
fill bad channel words
edm::ESHandle< CSCChamberTimeCorrections > theChamberTimingCorrections
bool useTimingCorrections() const
did we request reading timing correction info from db?
Definition: CSCConditions.h:88
const CSCDBNoiseMatrix::Item & noiseMatrix(const CSCDetId &detId, int channel) const
return raw noise matrix (unscaled short int elements)
void print() const
edm::ESHandle< CSCBadStrips > theBadStrips
float chipCorrection(const CSCDetId &detId, int channel) const
std::vector< std::bitset< 112 > > badWireWords
float chamberTimingCorrection(const CSCDetId &detId) const
void crossTalk(const CSCDetId &id, int channel, std::vector< float > &ct) const
fill vector (dim 4, must be allocated by caller) with crosstalk sl, il, sr, ir
edm::ESWatcher< CSCDBGainsRcd > gainsWatcher_
std::vector< std::bitset< 80 > > badStripWords
float anodeBXoffset(const CSCDetId &detId) const
bool useTimingCorrections_
float pedestal(const CSCDetId &detId, int channel) const
in ADC counts
bool readBadChannels() const
did we request reading bad channel info from db?
Definition: CSCConditions.h:82
CSCDetId detIdFromChamberIndex(IndexType ici) const
Definition: CSCIndexer.cc:62
edm::ESHandle< CSCBadChambers > theBadChambers
int j
Definition: DBlmapReader.cc:9
edm::ESHandle< CSCDBGains > theGains
Definition: CSCConditions.h:99
void initializeEvent(const edm::EventSetup &es)
fetch the maps from the database
const std::bitset< 112 > & badWireWord(const CSCDetId &id) const
float crosstalkSlope(const CSCDetId &detId, int channel, bool leftRight) const
float crosstalkIntercept(const CSCDetId &detId, int channel, bool leftRight) const
bool isInBadChamber(const CSCDetId &id) const
Is the gven chamber flagged as bad?
float gain(const CSCDetId &detId, int channel) const
channels count from 1
CSCConditions(const edm::ParameterSet &ps)
void fillBadWireWords()
edm::ESHandle< CSCDBChipSpeedCorrection > theChipCorrections
edm::ESHandle< CSCDBCrosstalk > theCrosstalk
bool readBadChambers() const
did we request reading bad chamber info from db?
Definition: CSCConditions.h:85
const T & get() const
Definition: EventSetup.h:55
string const
Definition: compareJSON.py:14
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:59
IndexType layerIndex(const CSCDetId &id) const
Definition: CSCIndexer.h:79
IndexType chipIndex(IndexType ie, IndexType is, IndexType ir, IndexType ic, IndexType il, IndexType ichip) const
Definition: CSCIndexer.h:236
float averageGain() const
average gain over entire CSC system (logically const although must be cached here).
float theAverageGain
bool isValid() const
Definition: ESHandle.h:37
edm::ESHandle< CSCDBPedestals > thePedestals