CMS 3D CMS Logo

EcalDccWeightBuilder.cc
Go to the documentation of this file.
1 /*
2  *
3  * authors: Ph. Gras (CEA/Saclay), F. Cavallari (INFN/Roma)
4  * some code copied from CalibCalorimetry/EcalTPGTools code
5  * written by P. Paganini and F. Cavallari
6  */
7 
8 #define DB_WRITE_SUPPORT
9 
11 
12 #include <limits>
13 #include <algorithm>
14 #include <fstream>
15 #include <iomanip>
16 #include <TFile.h>
17 #include <TTree.h>
18 
25 
29 
30 #ifdef DB_WRITE_SUPPORT
33 #endif //DB_WRITE_SUPPORT defined
34 
36 
37 using namespace std;
38 using namespace edm;
39 
40 const double EcalDccWeightBuilder::weightScale_ = 1024.;
41 
42 
43 //TODO: handling case of weight encoding saturation: weights shall be downscaled to prevent saturation
44 
46  dcc1stSample_(ps.getParameter<int>("dcc1stSample")),
47  sampleToSkip_(ps.getParameter<int>("sampleToSkip")),
48  nDccWeights_(ps.getParameter<int>("nDccWeights")),
49  inputWeights_(ps.getParameter<vector<double> >("inputWeights")),
50  mode_(ps.getParameter<string>("mode")),
51  dccWeightsWithIntercalib_(ps.getParameter<bool>("dccWeightsWithIntercalib")),
52  writeToDB_(ps.getParameter<bool>("writeToDB")),
53  writeToAsciiFile_(ps.getParameter<bool>("writeToAsciiFile")),
54  writeToRootFile_(ps.getParameter<bool>("writeToRootFile")),
55  asciiOutputFileName_(ps.getParameter<string>("asciiOutputFileName")),
56  rootOutputFileName_(ps.getParameter<string>("rootOutputFileName")),
57  dbSid_(ps.getParameter<string>("dbSid")),
58  dbUser_(ps.getParameter<string>("dbUser")),
59  dbPassword_(ps.getUntrackedParameter<string>("dbPassword","")),
60  dbTag_(ps.getParameter<string>("dbTag")),
61  dbVersion_(ps.getParameter<int>("dbVersion")),
62  sqlMode_(ps.getParameter<bool>("sqlMode")),
63  calibMap_(emptyCalibMap_)
64 {
65  if(mode_=="weightsFromConfig"){
67  if(inputWeights_.size()!=(unsigned)nDccWeights_){
68  throw cms::Exception("Config")
69  << "Inconsistent configuration. 'nDccWeights' parameters indicated "
70  << nDccWeights_ << " weights while parameter 'inputWeights_' contains "
71  << inputWeights_.size() << " weight values!\n";
72  }
73  } else if(mode_=="computeWeights"){
75  } else{
76  throw cms::Exception("Config")
77  << "Invalid value ('" << mode_ << "') for parameter mode. "
78  << "Valid values are: 'weightsFromConfig' and 'computeWeights'\n";
79  }
80 }
81 
82 
83 void
85  const edm::EventSetup& es){
86 
88  es.get<EcalMappingRcd>().get(handle);
89  ecalElectronicsMap_ = handle.product();
90 
91  // Retrieval of intercalib constants
94  es.get<EcalIntercalibConstantsRcd>().get(hIntercalib) ;
95  const EcalIntercalibConstants* intercalib = hIntercalib.product();
96  calibMap_ = intercalib->getMap();
97  }
98 
99  //gets geometry
100  es.get<CaloGeometryRecord>().get(geom_);
101 
102 
103  //computes the weights:
105 
106  //Writing out weights.
110 }
111 
112 void EcalDccWeightBuilder::computeAllWeights(bool withIntercalib){
113  const int nw = nDccWeights_;
114  int iSkip0_ = sampleToSkip_>=0?(sampleToSkip_-dcc1stSample_):-1;
115 
116  EcalSimParameterMap parameterMap;
117  const vector<DetId>& ebDetIds
119 
120  // cout << __FILE__ << ":" << __LINE__ << ": "
121  // << "Number of EB det IDs: " << ebDetIds.size() << "\n";
122 
123  const vector<DetId>& eeDetIds
125 
126  // cout << __FILE__ << ":" << __LINE__ << ": "
127  // << "Number of EE det IDs: " << eeDetIds.size() << "\n";
128 
129 
130  vector<DetId> detIds(ebDetIds.size()+eeDetIds.size());
131  copy(ebDetIds.begin(), ebDetIds.end(), detIds.begin());
132  copy(eeDetIds.begin(), eeDetIds.end(), detIds.begin()+ebDetIds.size());
133 
134  vector<double> baseWeights(nw); //weight obtained from signal shape
135  vector<double> w(nw); //weight*intercalib
136  vector<int> W(nw); //weight in hw encoding (integrer)
137  double prevPhase = numeric_limits<double>::min();
138 
139 
141  assert(inputWeights_.size()==baseWeights.size());
142  copy(inputWeights_.begin(), inputWeights_.end(), baseWeights.begin());
143  }
144 
145  for(vector<DetId>::const_iterator it = detIds.begin();
146  it != detIds.end(); ++it){
147 
148  double phase = parameterMap.simParameters(*it).timePhase();
149  int binOfMax = parameterMap.simParameters(*it).binOfMaximum();
150 
151 #if 0
152  //for debugging...
153  cout << __FILE__ << ":" << __LINE__ << ": ";
154  if(it->subdetId()==EcalBarrel){
155  cout << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
156  << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
157  } else if(it->subdetId()==EcalEndcap){
158  cout << "ix = " << setw(3) << ((EEDetId)(*it)).ix()
159  << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
160  << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
161  } else{
162  throw cms::Exception("EcalDccWeightBuilder")
163  << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
164  << "Got a detId which is neither tagged as ECAL Barrel "
165  << "not ECAL endcap while looping on ECAL cell detIds\n";
166  }
167  cout << " -> phase: " << phase << "\n";
168  cout << " -> binOfMax: " << binOfMax << "\n";
169 #endif
170 
171  try{
172  EBShape ebShape;
173  EEShape eeShape;
174  EcalShapeBase* pShape;
175 
176  if(it->subdetId()==EcalBarrel){
177  pShape = &ebShape;
178  } else if(it->subdetId()==EcalEndcap){
179  pShape = &eeShape;
180  } else{
181  throw cms::Exception("EcalDccWeightBuilder")
182  << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
183  << "Got a detId which is neither tagged as ECAL Barrel "
184  << "not ECAL endcap while looping on ECAL cell detIds\n";
185  }
186 
187  if(phase!=prevPhase){
188  if(imode_==COMPUTE_WEIGHTS){
189  if(it->subdetId()==EcalBarrel){
190  computeWeights(*pShape, binOfMax, phase,
191  dcc1stSample_-1, nDccWeights_, iSkip0_,
192  baseWeights);
193  }
194  prevPhase = phase;
195  }
196  }
197  for(int i = 0; i < nw; ++i){
198  w[i] = baseWeights[i];
199  if(withIntercalib) w[i]*= intercalib(*it);
200  }
201  unbiasWeights(w, &W);
202  encodedWeights_[*it] = W;
203  } catch(std::exception& e){
204  cout << __FILE__ << ":" << __LINE__ << ": ";
205  if(it->subdetId()==EcalBarrel){
206  cout << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
207  << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
208  } else if(it->subdetId()==EcalEndcap){
209  cout << "ix = " << setw(3) << ((EEDetId)(*it)).ix()
210  << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
211  << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
212  } else{
213  cout << "DetId " << (uint32_t) (*it);
214  }
215  cout << "phase: " << phase << "\n";
216  throw;
217  }
218  }
219 }
220 
221 void
223  int binOfMax,
224  double timePhase,
225  int iFirst,
226  int nWeights, int iSkip,
227  vector<double>& result){
228  double sum2 = 0.;
229  double sum = 0;
230  result.resize(nWeights);
231 
232  int nActualWeights = 0;
233 
234  const double tzero = -(binOfMax-1)*25+timePhase + shape.timeToRise();//ns
235 
236  for(int i=0; i<nWeights; ++i){
237  double t_ns = tzero+(iFirst+i)*25;
238  double s = shape(t_ns);
239  if(i==iSkip){
240  continue;
241  }
242  result[i] = s;
243  sum += s;
244  sum2 += s*s;
245  ++nActualWeights;
246  }
247  for(int i=0; i<nWeights; ++i){
248  if(i==iSkip){
249  result[i] = 0;
250  } else{
251  result[i] = (result[i]-sum/nActualWeights)/(sum2-sum*sum/nActualWeights);
252  }
253  }
254 }
255 
257  return lround(w * weightScale_);
258 }
259 
261  return ((double) W) / weightScale_;
262 }
263 
264 
265 template<class T>
266 void EcalDccWeightBuilder::sort(const std::vector<T>& a,
267  std::vector<int>& s,
268  bool decreasingOrder){
269  // cout << __FILE__ << ":" << __LINE__ << ": "
270  // << "sort input array:" ;
271  // for(unsigned i=0; i<a.size(); ++i){
272  // cout << "\t" << a[i];
273  // }
274  // cout << "\n";
275 
276  //performs a bubble sort: adjacent elements are successively swapped 2 by 2
277  //until the list is finally sorted.
278  bool changed = false;
279  s.resize(a.size());
280  for(unsigned i=0; i<a.size(); ++i) s[i] = i;
281  if(a.size() == 0) return;
282  do {
283  changed = false;
284  for(unsigned i = 0; i < a.size()-1; ++i){
285  const int j = s[i];
286  const int nextj = s[i+1];
287  if((decreasingOrder && (a[j] < a[nextj]))
288  || (!decreasingOrder && (a[j] > a[nextj]))){
289  std::swap(s[i], s[i+1]);
290  changed = true;
291  }
292  }
293  } while(changed);
294 
295  // cout << __FILE__ << ":" << __LINE__ << ": "
296  // << "sorted list of indices:" ;
297  // for(unsigned i=0; i < s.size(); ++i){
298  // cout << "\t" << s[i];
299  // }
300  // cout << "\n";
301 }
302 
304  std::vector<int>* encodedWeights){
305  const unsigned nw = weights.size();
306 
307  //computes integer weights, weights residuals and weight sum residual:
308  vector<double> dw(nw); //weight residuals due to interger encoding
309  vector<int> W(nw); //integer weights
310  int wsum = 0;
311  for(unsigned i = 0; i < nw; ++i){
312  W[i] = encodeWeight(weights[i]);
313  dw[i] = decodeWeight(W[i]) - weights[i];
314  wsum += W[i];
315  }
316 
317 // cout << __FILE__ << ":" << __LINE__ << ": "
318 // << "weights before bias correction: ";
319 // for(unsigned i=0; i<weights.size(); ++i){
320 // const double w = weights[i];
321 // cout << "\t" << encodeWeight(w) << "(" << w << ", dw = " << dw[i] << ")";
322 // }
323 // cout << "\t sum: " << wsum << "\n";
324 
325  //sorts weight residuals in decreasing order:
326  vector<int> iw(nw);
327  sort(dw, iw, true);
328 
329  //compensates weight sum residual by adding or substracting 1 to weights
330  //starting from:
331  // 1) the weight with the minimal signed residual if the correction
332  // is positive (wsum<0)
333  // 2) the weight with the maximal signed residual if the correction
334  // is negative (wsum>0)
335  int wsumSign = wsum>0?1:-1;
336  int i = wsum>0?0:(nw-1);
337  while(wsum!=0){
338  W[iw[i]] -= wsumSign;
339  wsum -= wsumSign;
340  i += wsumSign;
341  if(i<0 || i>=(int)nw){ //recompute the residuals if a second iteration is
342  // needed (in principle, it is not expected with usual input weights), :
343  for(unsigned i = 0; i < nw; ++i){
344  dw[i] = decodeWeight(W[i]) - weights[i];
345  sort(dw, iw, true);
346  }
347  }
348  if(i<0) i = nw-1;
349  if(i>=(int)nw) i = 0;
350  }
351 
352 // cout << __FILE__ << ":" << __LINE__ << ": "
353 // << "weights after bias correction: ";
354 // for(unsigned i=0; i<weights.size(); ++i){
355 // cout << "\t" << W[i] << "(" << decodeWeight(W[i]) << ", dw = "
356 // << (decodeWeight(W[i])-weights[i]) << ")";
357 // }
358 // cout << "\n";
359 
360  //copy result
361  if(encodedWeights!=0) encodedWeights->resize(nw);
362  for(unsigned i = 0; i < nw; ++i){
363  weights[i] = decodeWeight(W[i]);
364  if(encodedWeights) (*encodedWeights)[i] = W[i];
365  }
366 }
367 
369  // get current intercalibration coeff
370  double coef;
372  = calibMap_.find(detId.rawId());
373  if(itCalib != calibMap_.end()){
374  coef = (*itCalib);
375  } else{
376  coef = 1.;
377  std::cout << (uint32_t) detId
378  << " not found in EcalIntercalibConstantMap"<<std::endl ;
379  }
380 #if 0
381  cout << __FILE__ << ":" << __LINE__ << ": ";
382  if(detId.subdetId()==EcalBarrel){
383  cout << "ieta = " << ((EBDetId)detId).ieta()
384  << " iphi = " << ((EBDetId)detId).iphi();
385  } else{
386  cout << "ix = " << ((EEDetId)detId).ix()
387  << " iy = " << ((EEDetId)detId).iy()
388  << " iz = " << ((EEDetId)detId).zside();
389  }
390  cout << " coef = " << coef << "\n";
391 #endif
392  return coef;
393 }
394 
396  string fName = asciiOutputFileName_.size()!=0?
397  asciiOutputFileName_.c_str()
398  :"dccWeights.txt";
399  ofstream file(fName.c_str());
400  if(!file.good()){
401  throw cms::Exception("Output")
402  << "Failed to open file '"
403  << fName
404  << "'for writing DCC weights\n";
405  }
406 
407  const char* comment = sqlMode_?"-- ":"# ";
408 
409  file << comment << "List of weights for amplitude estimation to be used in DCC for\n"
410  << comment << "zero suppresssion.\n\n";
411  if(!sqlMode_){
412  file << comment << "Note: RU: trigger tower in EB, supercrystal in EE\n"
413  << comment << " xtl: crystal electronic channel id in RU, from 1 to 25\n\n"
414  << comment << " DetId SM FED RU xtl weights[0..5]...\n";
415  }
416 
417  if(sqlMode_){
418  file << "variable recid number;\n"
419  "exec select COND2CONF_INFO_SQ.NextVal into :recid from DUAL;\n"
420  "insert into weights_info (rec_id,tag,version) values (:recid,'"
421  << dbTag_ << "'," << dbVersion_ << ");\n";
422  file << "\n" << comment
423  << "index of first sample used in the weighting sum\n"
424  "begin\n"
425  " for fedid in " << ecalDccFedIdMin << ".." << ecalDccFedIdMax
426  << " loop\n"
427  " insert into dcc_weightsample_dat (rec_id, logic_id, sample_id, \n"
428  " weight_number)\n"
429  " values(:recid,fedid," << dcc1stSample_ << ",1);\n"
430  " end loop;\n"
431  "end;\n"
432  "/\n";
433  } else{
434  file << "1st DCC sample: " << dcc1stSample_ << "\n";
435  }
436 
437  file << "\n" << comment << "list of weights per crystal channel\n";
438 
439  for(map<DetId, std::vector<int32_t> >::const_iterator it
440  = encodedWeights_.begin();
441  it != encodedWeights_.end();
442  ++it){
443  const DetId& detId = it->first;
444 
445  int fedId;
446  int smId;
447  int ruId;
448  int xtalId;
449 
450  //detId -> fedId, smId, ruId, xtalId
451  dbId(detId, fedId, smId, ruId, xtalId);
452 
453  char delim = sqlMode_?',':' ';
454 
455  if(sqlMode_) file << "-- detId " << detId.rawId() << "\n"
456  << "insert into dcc_weights_dat(rec_id,sm_id,fed_id,"
457  "tt_id, cry_id,\n"
458  "weight_0,weight_1,weight_2,weight_3,weight_4,weight_5) \n"
459  "values ("
460  ":recid";
461 
462  const vector<int>& weights = it->second;
463  if(!sqlMode_) file << setw(10) << detId.rawId();
464  file << delim << setw(2) << smId;
465  file << delim << setw(3) << fedId;
466  file << delim << setw(2) << ruId;
467  file << delim << setw(2) << xtalId;
468 
469  for(unsigned i=0; i<weights.size(); ++i){
470  file << delim << setw(5) << weights[i];
471  }
472  if(sqlMode_) file << ");";
473  file << "\n";
474  }
475  if(!file.good()){
476  throw cms::Exception("Output") << "Error while writing DCC weights to '"
477  << fName << "' file.";
478  }
479 }
481  string fName = rootOutputFileName_.size()!=0?
482  rootOutputFileName_.c_str()
483  :"dccWeights.root";
484  TFile file(fName.c_str(), "RECREATE");
485  if(file.IsZombie()){
486  throw cms::Exception("Output")
487  << "Failed to open file '"
488  << fName
489  << "'for writing DCC weights\n";
490  }
491  TTree t("dccWeights", "Weights for DCC ZS filter");
492  const int nWeightMax = 20; //normally n_weights = 6. A different might be used
493  // used for test purposes.
494  struct {
495  Int_t detId;
496  Int_t fedId;
497  Int_t smId;
498  Int_t ruId;
499  Int_t xtalId;
500  Int_t n_weights;
501  Int_t weights[nWeightMax];
502  } buf;
503  t.Branch("weights", &buf,
504  "rawDetId/I:"
505  "feId/I:"
506  "smSlotId/I:"
507  "ruId/I:"
508  "xtalInRuId/I:"
509  "n_weights/I:"
510  "weights[n_weights]/I");
511  for(map<DetId, std::vector<int32_t> >::const_iterator it
512  = encodedWeights_.begin();
513  it != encodedWeights_.end();
514  ++it){
515  buf.detId = it->first.rawId();
516  buf.n_weights = it->second.size();
517 
518  //detId -> fedId, smId, ruId, xtalId
519  dbId(buf.detId, buf.fedId, buf.smId, buf.ruId, buf.xtalId);
520 
521  if(buf.n_weights>nWeightMax){
522  throw cms::Exception("EcalDccWeight")
523  << "Number of weights (" << buf.n_weights
524  << ") for DetId " << buf.detId
525  << " exceeded maximum limit (" << nWeightMax
526  << ") of root output format. ";
527  }
528  copy(it->second.begin(), it->second.end(), buf.weights);
529  t.Fill();
530  }
531  t.Write();
532  file.Close();
533 }
534 
535 #ifndef DB_WRITE_SUPPORT
537  throw cms::Exception("DccWeight")
538  << "Code was compiled without support for writing dcc weights directly "
539  " into configuration DB. Configurable writeToDB must be set to False. "
540  "sqlMode can be used to produce an SQL*PLUS script to fill the DB\n";
541 }
542 #else //DB_WRITE_SUPPORT defined
544  cout << "going to write to the online DB "<<dbSid_<<" user "<<dbUser_<<endl;;
545  EcalCondDBInterface* econn;
546 
547  try {
548  cout << "Making connection..." << flush;
549  const string& filePrefix = string("file:");
550  if(dbPassword_.find(filePrefix)==0){ //password must be read for a file
551  string fileName = dbPassword_.substr(filePrefix.size());
552  //substitute dbPassword_ value by the password read from the file
553  PasswordReader pr;
554  pr.readPassword(fileName, dbUser_, dbPassword_);
555  }
556 
557  // cout << __FILE__ << ":" << __LINE__ << ": "
558  // << "Password: " << dbPassword_ << "\n";
559 
561  cout << "Done." << endl;
562  } catch (runtime_error &e) {
563  cerr << e.what() << endl;
564  exit(-1);
565  }
566 
567  ODFEWeightsInfo weight_info;
568  weight_info.setConfigTag(dbTag_);
569  weight_info.setVersion(dbVersion_);
570  cout << "Inserting in DB..." << endl;
571 
572  econn->insertConfigSet(&weight_info);
573 
574  int weight_id=weight_info.getId();
575  cout << "WeightInfo inserted with ID "<< weight_id<< endl;
576 
577  vector<ODWeightsDat> datadel;
578  datadel.reserve(encodedWeights_.size());
579 
580  vector<ODWeightsSamplesDat> dcc1stSampleConfig(nDccs);
581  for(int i = ecalDccFedIdMin; i <= ecalDccFedIdMax; ++i){
582  dcc1stSampleConfig[i].setId(weight_id);
583  dcc1stSampleConfig[i].setFedId(601+i);
584  dcc1stSampleConfig[i].setSampleId(dcc1stSample_);
585  dcc1stSampleConfig[i].setWeightNumber(-1); //not used.
586  }
587  econn->insertConfigDataArraySet(dcc1stSampleConfig, &weight_info);
588 
589  for(map<DetId, std::vector<int32_t> >::const_iterator it
590  = encodedWeights_.begin();
591  it != encodedWeights_.end();
592  ++it){
593  const DetId& detId = it->first;
594  const unsigned nWeights = 6;
595  vector<int> weights(nWeights);
596 
597  for(unsigned i=0; i<weights.size(); ++i){
598  //completing the weight vector with zeros in case it has
599  //less than 6 elements:
600  const vector<int>& w = it->second;
601  weights[i] = i<w.size()?w[i]:0;
602  }
603 
604  ODWeightsDat one_dat;
605  one_dat.setId(weight_id);
606 
607  int fedId;
608  int smId;
609  int ruId;
610  int xtalId;
611 
612  //detId -> fedId, smId, ruId, xtalId
613  dbId(detId, fedId, smId, ruId, xtalId);
614 
615  one_dat.setSMId(smId);
616  one_dat.setFedId(fedId);
617  one_dat.setTTId(ruId);
618  one_dat.setCrystalId(xtalId);
619 
620  one_dat.setWeight0(weights[0]);
621  one_dat.setWeight1(weights[1]);
622  one_dat.setWeight2(weights[2]);
623  one_dat.setWeight3(weights[3]);
624  one_dat.setWeight4(weights[4]);
625  one_dat.setWeight5(weights[5]);
626 
627  datadel.push_back(one_dat);
628  }
629  econn->insertConfigDataArraySet(datadel,&weight_info);
630  std::cout<< " .. done insertion in DB "<< endl;
631  delete econn;
632  cout<< "closed DB connection ... done" << endl;
633 }
634 #endif //DB_WRITE_SUPPORT not defined
635 
636 
637 void EcalDccWeightBuilder::dbId(const DetId& detId, int& fedId, int& smId,
638  int& ruId,
639  int& xtalId) const{
640  const EcalElectronicsId& elecId
642 
643  fedId = 600 + elecId.dccId();
645 
646  if(detId.subdetId()==EcalBarrel) {
647  smId=((EBDetId)detId).ism();
648  } else{
649  smId = 10000-fedId; //no SM in EE. Use some unique value to satisfy
650  // current DB PK constraints.
651  }
652  const int stripLength = 5; //1 strip = 5 crystals in a row
653  xtalId = (elecId.stripId()-1)*stripLength + elecId.xtalId();
654 
655 #if 0
656  cout << __FILE__ << ":" << __LINE__ << ": FED ID "
657  << fedId << "\n";
658 
659  cout << __FILE__ << ":" << __LINE__ << ": SM logical ID "
660  << smId << "\n";
661 
662  cout << __FILE__ << ":" << __LINE__ << ": RU ID (TT or SC): "
663  << ruId << "\n";
664 
665  cout << __FILE__ << ":" << __LINE__ << ": strip:"
666  << elecId.stripId() << "\n";
667 
668  cout << __FILE__ << ":" << __LINE__ << ": xtal in strip: "
669  << elecId.xtalId() << "\n";
670 
671  cout << __FILE__ << ":" << __LINE__ << ": xtalId in RU: "
672  << xtalId << "\n";
673 #endif
674 }
void setId(int dac)
Definition: ODWeightsDat.h:21
void setWeight3(float x)
Definition: ODWeightsDat.h:39
void setWeight2(float x)
Definition: ODWeightsDat.h:38
int xtalId() const
get the channel id
const double w
Definition: UKUtility.cc:23
std::map< DetId, std::vector< int > > encodedWeights_
virtual const CaloSimParameters & simParameters(const DetId &id) const
return the sim parameters relative to the right subdet
int stripId() const
get the tower id
void setSMId(int dac)
Definition: ODWeightsDat.h:24
const self & getMap() const
double intercalib(const DetId &detId)
Ecal readout channel identification [32:20] Unused (so far) [19:13] DCC id [12:6] tower [5:3] strip [...
std::vector< double > inputWeights_
edm::ESHandle< CaloGeometry > geom_
static const double weightScale_
void dbId(const DetId &detId, int &fedId, int &smId, int &ruId, int &xtalId) const
int towerId() const
get the tower id
void setWeight4(float x)
Definition: ODWeightsDat.h:40
void computeWeights(const EcalShapeBase &shape, int binOfMax, double timePhase, int iFirst0, int nWeights, int iSkip0, std::vector< double > &result)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
double timePhase() const
the adjustment you need to apply to get the signal where you want it
void setVersion(int id)
void sort(const std::vector< T > &a, std::vector< int > &s, bool decreasingOrder=false)
EcalElectronicsId getElectronicsId(const DetId &id) const
Get the electronics id for this det id.
void setWeight1(float x)
Definition: ODWeightsDat.h:37
void insertConfigDataArraySet(const std::vector< DATT > &data, ICONF *iconf) noexcept(false)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
EcalIntercalibConstantMap & calibMap_
void setTTId(int dac)
Definition: ODWeightsDat.h:30
void insertConfigSet(ICONF *iconf) noexcept(false)
void setWeight0(float x)
Definition: ODWeightsDat.h:36
void setCrystalId(int dac)
Definition: ODWeightsDat.h:33
T min(T a, T b)
Definition: MathUtil.h:58
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
int dccId() const
get the DCC (Ecal Local DCC value not global one) id
Definition: EBShape.h:6
Definition: DetId.h:18
void readPassword(const std::string &fileName, const std::string &user, std::string &password)
const T & get() const
Definition: EventSetup.h:56
std::vector< Item >::const_iterator const_iterator
void unbiasWeights(std::vector< double > &weights, std::vector< int32_t > *encodedWeigths)
EcalDccWeightBuilder(edm::ParameterSet const &ps)
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:92
int getId() const
void setConfigTag(std::string x)
Definition: IODConfig.h:31
HLT enums.
void setWeight5(float x)
Definition: ODWeightsDat.h:41
double a
Definition: hdecay.h:121
void analyze(const edm::Event &event, const edm::EventSetup &es)
const_iterator find(uint32_t rawId) const
static const double tzero[3]
const_iterator end() const
static const int ecalDccFedIdMax
void computeAllWeights(bool withIntercalib)
void setFedId(int dac)
Definition: ODWeightsDat.h:27
int binOfMaximum() const
T const * product() const
Definition: ESHandle.h:86
virtual double timeToRise() const
const EcalElectronicsMapping * ecalElectronicsMap_
static const int ecalDccFedIdMin
Definition: event.py:1
Definition: EEShape.h:6
#define comment(par)
Definition: vmac.h:161