CMS 3D CMS Logo

GenWeightsTableProducer.cc
Go to the documentation of this file.
13 
14 #include <vector>
15 #include <unordered_map>
16 #include <iostream>
17 #include <regex>
18 
19 namespace {
21  struct Counter {
22  Counter() :
23  num(0), sumw(0), sumw2(0), sumPDF(), sumScale(), sumNamed() {}
24 
25  // the counters
26  long long num;
27  long double sumw;
28  long double sumw2;
29  std::vector<long double> sumPDF, sumScale, sumNamed;
30 
31  void clear() {
32  num = 0; sumw = 0; sumw2 = 0;
33  sumPDF.clear(); sumScale.clear(); sumNamed.clear();
34  }
35 
36  // inc the counters
37  void incGenOnly(double w) {
38  num++; sumw += w; sumw2 += (w*w);
39  }
40  void incLHE(double w0, const std::vector<double> & wScale, const std::vector<double> & wPDF, const std::vector<double> & wNamed) {
41  // add up weights
42  incGenOnly(w0);
43  // then add up variations
44  if (!wScale.empty()) {
45  if (sumScale.empty()) sumScale.resize(wScale.size(), 0);
46  for (unsigned int i = 0, n = wScale.size(); i < n; ++i) sumScale[i] += (w0 * wScale[i]);
47  }
48  if (!wPDF.empty()) {
49  if (sumPDF.empty()) sumPDF.resize(wPDF.size(), 0);
50  for (unsigned int i = 0, n = wPDF.size(); i < n; ++i) sumPDF[i] += (w0 * wPDF[i]);
51  }
52  if (!wNamed.empty()) {
53  if (sumNamed.empty()) sumNamed.resize(wNamed.size(), 0);
54  for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) sumNamed[i] += (w0 * wNamed[i]);
55  }
56  }
57 
58  void merge(const Counter & other) {
59  num += other.num; sumw += other.sumw; sumw2 += other.sumw2;
60  if (sumScale.empty() && !other.sumScale.empty()) sumScale.resize(other.sumScale.size(),0);
61  if (sumPDF.empty() && !other.sumPDF.empty()) sumPDF.resize(other.sumPDF.size(),0);
62  if (sumNamed.empty() && !other.sumNamed.empty()) sumNamed.resize(other.sumNamed.size(),0);
63  for (unsigned int i = 0, n = sumScale.size(); i < n; ++i) sumScale[i] += other.sumScale[i];
64  for (unsigned int i = 0, n = sumPDF.size(); i < n; ++i) sumPDF[i] += other.sumPDF[i];
65  for (unsigned int i = 0, n = sumNamed.size(); i < n; ++i) sumNamed[i] += other.sumNamed[i];
66  }
67  };
68 
70  struct DynamicWeightChoice {
71  // choice of LHE weights
72  // ---- scale ----
73  std::vector<std::string> scaleWeightIDs;
74  std::string scaleWeightsDoc;
75  // ---- pdf ----
76  std::vector<std::string> pdfWeightIDs;
77  std::string pdfWeightsDoc;
78  };
79 
80  float stof_fortrancomp(const std::string &str) {
81  std::string::size_type match = str.find("d");
82  if (match != std::string::npos) {
83  std::string pre = str.substr(0,match);
84  std::string post = str.substr(match+1);
85  return std::stof(pre) * std::pow(10.0f, std::stof(post));
86  } else {
87  return std::stof(str);
88  }
89  }
91  struct ScaleVarWeight {
92  std::string wid, label;
93  std::pair<float,float> scales;
94  ScaleVarWeight(const std::string & id, const std::string & text, const std::string & muR, const std::string & muF) :
95  wid(id), label(text), scales(stof_fortrancomp(muR), stof_fortrancomp(muF)) {}
96  bool operator<(const ScaleVarWeight & other) { return (scales == other.scales ? wid < other.wid : scales < other.scales); }
97  };
98  struct PDFSetWeights {
99  std::vector<std::string> wids;
100  std::pair<unsigned int,unsigned int> lhaIDs;
101  PDFSetWeights(const std::string & wid, unsigned int lhaID) : wids(1,wid), lhaIDs(lhaID,lhaID) {}
102  bool operator<(const PDFSetWeights & other) const { return lhaIDs < other.lhaIDs; }
103  void add(const std::string & wid, unsigned int lhaID) {
104  wids.push_back(wid);
105  lhaIDs.second = lhaID;
106  }
107  bool maybe_add(const std::string & wid, unsigned int lhaID) {
108  if (lhaID == lhaIDs.second+1) {
109  lhaIDs.second++;
110  wids.push_back(wid);
111  return true;
112  } else {
113  return false;
114  }
115  }
116  };
117 }
118 
119 class GenWeightsTableProducer : public edm::global::EDProducer<edm::StreamCache<Counter>, edm::RunCache<DynamicWeightChoice>, edm::RunSummaryCache<Counter>, edm::EndRunProducer> {
120  public:
122  genTag_(consumes<GenEventInfoProduct>(params.getParameter<edm::InputTag>("genEvent"))),
123  lheLabel_(params.getParameter<edm::InputTag>("lheInfo")),
124  lheTag_(consumes<LHEEventProduct>(lheLabel_)),
125  lheRunTag_(consumes<LHERunInfoProduct, edm::InRun>(lheLabel_)),
126  namedWeightIDs_(params.getParameter<std::vector<std::string>>("namedWeightIDs")),
127  namedWeightLabels_(params.getParameter<std::vector<std::string>>("namedWeightLabels")),
128  lheWeightPrecision_(params.getParameter<int32_t>("lheWeightPrecision")),
129  maxPdfWeights_(params.getParameter<uint32_t>("maxPdfWeights")),
130  debug_(params.getUntrackedParameter<bool>("debug",false)), debugRun_(debug_.load()),
131  hasIssuedWarning_(false)
132  {
133  produces<nanoaod::FlatTable>();
134  produces<nanoaod::FlatTable>("LHEScale");
135  produces<nanoaod::FlatTable>("LHEPdf");
136  produces<nanoaod::FlatTable>("LHENamed");
137  produces<nanoaod::MergeableCounterTable,edm::Transition::EndRun>();
138  if (namedWeightIDs_.size() != namedWeightLabels_.size()) {
139  throw cms::Exception("Configuration", "Size mismatch between namedWeightIDs & namedWeightLabels");
140  }
141  for (const edm::ParameterSet & pdfps : params.getParameter<std::vector<edm::ParameterSet>>("preferredPDFs")) {
142  const std::string & name = pdfps.getParameter<std::string>("name");
143  uint32_t lhaid = pdfps.getParameter<uint32_t>("lhaid");
144  preferredPDFLHAIDs_.push_back(lhaid);
145  lhaNameToID_[name] = lhaid;
146  lhaNameToID_[name+".LHgrid"] = lhaid;
147  }
148  }
149 
151 
152  void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
153  // get my counter for weights
154  Counter * counter = streamCache(id);
155 
156  // generator information (always available)
158  iEvent.getByToken(genTag_, genInfo);
159  double weight = genInfo->weight();
160 
161  // table for gen info, always available
162  auto out = std::make_unique<nanoaod::FlatTable>(1, "genWeight", true);
163  out->setDoc("generator weight");
164  out->addColumnValue<float>("", weight, "generator weight", nanoaod::FlatTable::FloatColumn);
165  iEvent.put(std::move(out));
166 
167  // tables for LHE weights, may not be filled
168  std::unique_ptr<nanoaod::FlatTable> lheScaleTab, lhePdfTab, lheNamedTab;
169 
171  if (iEvent.getByToken(lheTag_, lheInfo)) {
172  // get the dynamic choice of weights
173  const DynamicWeightChoice * weightChoice = runCache(iEvent.getRun().index());
174  // go fill tables
175  fillLHEWeightTables(counter, weightChoice, weight, *lheInfo, lheScaleTab, lhePdfTab, lheNamedTab);
176  } else {
177  // minimal book-keeping of weights
178  counter->incGenOnly(weight);
179  // make dummy values
180  lheScaleTab.reset(new nanoaod::FlatTable(1, "LHEScaleWeights", true));
181  lhePdfTab.reset(new nanoaod::FlatTable(1, "LHEPdfWeights", true));
182  lheNamedTab.reset(new nanoaod::FlatTable(1, "LHENamedWeights", true));
183  if (!hasIssuedWarning_.exchange(true)) {
184  edm::LogWarning("LHETablesProducer") << "No LHEEventProduct, so there will be no LHE Tables\n";
185  }
186  }
187 
188  iEvent.put(std::move(lheScaleTab), "LHEScale");
189  iEvent.put(std::move(lhePdfTab), "LHEPdf");
190  iEvent.put(std::move(lheNamedTab), "LHENamed");
191  }
192 
194  Counter * counter,
195  const DynamicWeightChoice * weightChoice,
196  double genWeight,
197  const LHEEventProduct & lheProd,
198  std::unique_ptr<nanoaod::FlatTable> & outScale,
199  std::unique_ptr<nanoaod::FlatTable> & outPdf,
200  std::unique_ptr<nanoaod::FlatTable> & outNamed ) const
201  {
202  bool lheDebug = debug_.exchange(false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
203 
204  const std::vector<std::string> & scaleWeightIDs = weightChoice->scaleWeightIDs;
205  const std::vector<std::string> & pdfWeightIDs = weightChoice->pdfWeightIDs;
206 
207  double w0 = lheProd.originalXWGTUP();
208 
209  std::vector<double> wScale(scaleWeightIDs.size(), 1), wPDF(pdfWeightIDs.size(), 1), wNamed(namedWeightIDs_.size(), 1);
210  for (auto & weight : lheProd.weights()) {
211  if (lheDebug) printf("Weight %+9.5f rel %+9.5f for id %s\n", weight.wgt, weight.wgt/w0, weight.id.c_str());
212  // now we do it slowly, can be optimized
213  auto mScale = std::find(scaleWeightIDs.begin(), scaleWeightIDs.end(), weight.id);
214  if (mScale != scaleWeightIDs.end()) wScale[mScale-scaleWeightIDs.begin()] = weight.wgt/w0;
215 
216  auto mPDF = std::find(pdfWeightIDs.begin(), pdfWeightIDs.end(), weight.id);
217  if (mPDF != pdfWeightIDs.end()) wPDF[mPDF-pdfWeightIDs.begin()] = weight.wgt/w0;
218 
219  auto mNamed = std::find(namedWeightIDs_.begin(), namedWeightIDs_.end(), weight.id);
220  if (mNamed != namedWeightIDs_.end()) wNamed[mNamed-namedWeightIDs_.begin()] = weight.wgt/w0;
221  }
222 
223  outScale.reset(new nanoaod::FlatTable(wScale.size(), "LHEScaleWeight", false));
224  outScale->addColumn<float>("", wScale, weightChoice->scaleWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
225 
226  outPdf.reset(new nanoaod::FlatTable(wPDF.size(), "LHEPdfWeight", false));
227  outPdf->addColumn<float>("", wPDF, weightChoice->pdfWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
228 
229  outNamed.reset(new nanoaod::FlatTable(1, "LHEWeight", true));
230  outNamed->addColumnValue<float>("originalXWGTUP", lheProd.originalXWGTUP(), "Nominal event weight in the LHE file", nanoaod::FlatTable::FloatColumn);
231  for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) {
232  outNamed->addColumnValue<float>(namedWeightLabels_[i], wNamed[i], "LHE weight for id "+namedWeightIDs_[i]+", relative to nominal", nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
233  }
234 
235  counter->incLHE(genWeight, wScale, wPDF, wNamed);
236  }
237 
238  // create an empty counter
239  std::shared_ptr<DynamicWeightChoice> globalBeginRun(edm::Run const& iRun, edm::EventSetup const&) const override {
241 
242  bool lheDebug = debugRun_.exchange(false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
243  auto weightChoice = std::make_shared<DynamicWeightChoice>();
244 
245  // getByToken throws since we're not in the endRun (see https://github.com/cms-sw/cmssw/pull/18499)
246  //if (iRun.getByToken(lheRunTag_, lheInfo)) {
247  if (iRun.getByLabel(lheLabel_, lheInfo)) {
248  std::vector<ScaleVarWeight> scaleVariationIDs;
249  std::vector<PDFSetWeights> pdfSetWeightIDs;
250 
251  std::regex weightgroup("<weightgroup\\s+combine=\"(.*)\"\\s+(?:name|type)=\"(.*)\"\\s*>");
252  std::regex endweightgroup("</weightgroup>");
253  std::regex scalew("<weight\\s+(?:.*\\s+)?id=\"(\\d+)\">\\s*(?:lhapdf=\\d+|dyn=\\s*-?\\d+)?\\s*((?:mu[rR]|renscfact)=(\\S+)\\s+(?:mu[Ff]|facscfact)=(\\S+)(\\s+.*)?)</weight>");
254  std::regex pdfw("<weight\\s+id=\"(\\d+)\">\\s*(?:PDF set|lhapdf|PDF|pdfset)\\s*=\\s*(\\d+)\\s*(?:\\s.*)?</weight>");
255  std::regex pdfwOld("<weight\\s+(?:.*\\s+)?id=\"(\\d+)\">\\s*Member \\s*(\\d+)\\s*(?:.*)</weight>");
256  std::smatch groups;
257  for (auto iter=lheInfo->headers_begin(), end = lheInfo->headers_end(); iter != end; ++iter) {
258  if (iter->tag() != "initrwgt") {
259  if (lheDebug) std::cout << "Skipping LHE header with tag" << iter->tag() << std::endl;
260  continue;
261  }
262  if (lheDebug) std::cout << "Found LHE header with tag" << iter->tag() << std::endl;
263  const std::vector<std::string> & lines = iter->lines();
264  for (unsigned int iLine = 0, nLines = lines.size(); iLine < nLines; ++iLine) {
265  if (lheDebug) std::cout << lines[iLine];
266  if (std::regex_search(lines[iLine], groups, weightgroup)) {
267  std::string groupname = groups.str(2);
268  if (lheDebug) std::cout << ">>> Looks like the beginning of a weight group for '" << groupname << "'" << std::endl;
269  if (groupname.find("scale_variation") == 0 || groupname == "Central scale variation") {
270  if (lheDebug) std::cout << ">>> Looks like scale variation for theory uncertainties" << std::endl;
271  for ( ++iLine; iLine < nLines; ++iLine) {
272  if (lheDebug) std::cout << " " << lines[iLine];
273  if (std::regex_search(lines[iLine], groups, scalew)) {
274  if (lheDebug) std::cout << " >>> Scale weight " << groups[1].str() << " for " << groups[3].str() << " , " << groups[4].str() << " , " << groups[5].str() << std::endl;
275  scaleVariationIDs.emplace_back(groups.str(1), groups.str(2), groups.str(3), groups.str(4));
276  } else if (std::regex_search(lines[iLine], endweightgroup)) {
277  if (lheDebug) std::cout << ">>> Looks like the end of a weight group" << std::endl;
278  break;
279  } else if (std::regex_search(lines[iLine], weightgroup)) {
280  if (lheDebug) std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end of the group." << std::endl;
281  --iLine; // rewind by one, and go back to the outer loop
282  break;
283  }
284  }
285  } else if (groupname == "PDF_variation" || groupname.find("PDF_variation ") == 0) {
286  if (lheDebug) std::cout << ">>> Looks like a new-style block of PDF weights for one or more pdfs" << std::endl;
287  for ( ++iLine; iLine < nLines; ++iLine) {
288  if (lheDebug) std::cout << " " << lines[iLine];
289  if (std::regex_search(lines[iLine], groups, pdfw)) {
290  unsigned int lhaID = std::stoi(groups.str(2));
291  if (lheDebug) std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID << std::endl;
292  if (pdfSetWeightIDs.empty() || ! pdfSetWeightIDs.back().maybe_add(groups.str(1),lhaID)) {
293  pdfSetWeightIDs.emplace_back(groups.str(1),lhaID);
294  }
295  } else if (std::regex_search(lines[iLine], endweightgroup)) {
296  if (lheDebug) std::cout << ">>> Looks like the end of a weight group" << std::endl;
297  break;
298  } else if (std::regex_search(lines[iLine], weightgroup)) {
299  if (lheDebug) std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end of the group." << std::endl;
300  --iLine; // rewind by one, and go back to the outer loop
301  break;
302  }
303  }
304  } else if (groupname == "PDF_variation1" || groupname == "PDF_variation2") {
305  if (lheDebug) std::cout << ">>> Looks like a new-style block of PDF weights for multiple pdfs" << std::endl;
306  unsigned int lastid = 0;
307  for ( ++iLine; iLine < nLines; ++iLine) {
308  if (lheDebug) std::cout << " " << lines[iLine];
309  if (std::regex_search(lines[iLine], groups, pdfw)) {
310  unsigned int id = std::stoi(groups.str(1));
311  unsigned int lhaID = std::stoi(groups.str(2));
312  if (lheDebug) std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID << std::endl;
313  if (id != (lastid+1) || pdfSetWeightIDs.empty()) {
314  pdfSetWeightIDs.emplace_back(groups.str(1),lhaID);
315  } else {
316  pdfSetWeightIDs.back().add(groups.str(1),lhaID);
317  }
318  lastid = id;
319  } else if (std::regex_search(lines[iLine], endweightgroup)) {
320  if (lheDebug) std::cout << ">>> Looks like the end of a weight group" << std::endl;
321  break;
322  } else if (std::regex_search(lines[iLine], weightgroup)) {
323  if (lheDebug) std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end of the group." << std::endl;
324  --iLine; // rewind by one, and go back to the outer loop
325  break;
326  }
327  }
328  } else if (lhaNameToID_.find(groupname) != lhaNameToID_.end()) {
329  if (lheDebug) std::cout << ">>> Looks like an old-style PDF weight for an individual pdf" << std::endl;
330  unsigned int firstLhaID = lhaNameToID_.find(groupname)->second;
331  bool first = true;
332  for ( ++iLine; iLine < nLines; ++iLine) {
333  if (lheDebug) std::cout << " " << lines[iLine];
334  if (std::regex_search(lines[iLine], groups, pdfwOld)) {
335  unsigned int member = std::stoi(groups.str(2));
336  unsigned int lhaID = member+firstLhaID;
337  if (lheDebug) std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID << std::endl;
338  //if (member == 0) continue; // let's keep also the central value for now
339  if (first) {
340  pdfSetWeightIDs.emplace_back(groups.str(1),lhaID);
341  first = false;
342  } else {
343  pdfSetWeightIDs.back().add(groups.str(1),lhaID);
344  }
345  } else if (std::regex_search(lines[iLine], endweightgroup)) {
346  if (lheDebug) std::cout << ">>> Looks like the end of a weight group" << std::endl;
347  break;
348  } else if (std::regex_search(lines[iLine], weightgroup)) {
349  if (lheDebug) std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end of the group." << std::endl;
350  --iLine; // rewind by one, and go back to the outer loop
351  break;
352  }
353  }
354  } else {
355  for ( ++iLine; iLine < nLines; ++iLine) {
356  if (lheDebug) std::cout << " " << lines[iLine];
357  if (std::regex_search(lines[iLine], groups, endweightgroup)) {
358  if (lheDebug) std::cout << ">>> Looks like the end of a weight group" << std::endl;
359  break;
360  } else if (std::regex_search(lines[iLine], weightgroup)) {
361  if (lheDebug) std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end of the group." << std::endl;
362  --iLine; // rewind by one, and go back to the outer loop
363  break;
364  }
365  }
366  }
367  }
368  }
369  //std::cout << "============= END [ " << iter->tag() << " ] ============ \n\n" << std::endl;
370 
371  // ----- SCALE VARIATIONS -----
372  std::sort(scaleVariationIDs.begin(), scaleVariationIDs.end());
373  if (lheDebug) std::cout << "Found " << scaleVariationIDs.size() << " scale variations: " << std::endl;
374  std::stringstream scaleDoc; scaleDoc << "LHE scale variation weights (w_var / w_nominal); ";
375  for (unsigned int isw = 0, nsw = scaleVariationIDs.size(); isw < nsw; ++isw) {
376  const auto & sw = scaleVariationIDs[isw];
377  if (isw) scaleDoc << "; ";
378  scaleDoc << "[" << isw << "] is " << sw.label;
379  weightChoice->scaleWeightIDs.push_back(sw.wid);
380  if (lheDebug) printf(" id %s: scales ren = % .2f fact = % .2f text = %s\n", sw.wid.c_str(), sw.scales.first, sw.scales.second, sw.label.c_str());
381  }
382  if (!scaleVariationIDs.empty()) weightChoice->scaleWeightsDoc = scaleDoc.str();
383 
384  // ------ PDF VARIATIONS (take the preferred one) -----
385  if (lheDebug) {
386  std::cout << "Found " << pdfSetWeightIDs.size() << " PDF set errors: " << std::endl;
387  for (const auto & pw : pdfSetWeightIDs) printf("lhaIDs %6d - %6d (%3lu weights: %s, ... )\n", pw.lhaIDs.first, pw.lhaIDs.second, pw.wids.size(), pw.wids.front().c_str());
388  }
389 
390  std::stringstream pdfDoc; pdfDoc << "LHE pdf variation weights (w_var / w_nominal) for LHA IDs ";
391  bool found = false;
392  for (uint32_t lhaid : preferredPDFLHAIDs_) {
393  for (const auto & pw : pdfSetWeightIDs) {
394  if (pw.lhaIDs.first != lhaid && pw.lhaIDs.first != (lhaid+1)) continue; // sometimes the first weight is not saved if that PDF is the nominal one for the sample
395  if (pw.wids.size() == 1) continue; // only consider error sets
396  pdfDoc << pw.lhaIDs.first << " - " << pw.lhaIDs.second;
397  weightChoice->pdfWeightIDs = pw.wids;
398  if (maxPdfWeights_ < pw.wids.size()) {
399  weightChoice->pdfWeightIDs.resize(maxPdfWeights_); // drop some replicas
400  pdfDoc << ", truncated to the first " << maxPdfWeights_ << " replicas";
401  }
402  weightChoice->pdfWeightsDoc = pdfDoc.str();
403  found = true; break;
404  }
405  if (found) break;
406  }
407  }
408  }
409  return weightChoice;
410  }
411 
412 
413  // create an empty counter
414  std::unique_ptr<Counter> beginStream(edm::StreamID) const override {
415  return std::make_unique<Counter>();
416  }
417  // inizialize to zero at begin run
418  void streamBeginRun(edm::StreamID id, edm::Run const&, edm::EventSetup const&) const override {
419  streamCache(id)->clear();
420  }
421  // create an empty counter
422  std::shared_ptr<Counter> globalBeginRunSummary(edm::Run const&, edm::EventSetup const&) const override {
423  return std::make_shared<Counter>();
424  }
425  // add this stream to the summary
426  void streamEndRunSummary(edm::StreamID id, edm::Run const&, edm::EventSetup const&, Counter* runCounter) const override {
427  runCounter->merge(*streamCache(id));
428  }
429  // nothing to do per se
430  void globalEndRunSummary(edm::Run const&, edm::EventSetup const&, Counter* runCounter) const override {
431  }
432  // write the total to the run
433  void globalEndRunProduce(edm::Run& iRun, edm::EventSetup const&, Counter const* runCounter) const override {
434  auto out = std::make_unique<nanoaod::MergeableCounterTable>();
435  out->addInt("genEventCount", "event count", runCounter->num);
436  out->addFloat("genEventSumw", "sum of gen weights", runCounter->sumw);
437  out->addFloat("genEventSumw2", "sum of gen (weight^2)", runCounter->sumw2);
438 
439  double norm = runCounter->sumw ? 1.0/runCounter->sumw : 1;
440  auto sumScales = runCounter->sumScale; for (auto & val : sumScales) val *= norm;
441  out->addVFloat("LHEScaleSumw", "Sum of genEventWeight * LHEScaleWeight[i], divided by genEventSumw", sumScales);
442  auto sumPDFs = runCounter->sumPDF; for (auto & val : sumPDFs) val *= norm;
443  out->addVFloat("LHEPdfSumw", "Sum of genEventWeight * LHEPdfWeight[i], divided by genEventSumw", sumPDFs);
444  if (!runCounter->sumNamed.empty()) { // it could be empty if there's no LHE info in the sample
445  for (unsigned int i = 0, n = namedWeightLabels_.size(); i < n; ++i) {
446  out->addFloat("LHESumw_"+namedWeightLabels_[i], "Sum of genEventWeight * LHEWeight_"+namedWeightLabels_[i]+", divided by genEventSumw", runCounter->sumNamed[i] * norm);
447  }
448  }
449  iRun.put(std::move(out));
450  }
451  // nothing to do here
452  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override { }
453 
454  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions) {
456  desc.add<edm::InputTag>("genEvent", edm::InputTag("generator"))->setComment("tag for the GenEventInfoProduct, to get the main weight");
457  desc.add<edm::InputTag>("lheInfo", edm::InputTag("externalLHEProducer"))->setComment("tag for the LHE information (LHEEventProduct and LHERunInfoProduct)");
458 
460  prefpdf.add<std::string>("name");
461  prefpdf.add<uint32_t>("lhaid");
462  desc.addVPSet("preferredPDFs", prefpdf, std::vector<edm::ParameterSet>())->setComment("LHA PDF Ids of the preferred PDF sets, in order of preference (the first matching one will be used)");
463  desc.add<std::vector<std::string>>("namedWeightIDs")->setComment("set of LHA weight IDs for named LHE weights");
464  desc.add<std::vector<std::string>>("namedWeightLabels")->setComment("output names for the namedWeightIDs (in the same order)");
465  desc.add<int32_t>("lheWeightPrecision")->setComment("Number of bits in the mantissa for LHE weights");
466  desc.add<uint32_t>("maxPdfWeights")->setComment("Maximum number of PDF weights to save (to crop NN replicas)");
467  desc.addOptionalUntracked<bool>("debug")->setComment("dump out all LHE information for one event");
468  descriptions.add("genWeightsTable", desc);
469  }
470 
471 
472  protected:
477 
478  std::vector<uint32_t> preferredPDFLHAIDs_;
479  std::unordered_map<std::string,uint32_t> lhaNameToID_;
480  std::vector<std::string> namedWeightIDs_;
481  std::vector<std::string> namedWeightLabels_;
483  unsigned int maxPdfWeights_;
484 
485  mutable std::atomic<bool> debug_, debugRun_, hasIssuedWarning_;
486 };
487 
490 
double originalXWGTUP() const
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:248
void setComment(std::string const &value)
void streamBeginRun(edm::StreamID id, edm::Run const &, edm::EventSetup const &) const override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
std::unordered_map< std::string, uint32_t > lhaNameToID_
const double w
Definition: UKUtility.cc:23
void globalEndRun(edm::Run const &, edm::EventSetup const &) const override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
void globalEndRunProduce(edm::Run &iRun, edm::EventSetup const &, Counter const *runCounter) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void streamEndRunSummary(edm::StreamID id, edm::Run const &, edm::EventSetup const &, Counter *runCounter) const override
const edm::EDGetTokenT< GenEventInfoProduct > genTag_
std::shared_ptr< DynamicWeightChoice > globalBeginRun(edm::Run const &iRun, edm::EventSetup const &) const override
Definition: weight.py:1
headers_const_iterator headers_end() const
Run const & getRun() const
Definition: Event.cc:114
double weight() const
std::atomic< bool > hasIssuedWarning_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
uint16_t size_type
std::unique_ptr< Counter > beginStream(edm::StreamID) const override
const std::vector< WGT > & weights() const
std::function< unsigned int(align::ID)> Counter
int iEvent
Definition: GenABIO.cc:230
bool operator<(const FedChannelConnection &, const FedChannelConnection &)
headers_const_iterator headers_begin() const
void clear(CLHEP::HepGenMatrix &m)
Helper function: Reset all elements of a matrix to 0.
Definition: matutil.cc:167
const edm::EDGetTokenT< LHEEventProduct > lheTag_
RunIndex index() const
Definition: Run.cc:24
double f[11][100]
#define end
Definition: vmac.h:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
std::vector< std::string > namedWeightLabels_
GenWeightsTableProducer(edm::ParameterSet const &params)
std::shared_ptr< Counter > globalBeginRunSummary(edm::Run const &, edm::EventSetup const &) const override
def load(fileName)
Definition: svgfig.py:546
void put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Run.h:114
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< LHERunInfoProduct > lheRunTag_
std::vector< uint32_t > preferredPDFLHAIDs_
HLT enums.
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
ParameterDescriptionBase * addOptionalUntracked(U const &iLabel, T const &value)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< std::string > namedWeightIDs_
void globalEndRunSummary(edm::Run const &, edm::EventSetup const &, Counter *runCounter) const override
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510
Definition: Run.h:43
void fillLHEWeightTables(Counter *counter, const DynamicWeightChoice *weightChoice, double genWeight, const LHEEventProduct &lheProd, std::unique_ptr< nanoaod::FlatTable > &outScale, std::unique_ptr< nanoaod::FlatTable > &outPdf, std::unique_ptr< nanoaod::FlatTable > &outNamed) const
def merge(dictlist, TELL=False)
Definition: MatrixUtil.py:193