CMS 3D CMS Logo

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