CMS 3D CMS Logo

GenWeightsTableProducer.cc
Go to the documentation of this file.
15 #include "boost/algorithm/string.hpp"
16 
17 #include <memory>
18 
19 #include <vector>
20 #include <unordered_map>
21 #include <iostream>
22 #include <regex>
23 
24 namespace {
26  struct Counter {
27  Counter() : num(0), sumw(0), sumw2(0), sumPDF(), sumScale(), sumRwgt(), sumNamed(), sumPS() {}
28 
29  // the counters
30  long long num;
31  long double sumw;
32  long double sumw2;
33  std::vector<long double> sumPDF, sumScale, sumRwgt, sumNamed, sumPS;
34 
35  void clear() {
36  num = 0;
37  sumw = 0;
38  sumw2 = 0;
39  sumPDF.clear();
40  sumScale.clear();
41  sumRwgt.clear();
42  sumNamed.clear(), sumPS.clear();
43  }
44 
45  // inc the counters
46  void incGenOnly(double w) {
47  num++;
48  sumw += w;
49  sumw2 += (w * w);
50  }
51 
52  void incPSOnly(double w0, const std::vector<double>& wPS) {
53  if (!wPS.empty()) {
54  if (sumPS.empty())
55  sumPS.resize(wPS.size(), 0);
56  for (unsigned int i = 0, n = wPS.size(); i < n; ++i)
57  sumPS[i] += (w0 * wPS[i]);
58  }
59  }
60 
61  void incLHE(double w0,
62  const std::vector<double>& wScale,
63  const std::vector<double>& wPDF,
64  const std::vector<double>& wRwgt,
65  const std::vector<double>& wNamed,
66  const std::vector<double>& wPS) {
67  // add up weights
68  incGenOnly(w0);
69  // then add up variations
70  if (!wScale.empty()) {
71  if (sumScale.empty())
72  sumScale.resize(wScale.size(), 0);
73  for (unsigned int i = 0, n = wScale.size(); i < n; ++i)
74  sumScale[i] += (w0 * wScale[i]);
75  }
76  if (!wPDF.empty()) {
77  if (sumPDF.empty())
78  sumPDF.resize(wPDF.size(), 0);
79  for (unsigned int i = 0, n = wPDF.size(); i < n; ++i)
80  sumPDF[i] += (w0 * wPDF[i]);
81  }
82  if (!wRwgt.empty()) {
83  if (sumRwgt.empty())
84  sumRwgt.resize(wRwgt.size(), 0);
85  for (unsigned int i = 0, n = wRwgt.size(); i < n; ++i)
86  sumRwgt[i] += (w0 * wRwgt[i]);
87  }
88  if (!wNamed.empty()) {
89  if (sumNamed.empty())
90  sumNamed.resize(wNamed.size(), 0);
91  for (unsigned int i = 0, n = wNamed.size(); i < n; ++i)
92  sumNamed[i] += (w0 * wNamed[i]);
93  }
94  incPSOnly(w0, wPS);
95  }
96 
97  void merge(const Counter& other) {
98  num += other.num;
99  sumw += other.sumw;
100  sumw2 += other.sumw2;
101  if (sumScale.empty() && !other.sumScale.empty())
102  sumScale.resize(other.sumScale.size(), 0);
103  if (sumPDF.empty() && !other.sumPDF.empty())
104  sumPDF.resize(other.sumPDF.size(), 0);
105  if (sumRwgt.empty() && !other.sumRwgt.empty())
106  sumRwgt.resize(other.sumRwgt.size(), 0);
107  if (sumNamed.empty() && !other.sumNamed.empty())
108  sumNamed.resize(other.sumNamed.size(), 0);
109  if (sumPS.empty() && !other.sumPS.empty())
110  sumPS.resize(other.sumPS.size(), 0);
111  if (!other.sumScale.empty())
112  for (unsigned int i = 0, n = sumScale.size(); i < n; ++i)
113  sumScale[i] += other.sumScale[i];
114  if (!other.sumPDF.empty())
115  for (unsigned int i = 0, n = sumPDF.size(); i < n; ++i)
116  sumPDF[i] += other.sumPDF[i];
117  if (!other.sumRwgt.empty())
118  for (unsigned int i = 0, n = sumRwgt.size(); i < n; ++i)
119  sumRwgt[i] += other.sumRwgt[i];
120  if (!other.sumNamed.empty())
121  for (unsigned int i = 0, n = sumNamed.size(); i < n; ++i)
122  sumNamed[i] += other.sumNamed[i];
123  if (!other.sumPS.empty())
124  for (unsigned int i = 0, n = sumPS.size(); i < n; ++i)
125  sumPS[i] += other.sumPS[i];
126  }
127  };
128 
129  struct CounterMap {
130  std::map<std::string, Counter> countermap;
131  Counter* active_el = nullptr;
132  std::string active_label = "";
133  void merge(const CounterMap& other) {
134  for (const auto& y : other.countermap)
135  countermap[y.first].merge(y.second);
136  active_el = nullptr;
137  }
138  void clear() {
139  for (auto x : countermap)
140  x.second.clear();
141  active_el = nullptr;
142  active_label = "";
143  }
144  void setLabel(std::string label) {
145  active_el = &(countermap[label]);
146  active_label = label;
147  }
148  void checkLabelSet() {
149  if (!active_el)
150  throw cms::Exception("LogicError", "Called CounterMap::get() before setting the active label\n");
151  }
152  Counter* get() {
153  checkLabelSet();
154  return active_el;
155  }
156  std::string& getLabel() {
157  checkLabelSet();
158  return active_label;
159  }
160  };
161 
163  struct DynamicWeightChoice {
164  // choice of LHE weights
165  // ---- scale ----
166  std::vector<std::string> scaleWeightIDs;
167  std::string scaleWeightsDoc;
168  // ---- pdf ----
169  std::vector<std::string> pdfWeightIDs;
170  std::string pdfWeightsDoc;
171  // ---- rwgt ----
172  std::vector<std::string> rwgtIDs;
173  std::string rwgtWeightDoc;
174  };
175 
176  struct DynamicWeightChoiceGenInfo {
177  // choice of LHE weights
178  // ---- scale ----
179  std::vector<unsigned int> scaleWeightIDs;
180  std::string scaleWeightsDoc;
181  // ---- pdf ----
182  std::vector<unsigned int> pdfWeightIDs;
183  std::string pdfWeightsDoc;
184  // ---- ps ----
185  std::vector<unsigned int> psWeightIDs;
186  };
187 
188  struct LumiCacheInfoHolder {
189  CounterMap countermap;
190  DynamicWeightChoiceGenInfo weightChoice;
191  void clear() {
192  countermap.clear();
193  weightChoice = DynamicWeightChoiceGenInfo();
194  }
195  };
196 
197  float stof_fortrancomp(const std::string& str) {
198  std::string::size_type match = str.find('d');
199  if (match != std::string::npos) {
200  std::string pre = str.substr(0, match);
201  std::string post = str.substr(match + 1);
202  return std::stof(pre) * std::pow(10.0f, std::stof(post));
203  } else {
204  return std::stof(str);
205  }
206  }
208  struct ScaleVarWeight {
209  std::string wid, label;
210  std::pair<float, float> scales;
211  ScaleVarWeight(const std::string& id, const std::string& text, const std::string& muR, const std::string& muF)
212  : wid(id), label(text), scales(stof_fortrancomp(muR), stof_fortrancomp(muF)) {}
213  bool operator<(const ScaleVarWeight& other) {
214  return (scales == other.scales ? wid < other.wid : scales < other.scales);
215  }
216  };
217  struct PDFSetWeights {
218  std::vector<std::string> wids;
219  std::pair<unsigned int, unsigned int> lhaIDs;
220  PDFSetWeights(const std::string& wid, unsigned int lhaID) : wids(1, wid), lhaIDs(lhaID, lhaID) {}
221  bool operator<(const PDFSetWeights& other) const { return lhaIDs < other.lhaIDs; }
222  void add(const std::string& wid, unsigned int lhaID) {
223  wids.push_back(wid);
224  lhaIDs.second = lhaID;
225  }
226  bool maybe_add(const std::string& wid, unsigned int lhaID) {
227  if (lhaID == lhaIDs.second + 1) {
228  lhaIDs.second++;
229  wids.push_back(wid);
230  return true;
231  } else {
232  return false;
233  }
234  }
235  };
236 } // namespace
237 
238 class GenWeightsTableProducer : public edm::global::EDProducer<edm::StreamCache<LumiCacheInfoHolder>,
239  edm::RunCache<DynamicWeightChoice>,
240  edm::RunSummaryCache<CounterMap>,
241  edm::EndRunProducer> {
242 public:
244  : genTag_(consumes<GenEventInfoProduct>(params.getParameter<edm::InputTag>("genEvent"))),
245  lheLabel_(params.getParameter<std::vector<edm::InputTag>>("lheInfo")),
246  lheTag_(edm::vector_transform(lheLabel_,
247  [this](const edm::InputTag& tag) { return mayConsume<LHEEventProduct>(tag); })),
248  lheRunTag_(edm::vector_transform(
249  lheLabel_, [this](const edm::InputTag& tag) { return mayConsume<LHERunInfoProduct, edm::InRun>(tag); })),
250  genLumiInfoHeadTag_(
251  mayConsume<GenLumiInfoHeader, edm::InLumi>(params.getParameter<edm::InputTag>("genLumiInfoHeader"))),
252  namedWeightIDs_(params.getParameter<std::vector<std::string>>("namedWeightIDs")),
253  namedWeightLabels_(params.getParameter<std::vector<std::string>>("namedWeightLabels")),
254  lheWeightPrecision_(params.getParameter<int32_t>("lheWeightPrecision")),
255  maxPdfWeights_(params.getParameter<uint32_t>("maxPdfWeights")),
256  keepAllPSWeights_(params.getParameter<bool>("keepAllPSWeights")),
257  debug_(params.getUntrackedParameter<bool>("debug", false)),
258  debugRun_(debug_.load()),
259  hasIssuedWarning_(false) {
260  produces<nanoaod::FlatTable>();
261  produces<std::string>("genModel");
262  produces<nanoaod::FlatTable>("LHEScale");
263  produces<nanoaod::FlatTable>("LHEPdf");
264  produces<nanoaod::FlatTable>("LHEReweighting");
265  produces<nanoaod::FlatTable>("LHENamed");
266  produces<nanoaod::FlatTable>("PS");
267  produces<nanoaod::MergeableCounterTable, edm::Transition::EndRun>();
268  if (namedWeightIDs_.size() != namedWeightLabels_.size()) {
269  throw cms::Exception("Configuration", "Size mismatch between namedWeightIDs & namedWeightLabels");
270  }
271  for (const edm::ParameterSet& pdfps : params.getParameter<std::vector<edm::ParameterSet>>("preferredPDFs")) {
272  const std::string& name = pdfps.getParameter<std::string>("name");
273  uint32_t lhaid = pdfps.getParameter<uint32_t>("lhaid");
274  preferredPDFLHAIDs_.push_back(lhaid);
275  lhaNameToID_[name] = lhaid;
276  lhaNameToID_[name + ".LHgrid"] = lhaid;
277  }
278  }
279 
281 
282  void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
283  // get my counter for weights
284  Counter* counter = streamCache(id)->countermap.get();
285 
286  // generator information (always available)
288  iEvent.getByToken(genTag_, genInfo);
289  double weight = genInfo->weight();
290 
291  // table for gen info, always available
292  auto out = std::make_unique<nanoaod::FlatTable>(1, "genWeight", true);
293  out->setDoc("generator weight");
294  out->addColumnValue<float>("", weight, "generator weight", nanoaod::FlatTable::FloatColumn);
295  iEvent.put(std::move(out));
296 
297  std::string model_label = streamCache(id)->countermap.getLabel();
298  auto outM = std::make_unique<std::string>((!model_label.empty()) ? std::string("GenModel_") + model_label : "");
299  iEvent.put(std::move(outM), "genModel");
300  bool getLHEweightsFromGenInfo = !model_label.empty();
301 
302  // tables for LHE weights, may not be filled
303  std::unique_ptr<nanoaod::FlatTable> lheScaleTab, lhePdfTab, lheRwgtTab, lheNamedTab;
304  std::unique_ptr<nanoaod::FlatTable> genPSTab;
305 
307  for (const auto& lheTag : lheTag_) {
308  iEvent.getByToken(lheTag, lheInfo);
309  if (lheInfo.isValid()) {
310  break;
311  }
312  }
313  if (lheInfo.isValid()) {
314  if (getLHEweightsFromGenInfo)
315  edm::LogWarning("LHETablesProducer")
316  << "Found both a LHEEventProduct and a GenLumiInfoHeader: will only save weights from LHEEventProduct.\n";
317  // get the dynamic choice of weights
318  const DynamicWeightChoice* weightChoice = runCache(iEvent.getRun().index());
319  // go fill tables
320  fillLHEWeightTables(
321  counter, weightChoice, weight, *lheInfo, *genInfo, lheScaleTab, lhePdfTab, lheRwgtTab, lheNamedTab, genPSTab);
322  } else if (getLHEweightsFromGenInfo) {
323  const DynamicWeightChoiceGenInfo* weightChoice = &(streamCache(id)->weightChoice);
324  fillLHEPdfWeightTablesFromGenInfo(
325  counter, weightChoice, weight, *genInfo, lheScaleTab, lhePdfTab, lheNamedTab, genPSTab);
326  lheRwgtTab = std::make_unique<nanoaod::FlatTable>(1, "LHEReweightingWeights", true);
327  //lheNamedTab = std::make_unique<nanoaod::FlatTable>(1, "LHENamedWeights", true);
328  //genPSTab = std::make_unique<nanoaod::FlatTable>(1, "PSWeight", true);
329  } else {
330  // Still try to add the PS weights
331  fillOnlyPSWeightTable(counter, weight, *genInfo, genPSTab);
332  // make dummy values
333  lheScaleTab = std::make_unique<nanoaod::FlatTable>(1, "LHEScaleWeights", true);
334  lhePdfTab = std::make_unique<nanoaod::FlatTable>(1, "LHEPdfWeights", true);
335  lheRwgtTab = std::make_unique<nanoaod::FlatTable>(1, "LHEReweightingWeights", true);
336  lheNamedTab = std::make_unique<nanoaod::FlatTable>(1, "LHENamedWeights", true);
337  if (!hasIssuedWarning_.exchange(true)) {
338  edm::LogWarning("LHETablesProducer") << "No LHEEventProduct, so there will be no LHE Tables\n";
339  }
340  }
341 
342  iEvent.put(std::move(lheScaleTab), "LHEScale");
343  iEvent.put(std::move(lhePdfTab), "LHEPdf");
344  iEvent.put(std::move(lheRwgtTab), "LHEReweighting");
345  iEvent.put(std::move(lheNamedTab), "LHENamed");
346  iEvent.put(std::move(genPSTab), "PS");
347  }
348 
350  const DynamicWeightChoice* weightChoice,
351  double genWeight,
352  const LHEEventProduct& lheProd,
353  const GenEventInfoProduct& genProd,
354  std::unique_ptr<nanoaod::FlatTable>& outScale,
355  std::unique_ptr<nanoaod::FlatTable>& outPdf,
356  std::unique_ptr<nanoaod::FlatTable>& outRwgt,
357  std::unique_ptr<nanoaod::FlatTable>& outNamed,
358  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
359  bool lheDebug = debug_.exchange(
360  false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
361 
362  const std::vector<std::string>& scaleWeightIDs = weightChoice->scaleWeightIDs;
363  const std::vector<std::string>& pdfWeightIDs = weightChoice->pdfWeightIDs;
364  const std::vector<std::string>& rwgtWeightIDs = weightChoice->rwgtIDs;
365 
366  double w0 = lheProd.originalXWGTUP();
367 
368  std::vector<double> wScale(scaleWeightIDs.size(), 1), wPDF(pdfWeightIDs.size(), 1), wRwgt(rwgtWeightIDs.size(), 1),
369  wNamed(namedWeightIDs_.size(), 1);
370  for (auto& weight : lheProd.weights()) {
371  if (lheDebug)
372  printf("Weight %+9.5f rel %+9.5f for id %s\n", weight.wgt, weight.wgt / w0, weight.id.c_str());
373  // now we do it slowly, can be optimized
374  auto mScale = std::find(scaleWeightIDs.begin(), scaleWeightIDs.end(), weight.id);
375  if (mScale != scaleWeightIDs.end())
376  wScale[mScale - scaleWeightIDs.begin()] = weight.wgt / w0;
377 
378  auto mPDF = std::find(pdfWeightIDs.begin(), pdfWeightIDs.end(), weight.id);
379  if (mPDF != pdfWeightIDs.end())
380  wPDF[mPDF - pdfWeightIDs.begin()] = weight.wgt / w0;
381 
382  auto mRwgt = std::find(rwgtWeightIDs.begin(), rwgtWeightIDs.end(), weight.id);
383  if (mRwgt != rwgtWeightIDs.end())
384  wRwgt[mRwgt - rwgtWeightIDs.begin()] = weight.wgt / w0;
385 
386  auto mNamed = std::find(namedWeightIDs_.begin(), namedWeightIDs_.end(), weight.id);
387  if (mNamed != namedWeightIDs_.end())
388  wNamed[mNamed - namedWeightIDs_.begin()] = weight.wgt / w0;
389  }
390 
391  std::size_t vectorSize =
392  (genProd.weights().size() > 2)
393  ? (keepAllPSWeights_ ? (genProd.weights().size() - 2)
394  : ((genProd.weights().size() == 14 || genProd.weights().size() == 46) ? 4 : 1))
395  : 1;
396  std::vector<double> wPS(vectorSize, 1);
397  std::string psWeightDocStr;
398  if (vectorSize > 1) {
399  double nominal = genProd.weights()[1]; // Called 'Baseline' in GenLumiInfoHeader
400  if (keepAllPSWeights_) {
401  for (std::size_t i = 0; i < vectorSize; i++) {
402  wPS[i] = (genProd.weights()[i + 2]) / nominal;
403  }
404  psWeightDocStr = "All PS weights (w_var / w_nominal)";
405  } else {
406  for (std::size_t i = 6; i < 10; i++) {
407  wPS[i - 6] = (genProd.weights()[i]) / nominal;
408  }
409  psWeightDocStr =
410  "PS weights (w_var / w_nominal); [0] is ISR=0.5 FSR=1; [1] is ISR=1 "
411  "FSR=0.5; [2] is ISR=2 FSR=1; [3] is ISR=1 FSR=2 ";
412  }
413  } else {
414  psWeightDocStr = "dummy PS weight (1.0) ";
415  }
416  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
417  outPS->addColumn<float>("", wPS, psWeightDocStr, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
418 
419  outScale = std::make_unique<nanoaod::FlatTable>(wScale.size(), "LHEScaleWeight", false);
420  outScale->addColumn<float>(
421  "", wScale, weightChoice->scaleWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
422 
423  outPdf = std::make_unique<nanoaod::FlatTable>(wPDF.size(), "LHEPdfWeight", false);
424  outPdf->addColumn<float>(
425  "", wPDF, weightChoice->pdfWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
426 
427  outRwgt = std::make_unique<nanoaod::FlatTable>(wRwgt.size(), "LHEReweightingWeight", false);
428  outRwgt->addColumn<float>(
429  "", wRwgt, weightChoice->rwgtWeightDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
430 
431  outNamed = std::make_unique<nanoaod::FlatTable>(1, "LHEWeight", true);
432  outNamed->addColumnValue<float>("originalXWGTUP",
433  lheProd.originalXWGTUP(),
434  "Nominal event weight in the LHE file",
436  for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) {
437  outNamed->addColumnValue<float>(namedWeightLabels_[i],
438  wNamed[i],
439  "LHE weight for id " + namedWeightIDs_[i] + ", relative to nominal",
441  lheWeightPrecision_);
442  }
443 
444  counter->incLHE(genWeight, wScale, wPDF, wRwgt, wNamed, wPS);
445  }
446 
448  const DynamicWeightChoiceGenInfo* weightChoice,
449  double genWeight,
450  const GenEventInfoProduct& genProd,
451  std::unique_ptr<nanoaod::FlatTable>& outScale,
452  std::unique_ptr<nanoaod::FlatTable>& outPdf,
453  std::unique_ptr<nanoaod::FlatTable>& outNamed,
454  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
455  const std::vector<unsigned int>& scaleWeightIDs = weightChoice->scaleWeightIDs;
456  const std::vector<unsigned int>& pdfWeightIDs = weightChoice->pdfWeightIDs;
457  const std::vector<unsigned int>& psWeightIDs = weightChoice->psWeightIDs;
458 
459  auto weights = genProd.weights();
460  double w0 = (weights.size() > 1) ? weights.at(1) : 1.;
461  double originalXWGTUP = (weights.size() > 1) ? weights.at(1) : 1.;
462 
463  std::vector<double> wScale, wPDF, wPS;
464  for (auto id : scaleWeightIDs)
465  wScale.push_back(weights.at(id) / w0);
466  for (auto id : pdfWeightIDs) {
467  wPDF.push_back(weights.at(id) / w0);
468  }
469  if (!psWeightIDs.empty()) {
470  double psNom =
471  weights.at(psWeightIDs.at(0)); // normalise PS weights by "Baseline", which should be the first entry
472  for (std::size_t i = 1; i < psWeightIDs.size(); i++)
473  wPS.push_back(weights.at(psWeightIDs.at(i)) / psNom);
474  } else
475  wPS.push_back(1.0);
476  std::string psWeightDocStr;
477  if (keepAllPSWeights_) {
478  psWeightDocStr = "All PS weights (w_var / w_nominal)";
479  } else {
480  psWeightDocStr =
481  "PS weights (w_var / w_nominal); [0] is ISR=0.5 FSR=1; [1] is ISR=1 "
482  "FSR=0.5; [2] is ISR=2 FSR=1; [3] is ISR=1 FSR=2 ";
483  }
484 
485  outScale = std::make_unique<nanoaod::FlatTable>(wScale.size(), "LHEScaleWeight", false);
486  outScale->addColumn<float>(
487  "", wScale, weightChoice->scaleWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
488 
489  outPdf = std::make_unique<nanoaod::FlatTable>(wPDF.size(), "LHEPdfWeight", false);
490  outPdf->addColumn<float>(
491  "", wPDF, weightChoice->pdfWeightsDoc, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
492 
493  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
494  outPS->addColumn<float>("",
495  wPS,
496  wPS.size() > 1 ? psWeightDocStr : "dummy PS weight (1.0) ",
498  lheWeightPrecision_);
499 
500  outNamed = std::make_unique<nanoaod::FlatTable>(1, "LHEWeight", true);
501  outNamed->addColumnValue<float>(
502  "originalXWGTUP", originalXWGTUP, "Nominal event weight in the LHE file", nanoaod::FlatTable::FloatColumn);
503  /*for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) {
504  outNamed->addColumnValue<float>(namedWeightLabels_[i], wNamed[i], "LHE weight for id "+namedWeightIDs_[i]+", relative to nominal", nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
505  }*/
506 
507  counter->incLHE(genWeight, wScale, wPDF, std::vector<double>(), std::vector<double>(), wPS);
508  }
509 
511  double genWeight,
512  const GenEventInfoProduct& genProd,
513  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
514  std::size_t vectorSize =
515  (genProd.weights().size() > 2)
516  ? (keepAllPSWeights_ ? (genProd.weights().size() - 2)
517  : ((genProd.weights().size() == 14 || genProd.weights().size() == 46) ? 4 : 1))
518  : 1;
519  std::vector<double> wPS(vectorSize, 1);
520  std::string psWeightDocStr;
521  if (vectorSize > 1) {
522  double nominal = genProd.weights()[1]; // Called 'Baseline' in GenLumiInfoHeader
523  if (keepAllPSWeights_) {
524  for (std::size_t i = 0; i < vectorSize; i++) {
525  wPS[i] = (genProd.weights()[i + 2]) / nominal;
526  }
527  psWeightDocStr = "All PS weights (w_var / w_nominal)";
528  } else {
529  for (std::size_t i = 6; i < 10; i++) {
530  wPS[i - 6] = (genProd.weights()[i]) / nominal;
531  }
532  psWeightDocStr =
533  "PS weights (w_var / w_nominal); [0] is ISR=0.5 FSR=1; [1] is ISR=1 "
534  "FSR=0.5; [2] is ISR=2 FSR=1; [3] is ISR=1 FSR=2 ";
535  }
536  } else {
537  psWeightDocStr = "dummy PS weight (1.0) ";
538  }
539  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
540  outPS->addColumn<float>("", wPS, psWeightDocStr, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);
541 
542  counter->incGenOnly(genWeight);
543  counter->incPSOnly(genWeight, wPS);
544  }
545 
546  // create an empty counter
547  std::shared_ptr<DynamicWeightChoice> globalBeginRun(edm::Run const& iRun, edm::EventSetup const&) const override {
549 
550  bool lheDebug = debugRun_.exchange(
551  false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
552  auto weightChoice = std::make_shared<DynamicWeightChoice>();
553 
554  // getByToken throws since we're not in the endRun (see https://github.com/cms-sw/cmssw/pull/18499)
555  //if (iRun.getByToken(lheRunTag_, lheInfo)) {
556  for (const auto& lheLabel : lheLabel_) {
557  iRun.getByLabel(lheLabel, lheInfo);
558  if (lheInfo.isValid()) {
559  break;
560  }
561  }
562  if (lheInfo.isValid()) {
563  std::vector<ScaleVarWeight> scaleVariationIDs;
564  std::vector<PDFSetWeights> pdfSetWeightIDs;
565  std::vector<std::string> lheReweighingIDs;
566 
567  std::regex weightgroupmg26x("<weightgroup\\s+(?:name|type)=\"(.*)\"\\s+combine=\"(.*)\"\\s*>");
568  std::regex weightgroup("<weightgroup\\s+combine=\"(.*)\"\\s+(?:name|type)=\"(.*)\"\\s*>");
569  std::regex weightgroupRwgt("<weightgroup\\s+(?:name|type)=\"(.*)\"\\s*>");
570  std::regex endweightgroup("</weightgroup>");
571  std::regex scalewmg26x(
572  "<weight\\s+(?:.*\\s+)?id=\"(\\d+)\"\\s*(?:lhapdf=\\d+|dyn=\\s*-?\\d+)?\\s*((?:[mM][uU][rR]|renscfact)=\"("
573  "\\S+)\"\\s+(?:[mM][uU][Ff]|facscfact)=\"(\\S+)\")(\\s+.*)?</weight>");
574  std::regex scalew(
575  "<weight\\s+(?:.*\\s+)?id=\"(\\d+)\">\\s*(?:lhapdf=\\d+|dyn=\\s*-?\\d+)?\\s*((?:mu[rR]|renscfact)=(\\S+)\\s+("
576  "?:mu[Ff]|facscfact)=(\\S+)(\\s+.*)?)</weight>");
577  std::regex pdfw(
578  "<weight\\s+id=\"(\\d+)\">\\s*(?:PDF set|lhapdf|PDF|pdfset)\\s*=\\s*(\\d+)\\s*(?:\\s.*)?</weight>");
579  std::regex pdfwOld("<weight\\s+(?:.*\\s+)?id=\"(\\d+)\">\\s*Member \\s*(\\d+)\\s*(?:.*)</weight>");
580  std::regex pdfwmg26x(
581  "<weight\\s+id=\"(\\d+)\"\\s*MUR=\"(?:\\S+)\"\\s*MUF=\"(?:\\S+)\"\\s*(?:PDF "
582  "set|lhapdf|PDF|pdfset)\\s*=\\s*\"(\\d+)\"\\s*>\\s*(?:PDF=(\\d+)\\s*MemberID=(\\d+))?\\s*(?:\\s.*)?</"
583  "weight>");
584  std::regex rwgt("<weight\\s+id=\"(.+)\">(.+)?(</weight>)?");
585  std::smatch groups;
586  for (auto iter = lheInfo->headers_begin(), end = lheInfo->headers_end(); iter != end; ++iter) {
587  if (iter->tag() != "initrwgt") {
588  if (lheDebug)
589  std::cout << "Skipping LHE header with tag" << iter->tag() << std::endl;
590  continue;
591  }
592  if (lheDebug)
593  std::cout << "Found LHE header with tag" << iter->tag() << std::endl;
594  std::vector<std::string> lines = iter->lines();
595  bool missed_weightgroup =
596  false; //Needed because in some of the samples ( produced with MG26X ) a small part of the header info is ordered incorrectly
597  bool ismg26x = false;
598  for (unsigned int iLine = 0, nLines = lines.size(); iLine < nLines;
599  ++iLine) { //First start looping through the lines to see which weightgroup pattern is matched
600  boost::replace_all(lines[iLine], "&lt;", "<");
601  boost::replace_all(lines[iLine], "&gt;", ">");
602  if (std::regex_search(lines[iLine], groups, weightgroupmg26x)) {
603  ismg26x = true;
604  }
605  }
606  for (unsigned int iLine = 0, nLines = lines.size(); iLine < nLines; ++iLine) {
607  if (lheDebug)
608  std::cout << lines[iLine];
609  if (std::regex_search(lines[iLine], groups, ismg26x ? weightgroupmg26x : weightgroup)) {
610  std::string groupname = groups.str(2);
611  if (ismg26x)
612  groupname = groups.str(1);
613  if (lheDebug)
614  std::cout << ">>> Looks like the beginning of a weight group for '" << groupname << "'" << std::endl;
615  if (groupname.find("scale_variation") == 0 || groupname == "Central scale variation") {
616  if (lheDebug)
617  std::cout << ">>> Looks like scale variation for theory uncertainties" << std::endl;
618  for (++iLine; iLine < nLines; ++iLine) {
619  if (lheDebug)
620  std::cout << " " << lines[iLine];
621  if (std::regex_search(lines[iLine], groups, ismg26x ? scalewmg26x : scalew)) {
622  if (lheDebug)
623  std::cout << " >>> Scale weight " << groups[1].str() << " for " << groups[3].str() << " , "
624  << groups[4].str() << " , " << groups[5].str() << std::endl;
625  scaleVariationIDs.emplace_back(groups.str(1), groups.str(2), groups.str(3), groups.str(4));
626  } else if (std::regex_search(lines[iLine], endweightgroup)) {
627  if (lheDebug)
628  std::cout << ">>> Looks like the end of a weight group" << std::endl;
629  if (!missed_weightgroup) {
630  break;
631  } else
632  missed_weightgroup = false;
633  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
634  if (lheDebug)
635  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
636  "of the group."
637  << std::endl;
638  if (ismg26x)
639  missed_weightgroup = true;
640  --iLine; // rewind by one, and go back to the outer loop
641  break;
642  }
643  }
644  } else if (groupname == "PDF_variation" || groupname.find("PDF_variation ") == 0) {
645  if (lheDebug)
646  std::cout << ">>> Looks like a new-style block of PDF weights for one or more pdfs" << std::endl;
647  for (++iLine; iLine < nLines; ++iLine) {
648  if (lheDebug)
649  std::cout << " " << lines[iLine];
650  if (std::regex_search(lines[iLine], groups, pdfw)) {
651  unsigned int lhaID = std::stoi(groups.str(2));
652  if (lheDebug)
653  std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID
654  << std::endl;
655  if (pdfSetWeightIDs.empty() || !pdfSetWeightIDs.back().maybe_add(groups.str(1), lhaID)) {
656  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
657  }
658  } else if (std::regex_search(lines[iLine], endweightgroup)) {
659  if (lheDebug)
660  std::cout << ">>> Looks like the end of a weight group" << std::endl;
661  if (!missed_weightgroup) {
662  break;
663  } else
664  missed_weightgroup = false;
665  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
666  if (lheDebug)
667  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
668  "of the group."
669  << std::endl;
670  if (ismg26x)
671  missed_weightgroup = true;
672  --iLine; // rewind by one, and go back to the outer loop
673  break;
674  }
675  }
676  } else if (groupname == "PDF_variation1" || groupname == "PDF_variation2") {
677  if (lheDebug)
678  std::cout << ">>> Looks like a new-style block of PDF weights for multiple pdfs" << std::endl;
679  unsigned int lastid = 0;
680  for (++iLine; iLine < nLines; ++iLine) {
681  if (lheDebug)
682  std::cout << " " << lines[iLine];
683  if (std::regex_search(lines[iLine], groups, pdfw)) {
684  unsigned int id = std::stoi(groups.str(1));
685  unsigned int lhaID = std::stoi(groups.str(2));
686  if (lheDebug)
687  std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID
688  << std::endl;
689  if (id != (lastid + 1) || pdfSetWeightIDs.empty()) {
690  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
691  } else {
692  pdfSetWeightIDs.back().add(groups.str(1), lhaID);
693  }
694  lastid = id;
695  } else if (std::regex_search(lines[iLine], endweightgroup)) {
696  if (lheDebug)
697  std::cout << ">>> Looks like the end of a weight group" << std::endl;
698  if (!missed_weightgroup) {
699  break;
700  } else
701  missed_weightgroup = false;
702  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
703  if (lheDebug)
704  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
705  "of the group."
706  << std::endl;
707  if (ismg26x)
708  missed_weightgroup = true;
709  --iLine; // rewind by one, and go back to the outer loop
710  break;
711  }
712  }
713  } else if (lhaNameToID_.find(groupname) != lhaNameToID_.end()) {
714  if (lheDebug)
715  std::cout << ">>> Looks like an old-style PDF weight for an individual pdf" << std::endl;
716  unsigned int firstLhaID = lhaNameToID_.find(groupname)->second;
717  bool first = true;
718  for (++iLine; iLine < nLines; ++iLine) {
719  if (lheDebug)
720  std::cout << " " << lines[iLine];
721  if (std::regex_search(lines[iLine], groups, ismg26x ? pdfwmg26x : pdfwOld)) {
722  unsigned int member = 0;
723  if (ismg26x == 0) {
724  member = std::stoi(groups.str(2));
725  } else {
726  if (!groups.str(4).empty()) {
727  member = std::stoi(groups.str(4));
728  }
729  }
730  unsigned int lhaID = member + firstLhaID;
731  if (lheDebug)
732  std::cout << " >>> PDF weight " << groups.str(1) << " for " << member << " = " << lhaID
733  << std::endl;
734  //if (member == 0) continue; // let's keep also the central value for now
735  if (first) {
736  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
737  first = false;
738  } else {
739  pdfSetWeightIDs.back().add(groups.str(1), lhaID);
740  }
741  } else if (std::regex_search(lines[iLine], endweightgroup)) {
742  if (lheDebug)
743  std::cout << ">>> Looks like the end of a weight group" << std::endl;
744  if (!missed_weightgroup) {
745  break;
746  } else
747  missed_weightgroup = false;
748  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
749  if (lheDebug)
750  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
751  "of the group."
752  << std::endl;
753  if (ismg26x)
754  missed_weightgroup = true;
755  --iLine; // rewind by one, and go back to the outer loop
756  break;
757  }
758  }
759  } else {
760  for (++iLine; iLine < nLines; ++iLine) {
761  if (lheDebug)
762  std::cout << " " << lines[iLine];
763  if (std::regex_search(lines[iLine], groups, endweightgroup)) {
764  if (lheDebug)
765  std::cout << ">>> Looks like the end of a weight group" << std::endl;
766  if (!missed_weightgroup) {
767  break;
768  } else
769  missed_weightgroup = false;
770  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
771  if (lheDebug)
772  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
773  "of the group."
774  << std::endl;
775  if (ismg26x)
776  missed_weightgroup = true;
777  --iLine; // rewind by one, and go back to the outer loop
778  break;
779  }
780  }
781  }
782  } else if (std::regex_search(lines[iLine], groups, weightgroupRwgt)) {
783  std::string groupname = groups.str(1);
784  if (groupname == "mg_reweighting") {
785  if (lheDebug)
786  std::cout << ">>> Looks like a LHE weights for reweighting" << std::endl;
787  for (++iLine; iLine < nLines; ++iLine) {
788  if (lheDebug)
789  std::cout << " " << lines[iLine];
790  if (std::regex_search(lines[iLine], groups, rwgt)) {
791  std::string rwgtID = groups.str(1);
792  if (lheDebug)
793  std::cout << " >>> LHE reweighting weight: " << rwgtID << std::endl;
794  if (std::find(lheReweighingIDs.begin(), lheReweighingIDs.end(), rwgtID) == lheReweighingIDs.end()) {
795  // we're only interested in the beggining of the block
796  lheReweighingIDs.emplace_back(rwgtID);
797  }
798  } else if (std::regex_search(lines[iLine], endweightgroup)) {
799  if (lheDebug)
800  std::cout << ">>> Looks like the end of a weight group" << std::endl;
801  if (!missed_weightgroup) {
802  break;
803  } else
804  missed_weightgroup = false;
805  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
806  if (lheDebug)
807  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
808  "of the group."
809  << std::endl;
810  if (ismg26x)
811  missed_weightgroup = true;
812  --iLine; // rewind by one, and go back to the outer loop
813  break;
814  }
815  }
816  }
817  }
818  }
819  //std::cout << "============= END [ " << iter->tag() << " ] ============ \n\n" << std::endl;
820 
821  // ----- SCALE VARIATIONS -----
822  std::sort(scaleVariationIDs.begin(), scaleVariationIDs.end());
823  if (lheDebug)
824  std::cout << "Found " << scaleVariationIDs.size() << " scale variations: " << std::endl;
825  std::stringstream scaleDoc;
826  scaleDoc << "LHE scale variation weights (w_var / w_nominal); ";
827  for (unsigned int isw = 0, nsw = scaleVariationIDs.size(); isw < nsw; ++isw) {
828  const auto& sw = scaleVariationIDs[isw];
829  if (isw)
830  scaleDoc << "; ";
831  scaleDoc << "[" << isw << "] is " << sw.label;
832  weightChoice->scaleWeightIDs.push_back(sw.wid);
833  if (lheDebug)
834  printf(" id %s: scales ren = % .2f fact = % .2f text = %s\n",
835  sw.wid.c_str(),
836  sw.scales.first,
837  sw.scales.second,
838  sw.label.c_str());
839  }
840  if (!scaleVariationIDs.empty())
841  weightChoice->scaleWeightsDoc = scaleDoc.str();
842 
843  // ------ PDF VARIATIONS (take the preferred one) -----
844  if (lheDebug) {
845  std::cout << "Found " << pdfSetWeightIDs.size() << " PDF set errors: " << std::endl;
846  for (const auto& pw : pdfSetWeightIDs)
847  printf("lhaIDs %6d - %6d (%3lu weights: %s, ... )\n",
848  pw.lhaIDs.first,
849  pw.lhaIDs.second,
850  pw.wids.size(),
851  pw.wids.front().c_str());
852  }
853 
854  // ------ LHE REWEIGHTING -------
855  if (lheDebug) {
856  std::cout << "Found " << lheReweighingIDs.size() << " reweighting weights" << std::endl;
857  }
858  std::copy(lheReweighingIDs.begin(), lheReweighingIDs.end(), std::back_inserter(weightChoice->rwgtIDs));
859 
860  std::stringstream pdfDoc;
861  pdfDoc << "LHE pdf variation weights (w_var / w_nominal) for LHA IDs ";
862  bool found = false;
863  for (uint32_t lhaid : preferredPDFLHAIDs_) {
864  for (const auto& pw : pdfSetWeightIDs) {
865  if (pw.lhaIDs.first != lhaid && pw.lhaIDs.first != (lhaid + 1))
866  continue; // sometimes the first weight is not saved if that PDF is the nominal one for the sample
867  if (pw.wids.size() == 1)
868  continue; // only consider error sets
869  pdfDoc << pw.lhaIDs.first << " - " << pw.lhaIDs.second;
870  weightChoice->pdfWeightIDs = pw.wids;
871  if (maxPdfWeights_ < pw.wids.size()) {
872  weightChoice->pdfWeightIDs.resize(maxPdfWeights_); // drop some replicas
873  pdfDoc << ", truncated to the first " << maxPdfWeights_ << " replicas";
874  }
875  weightChoice->pdfWeightsDoc = pdfDoc.str();
876  found = true;
877  break;
878  }
879  if (found)
880  break;
881  }
882  }
883  }
884  return weightChoice;
885  }
886 
887  // create an empty counter
888  std::unique_ptr<LumiCacheInfoHolder> beginStream(edm::StreamID) const override {
889  return std::make_unique<LumiCacheInfoHolder>();
890  }
891  // inizialize to zero at begin run
892  void streamBeginRun(edm::StreamID id, edm::Run const&, edm::EventSetup const&) const override {
893  streamCache(id)->clear();
894  }
896  edm::LuminosityBlock const& lumiBlock,
897  edm::EventSetup const& eventSetup) const override {
898  auto counterMap = &(streamCache(id)->countermap);
899  edm::Handle<GenLumiInfoHeader> genLumiInfoHead;
900  lumiBlock.getByToken(genLumiInfoHeadTag_, genLumiInfoHead);
901  if (!genLumiInfoHead.isValid())
902  edm::LogWarning("LHETablesProducer")
903  << "No GenLumiInfoHeader product found, will not fill generator model string.\n";
905  if (genLumiInfoHead.isValid()) {
906  label = genLumiInfoHead->configDescription();
907  boost::replace_all(label, "-", "_");
908  boost::replace_all(label, "/", "_");
909  }
910  counterMap->setLabel(label);
911 
912  if (genLumiInfoHead.isValid()) {
913  auto weightChoice = &(streamCache(id)->weightChoice);
914 
915  std::vector<ScaleVarWeight> scaleVariationIDs;
916  std::vector<PDFSetWeights> pdfSetWeightIDs;
917  weightChoice->psWeightIDs.clear();
918 
919  std::regex scalew("LHE,\\s+id\\s+=\\s+(\\d+),\\s+(.+)\\,\\s+mur=(\\S+)\\smuf=(\\S+)");
920  std::regex pdfw("LHE,\\s+id\\s+=\\s+(\\d+),\\s+(.+),\\s+Member\\s+(\\d+)\\s+of\\ssets\\s+(\\w+\\b)");
921  std::smatch groups;
922  auto weightNames = genLumiInfoHead->weightNames();
923  std::unordered_map<std::string, uint32_t> knownPDFSetsFromGenInfo_;
924  unsigned int weightIter = 0;
925  for (const auto& line : weightNames) {
926  if (std::regex_search(line, groups, scalew)) { // scale variation
927  auto id = groups.str(1);
928  auto group = groups.str(2);
929  auto mur = groups.str(3);
930  auto muf = groups.str(4);
931  if (group.find("Central scale variation") != std::string::npos)
932  scaleVariationIDs.emplace_back(groups.str(1), groups.str(2), groups.str(3), groups.str(4));
933  } else if (std::regex_search(line, groups, pdfw)) { // PDF variation
934  auto id = groups.str(1);
935  auto group = groups.str(2);
936  auto memberid = groups.str(3);
937  auto pdfset = groups.str(4);
938  if (group.find(pdfset) != std::string::npos) {
939  if (knownPDFSetsFromGenInfo_.find(pdfset) == knownPDFSetsFromGenInfo_.end()) {
940  knownPDFSetsFromGenInfo_[pdfset] = std::atoi(id.c_str());
941  pdfSetWeightIDs.emplace_back(id, std::atoi(id.c_str()));
942  } else
943  pdfSetWeightIDs.back().add(id, std::atoi(id.c_str()));
944  }
945  } else if (line.find("Baseline") != std::string::npos || line.find("isr") != std::string::npos ||
946  line.find("fsr") != std::string::npos) {
947  if (keepAllPSWeights_ || line.find("Baseline") != std::string::npos ||
948  line.find("Def") != std::string::npos) {
949  weightChoice->psWeightIDs.push_back(weightIter); // PS variations
950  }
951  }
952  weightIter++;
953  }
954 
955  weightChoice->scaleWeightIDs.clear();
956  weightChoice->pdfWeightIDs.clear();
957 
958  std::sort(scaleVariationIDs.begin(), scaleVariationIDs.end());
959  std::stringstream scaleDoc;
960  scaleDoc << "LHE scale variation weights (w_var / w_nominal); ";
961  for (unsigned int isw = 0, nsw = scaleVariationIDs.size(); isw < nsw; ++isw) {
962  const auto& sw = scaleVariationIDs[isw];
963  if (isw)
964  scaleDoc << "; ";
965  scaleDoc << "[" << isw << "] is " << sw.label;
966  weightChoice->scaleWeightIDs.push_back(std::atoi(sw.wid.c_str()));
967  }
968  if (!scaleVariationIDs.empty())
969  weightChoice->scaleWeightsDoc = scaleDoc.str();
970  std::stringstream pdfDoc;
971  pdfDoc << "LHE pdf variation weights (w_var / w_nominal) for LHA names ";
972  bool found = false;
973  for (const auto& pw : pdfSetWeightIDs) {
974  if (pw.wids.size() == 1)
975  continue; // only consider error sets
976  for (const auto& wantedpdf : lhaNameToID_) {
977  auto pdfname = wantedpdf.first;
978  if (knownPDFSetsFromGenInfo_.find(pdfname) == knownPDFSetsFromGenInfo_.end())
979  continue;
980  uint32_t lhaid = knownPDFSetsFromGenInfo_.at(pdfname);
981  if (pw.lhaIDs.first != lhaid)
982  continue;
983  pdfDoc << pdfname;
984  for (const auto& x : pw.wids)
985  weightChoice->pdfWeightIDs.push_back(std::atoi(x.c_str()));
986  if (maxPdfWeights_ < pw.wids.size()) {
987  weightChoice->pdfWeightIDs.resize(maxPdfWeights_); // drop some replicas
988  pdfDoc << ", truncated to the first " << maxPdfWeights_ << " replicas";
989  }
990  weightChoice->pdfWeightsDoc = pdfDoc.str();
991  found = true;
992  break;
993  }
994  if (found)
995  break;
996  }
997  }
998  }
999  // create an empty counter
1000  std::shared_ptr<CounterMap> globalBeginRunSummary(edm::Run const&, edm::EventSetup const&) const override {
1001  return std::make_shared<CounterMap>();
1002  }
1003  // add this stream to the summary
1005  edm::Run const&,
1006  edm::EventSetup const&,
1007  CounterMap* runCounterMap) const override {
1008  runCounterMap->merge(streamCache(id)->countermap);
1009  }
1010  // nothing to do per se
1011  void globalEndRunSummary(edm::Run const&, edm::EventSetup const&, CounterMap* runCounterMap) const override {}
1012  // write the total to the run
1013  void globalEndRunProduce(edm::Run& iRun, edm::EventSetup const&, CounterMap const* runCounterMap) const override {
1014  auto out = std::make_unique<nanoaod::MergeableCounterTable>();
1015 
1016  for (const auto& x : runCounterMap->countermap) {
1017  auto runCounter = &(x.second);
1018  std::string label = (!x.first.empty()) ? (std::string("_") + x.first) : "";
1019  std::string doclabel = (!x.first.empty()) ? (std::string(", for model label ") + x.first) : "";
1020 
1021  out->addInt("genEventCount" + label, "event count" + doclabel, runCounter->num);
1022  out->addFloat("genEventSumw" + label, "sum of gen weights" + doclabel, runCounter->sumw);
1023  out->addFloat("genEventSumw2" + label, "sum of gen (weight^2)" + doclabel, runCounter->sumw2);
1024 
1025  double norm = runCounter->sumw ? 1.0 / runCounter->sumw : 1;
1026  auto sumScales = runCounter->sumScale;
1027  for (auto& val : sumScales)
1028  val *= norm;
1029  out->addVFloatWithNorm("LHEScaleSumw" + label,
1030  "Sum of genEventWeight * LHEScaleWeight[i], divided by genEventSumw" + doclabel,
1031  sumScales,
1032  runCounter->sumw);
1033  auto sumPDFs = runCounter->sumPDF;
1034  for (auto& val : sumPDFs)
1035  val *= norm;
1036  out->addVFloatWithNorm("LHEPdfSumw" + label,
1037  "Sum of genEventWeight * LHEPdfWeight[i], divided by genEventSumw" + doclabel,
1038  sumPDFs,
1039  runCounter->sumw);
1040  if (!runCounter->sumRwgt.empty()) {
1041  auto sumRwgts = runCounter->sumRwgt;
1042  for (auto& val : sumRwgts)
1043  val *= norm;
1044  out->addVFloatWithNorm("LHEReweightingSumw" + label,
1045  "Sum of genEventWeight * LHEReweightingWeight[i], divided by genEventSumw" + doclabel,
1046  sumRwgts,
1047  runCounter->sumw);
1048  }
1049  if (!runCounter->sumNamed.empty()) { // it could be empty if there's no LHE info in the sample
1050  for (unsigned int i = 0, n = namedWeightLabels_.size(); i < n; ++i) {
1051  out->addFloatWithNorm(
1052  "LHESumw_" + namedWeightLabels_[i] + label,
1053  "Sum of genEventWeight * LHEWeight_" + namedWeightLabels_[i] + ", divided by genEventSumw" + doclabel,
1054  runCounter->sumNamed[i] * norm,
1055  runCounter->sumw);
1056  }
1057  }
1058  }
1059  iRun.put(std::move(out));
1060  }
1061  // nothing to do here
1062  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override {}
1063 
1066  desc.add<edm::InputTag>("genEvent", edm::InputTag("generator"))
1067  ->setComment("tag for the GenEventInfoProduct, to get the main weight");
1068  desc.add<edm::InputTag>("genLumiInfoHeader", edm::InputTag("generator"))
1069  ->setComment("tag for the GenLumiInfoProduct, to get the model string");
1070  desc.add<std::vector<edm::InputTag>>("lheInfo", std::vector<edm::InputTag>{{"externalLHEProducer"}, {"source"}})
1071  ->setComment("tag(s) for the LHE information (LHEEventProduct and LHERunInfoProduct)");
1072 
1074  prefpdf.add<std::string>("name");
1075  prefpdf.add<uint32_t>("lhaid");
1076  desc.addVPSet("preferredPDFs", prefpdf, std::vector<edm::ParameterSet>())
1077  ->setComment(
1078  "LHA PDF Ids of the preferred PDF sets, in order of preference (the first matching one will be used)");
1079  desc.add<std::vector<std::string>>("namedWeightIDs")->setComment("set of LHA weight IDs for named LHE weights");
1080  desc.add<std::vector<std::string>>("namedWeightLabels")
1081  ->setComment("output names for the namedWeightIDs (in the same order)");
1082  desc.add<int32_t>("lheWeightPrecision")->setComment("Number of bits in the mantissa for LHE weights");
1083  desc.add<uint32_t>("maxPdfWeights")->setComment("Maximum number of PDF weights to save (to crop NN replicas)");
1084  desc.add<bool>("keepAllPSWeights")->setComment("Store all PS weights found");
1085  desc.addOptionalUntracked<bool>("debug")->setComment("dump out all LHE information for one event");
1086  descriptions.add("genWeightsTable", desc);
1087  }
1088 
1089 protected:
1091  const std::vector<edm::InputTag> lheLabel_;
1092  const std::vector<edm::EDGetTokenT<LHEEventProduct>> lheTag_;
1093  const std::vector<edm::EDGetTokenT<LHERunInfoProduct>> lheRunTag_;
1095 
1096  std::vector<uint32_t> preferredPDFLHAIDs_;
1097  std::unordered_map<std::string, uint32_t> lhaNameToID_;
1098  std::vector<std::string> namedWeightIDs_;
1099  std::vector<std::string> namedWeightLabels_;
1101  unsigned int maxPdfWeights_;
1103 
1104  mutable std::atomic<bool> debug_, debugRun_, hasIssuedWarning_;
1105 };
1106 
const std::vector< edm::EDGetTokenT< LHERunInfoProduct > > lheRunTag_
double originalXWGTUP() const
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:303
void setComment(std::string const &value)
std::unique_ptr< LumiCacheInfoHolder > beginStream(edm::StreamID) const override
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:137
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
void streamBeginLuminosityBlock(edm::StreamID id, edm::LuminosityBlock const &lumiBlock, edm::EventSetup const &eventSetup) const override
const double w
Definition: UKUtility.cc:23
def copy(args, dbName)
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:579
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::EDGetTokenT< GenEventInfoProduct > genTag_
const std::vector< std::string > & weightNames() const
void fillLHEPdfWeightTablesFromGenInfo(Counter *counter, const DynamicWeightChoiceGenInfo *weightChoice, double genWeight, const GenEventInfoProduct &genProd, std::unique_ptr< nanoaod::FlatTable > &outScale, std::unique_ptr< nanoaod::FlatTable > &outPdf, std::unique_ptr< nanoaod::FlatTable > &outNamed, std::unique_ptr< nanoaod::FlatTable > &outPS) const
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_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::unordered_map< std::string, uint32_t > lhaNameToID_
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
void fillLHEWeightTables(Counter *counter, const DynamicWeightChoice *weightChoice, double genWeight, const LHEEventProduct &lheProd, const GenEventInfoProduct &genProd, std::unique_ptr< nanoaod::FlatTable > &outScale, std::unique_ptr< nanoaod::FlatTable > &outPdf, std::unique_ptr< nanoaod::FlatTable > &outRwgt, std::unique_ptr< nanoaod::FlatTable > &outNamed, std::unique_ptr< nanoaod::FlatTable > &outPS) const
void globalEndRunSummary(edm::Run const &, edm::EventSetup const &, CounterMap *runCounterMap) const override
const std::vector< edm::InputTag > lheLabel_
const std::vector< WGT > & weights() const
std::function< unsigned int(align::ID)> Counter
std::shared_ptr< CounterMap > globalBeginRunSummary(edm::Run const &, edm::EventSetup const &) const override
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
RunIndex index() const
Definition: Run.cc:24
double f[11][100]
#define end
Definition: vmac.h:39
void fillOnlyPSWeightTable(Counter *counter, double genWeight, const GenEventInfoProduct &genProd, std::unique_ptr< nanoaod::FlatTable > &outPS) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
const std::string & configDescription() const
std::vector< std::string > namedWeightLabels_
GenWeightsTableProducer(edm::ParameterSet const &params)
void put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Run.h:115
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::vector< edm::EDGetTokenT< LHEEventProduct > > lheTag_
void streamEndRunSummary(edm::StreamID id, edm::Run const &, edm::EventSetup const &, CounterMap *runCounterMap) const override
void globalEndRunProduce(edm::Run &iRun, edm::EventSetup const &, CounterMap const *runCounterMap) const override
std::vector< uint32_t > preferredPDFLHAIDs_
HLT enums.
std::vector< double > & weights()
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_
#define str(s)
const edm::EDGetTokenT< GenLumiInfoHeader > genLumiInfoHeadTag_
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:44
def merge(dictlist, TELL=False)
Definition: MatrixUtil.py:193