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> defPSWeightIDs = {6, 7, 8, 9};
186  std::vector<unsigned int> defPSWeightIDs_alt = {27, 5, 26, 4};
187  bool matchPS_alt = false;
188  std::vector<unsigned int> psWeightIDs;
189  unsigned int psBaselineID = 1;
190  std::string psWeightsDoc;
191 
192  void setMissingWeight(int idx) { psWeightIDs[idx] = (matchPS_alt) ? defPSWeightIDs_alt[idx] : defPSWeightIDs[idx]; }
193 
194  bool empty() const { return scaleWeightIDs.empty() && pdfWeightIDs.empty() && psWeightIDs.empty(); }
195  };
196 
197  struct LumiCacheInfoHolder {
198  CounterMap countermap;
199  DynamicWeightChoiceGenInfo weightChoice;
200  void clear() {
201  countermap.clear();
202  weightChoice = DynamicWeightChoiceGenInfo();
203  }
204  };
205 
206  float stof_fortrancomp(const std::string& str) {
207  std::string::size_type match = str.find('d');
208  if (match != std::string::npos) {
209  std::string pre = str.substr(0, match);
210  std::string post = str.substr(match + 1);
211  return std::stof(pre) * std::pow(10.0f, std::stof(post));
212  } else {
213  return std::stof(str);
214  }
215  }
217  struct ScaleVarWeight {
218  std::string wid, label;
219  std::pair<float, float> scales;
220  ScaleVarWeight(const std::string& id, const std::string& text, const std::string& muR, const std::string& muF)
221  : wid(id), label(text), scales(stof_fortrancomp(muR), stof_fortrancomp(muF)) {}
222  bool operator<(const ScaleVarWeight& other) {
223  return (scales == other.scales ? wid < other.wid : scales < other.scales);
224  }
225  };
226  struct PDFSetWeights {
227  std::vector<std::string> wids;
228  std::pair<unsigned int, unsigned int> lhaIDs;
229  PDFSetWeights(const std::string& wid, unsigned int lhaID) : wids(1, wid), lhaIDs(lhaID, lhaID) {}
230  bool operator<(const PDFSetWeights& other) const { return lhaIDs < other.lhaIDs; }
231  void add(const std::string& wid, unsigned int lhaID) {
232  wids.push_back(wid);
233  lhaIDs.second = lhaID;
234  }
235  bool maybe_add(const std::string& wid, unsigned int lhaID) {
236  if (lhaID == lhaIDs.second + 1) {
237  lhaIDs.second++;
238  wids.push_back(wid);
239  return true;
240  } else {
241  return false;
242  }
243  }
244  };
245 } // namespace
246 
247 class GenWeightsTableProducer : public edm::global::EDProducer<edm::StreamCache<LumiCacheInfoHolder>,
248  edm::RunCache<DynamicWeightChoice>,
249  edm::RunSummaryCache<CounterMap>,
250  edm::EndRunProducer> {
251 public:
253  : genTag_(consumes<GenEventInfoProduct>(params.getParameter<edm::InputTag>("genEvent"))),
254  lheLabel_(params.getParameter<std::vector<edm::InputTag>>("lheInfo")),
256  [this](const edm::InputTag& tag) { return mayConsume<LHEEventProduct>(tag); })),
258  lheLabel_, [this](const edm::InputTag& tag) { return mayConsume<LHERunInfoProduct, edm::InRun>(tag); })),
260  mayConsume<GenLumiInfoHeader, edm::InLumi>(params.getParameter<edm::InputTag>("genLumiInfoHeader"))),
261  namedWeightIDs_(params.getParameter<std::vector<std::string>>("namedWeightIDs")),
262  namedWeightLabels_(params.getParameter<std::vector<std::string>>("namedWeightLabels")),
263  lheWeightPrecision_(params.getParameter<int32_t>("lheWeightPrecision")),
264  maxPdfWeights_(params.getParameter<uint32_t>("maxPdfWeights")),
265  keepAllPSWeights_(params.getParameter<bool>("keepAllPSWeights")),
266  debug_(params.getUntrackedParameter<bool>("debug", false)),
267  debugRun_(debug_.load()),
268  hasIssuedWarning_(false),
269  psWeightWarning_(false) {
270  produces<nanoaod::FlatTable>();
271  produces<std::string>("genModel");
272  produces<nanoaod::FlatTable>("LHEScale");
273  produces<nanoaod::FlatTable>("LHEPdf");
274  produces<nanoaod::FlatTable>("LHEReweighting");
275  produces<nanoaod::FlatTable>("LHENamed");
276  produces<nanoaod::FlatTable>("PS");
277  produces<nanoaod::MergeableCounterTable, edm::Transition::EndRun>();
278  if (namedWeightIDs_.size() != namedWeightLabels_.size()) {
279  throw cms::Exception("Configuration", "Size mismatch between namedWeightIDs & namedWeightLabels");
280  }
281  for (const edm::ParameterSet& pdfps : params.getParameter<std::vector<edm::ParameterSet>>("preferredPDFs")) {
282  const std::string& name = pdfps.getParameter<std::string>("name");
283  uint32_t lhaid = pdfps.getParameter<uint32_t>("lhaid");
284  preferredPDFLHAIDs_.push_back(lhaid);
286  lhaNameToID_[name + ".LHgrid"] = lhaid;
287  }
288  }
289 
291 
292  void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
293  // get my counter for weights
294  Counter* counter = streamCache(id)->countermap.get();
295 
296  // generator information (always available)
298  iEvent.getByToken(genTag_, genInfo);
299  double weight = genInfo->weight();
300 
301  // table for gen info, always available
302  auto out = std::make_unique<nanoaod::FlatTable>(1, "genWeight", true);
303  out->setDoc("generator weight");
304  out->addColumnValue<float>("", weight, "generator weight");
305  iEvent.put(std::move(out));
306 
307  std::string model_label = streamCache(id)->countermap.getLabel();
308  auto outM = std::make_unique<std::string>((!model_label.empty()) ? std::string("GenModel_") + model_label : "");
309  iEvent.put(std::move(outM), "genModel");
310  bool getLHEweightsFromGenInfo = !model_label.empty();
311 
312  // tables for LHE weights, may not be filled
313  std::unique_ptr<nanoaod::FlatTable> lheScaleTab, lhePdfTab, lheRwgtTab, lheNamedTab;
314  std::unique_ptr<nanoaod::FlatTable> genPSTab;
315 
317  for (const auto& lheTag : lheTag_) {
318  iEvent.getByToken(lheTag, lheInfo);
319  if (lheInfo.isValid()) {
320  break;
321  }
322  }
323 
324  const auto genWeightChoice = &(streamCache(id)->weightChoice);
325  if (lheInfo.isValid()) {
326  if (getLHEweightsFromGenInfo && !hasIssuedWarning_.exchange(true))
327  edm::LogWarning("LHETablesProducer")
328  << "Found both a LHEEventProduct and a GenLumiInfoHeader: will only save weights from LHEEventProduct.\n";
329  // get the dynamic choice of weights
330  const DynamicWeightChoice* weightChoice = runCache(iEvent.getRun().index());
331  // go fill tables
333  weightChoice,
334  genWeightChoice,
335  weight,
336  *lheInfo,
337  *genInfo,
338  lheScaleTab,
339  lhePdfTab,
340  lheRwgtTab,
341  lheNamedTab,
342  genPSTab);
343  } else if (getLHEweightsFromGenInfo) {
345  counter, genWeightChoice, weight, *genInfo, lheScaleTab, lhePdfTab, lheNamedTab, genPSTab);
346  lheRwgtTab = std::make_unique<nanoaod::FlatTable>(1, "LHEReweightingWeights", true);
347  //lheNamedTab.reset(new nanoaod::FlatTable(1, "LHENamedWeights", true));
348  //genPSTab.reset(new nanoaod::FlatTable(1, "PSWeight", true));
349  } else {
350  // Still try to add the PS weights
351  fillOnlyPSWeightTable(counter, genWeightChoice, weight, *genInfo, genPSTab);
352  // make dummy values
353  lheScaleTab = std::make_unique<nanoaod::FlatTable>(1, "LHEScaleWeights", true);
354  lhePdfTab = std::make_unique<nanoaod::FlatTable>(1, "LHEPdfWeights", true);
355  lheRwgtTab = std::make_unique<nanoaod::FlatTable>(1, "LHEReweightingWeights", true);
356  lheNamedTab = std::make_unique<nanoaod::FlatTable>(1, "LHENamedWeights", true);
357  if (!hasIssuedWarning_.exchange(true)) {
358  edm::LogWarning("LHETablesProducer") << "No LHEEventProduct, so there will be no LHE Tables\n";
359  }
360  }
361 
362  iEvent.put(std::move(lheScaleTab), "LHEScale");
363  iEvent.put(std::move(lhePdfTab), "LHEPdf");
364  iEvent.put(std::move(lheRwgtTab), "LHEReweighting");
365  iEvent.put(std::move(lheNamedTab), "LHENamed");
366  iEvent.put(std::move(genPSTab), "PS");
367  }
368 
370  const DynamicWeightChoice* weightChoice,
371  const DynamicWeightChoiceGenInfo* genWeightChoice,
372  double genWeight,
373  const LHEEventProduct& lheProd,
374  const GenEventInfoProduct& genProd,
375  std::unique_ptr<nanoaod::FlatTable>& outScale,
376  std::unique_ptr<nanoaod::FlatTable>& outPdf,
377  std::unique_ptr<nanoaod::FlatTable>& outRwgt,
378  std::unique_ptr<nanoaod::FlatTable>& outNamed,
379  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
380  bool lheDebug = debug_.exchange(
381  false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
382 
383  const std::vector<std::string>& scaleWeightIDs = weightChoice->scaleWeightIDs;
384  const std::vector<std::string>& pdfWeightIDs = weightChoice->pdfWeightIDs;
385  const std::vector<std::string>& rwgtWeightIDs = weightChoice->rwgtIDs;
386 
387  double w0 = lheProd.originalXWGTUP();
388 
389  std::vector<double> wScale(scaleWeightIDs.size(), 1), wPDF(pdfWeightIDs.size(), 1), wRwgt(rwgtWeightIDs.size(), 1),
390  wNamed(namedWeightIDs_.size(), 1);
391  for (auto& weight : lheProd.weights()) {
392  if (lheDebug)
393  printf("Weight %+9.5f rel %+9.5f for id %s\n", weight.wgt, weight.wgt / w0, weight.id.c_str());
394  // now we do it slowly, can be optimized
395  auto mScale = std::find(scaleWeightIDs.begin(), scaleWeightIDs.end(), weight.id);
396  if (mScale != scaleWeightIDs.end())
397  wScale[mScale - scaleWeightIDs.begin()] = weight.wgt / w0;
398 
399  auto mPDF = std::find(pdfWeightIDs.begin(), pdfWeightIDs.end(), weight.id);
400  if (mPDF != pdfWeightIDs.end())
401  wPDF[mPDF - pdfWeightIDs.begin()] = weight.wgt / w0;
402 
403  auto mRwgt = std::find(rwgtWeightIDs.begin(), rwgtWeightIDs.end(), weight.id);
404  if (mRwgt != rwgtWeightIDs.end())
405  wRwgt[mRwgt - rwgtWeightIDs.begin()] = weight.wgt / w0;
406 
407  auto mNamed = std::find(namedWeightIDs_.begin(), namedWeightIDs_.end(), weight.id);
408  if (mNamed != namedWeightIDs_.end())
409  wNamed[mNamed - namedWeightIDs_.begin()] = weight.wgt / w0;
410  }
411 
412  std::vector<double> wPS;
413  std::string psWeightDocStr;
414  setPSWeightInfo(genProd.weights(), genWeightChoice, wPS, psWeightDocStr);
415 
416  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
417  outPS->addColumn<float>("", wPS, psWeightDocStr, lheWeightPrecision_);
418 
419  outScale = std::make_unique<nanoaod::FlatTable>(wScale.size(), "LHEScaleWeight", false);
420  outScale->addColumn<float>("", wScale, weightChoice->scaleWeightsDoc, lheWeightPrecision_);
421 
422  outPdf = std::make_unique<nanoaod::FlatTable>(wPDF.size(), "LHEPdfWeight", false);
423  outPdf->addColumn<float>("", wPDF, weightChoice->pdfWeightsDoc, lheWeightPrecision_);
424 
425  outRwgt = std::make_unique<nanoaod::FlatTable>(wRwgt.size(), "LHEReweightingWeight", false);
426  outRwgt->addColumn<float>("", wRwgt, weightChoice->rwgtWeightDoc, lheWeightPrecision_);
427 
428  outNamed = std::make_unique<nanoaod::FlatTable>(1, "LHEWeight", true);
429  outNamed->addColumnValue<float>("originalXWGTUP", lheProd.originalXWGTUP(), "Nominal event weight in the LHE file");
430  for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) {
431  outNamed->addColumnValue<float>(namedWeightLabels_[i],
432  wNamed[i],
433  "LHE weight for id " + namedWeightIDs_[i] + ", relative to nominal",
435  }
436 
437  counter->incLHE(genWeight, wScale, wPDF, wRwgt, wNamed, wPS);
438  }
439 
441  const DynamicWeightChoiceGenInfo* weightChoice,
442  double genWeight,
443  const GenEventInfoProduct& genProd,
444  std::unique_ptr<nanoaod::FlatTable>& outScale,
445  std::unique_ptr<nanoaod::FlatTable>& outPdf,
446  std::unique_ptr<nanoaod::FlatTable>& outNamed,
447  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
448  const std::vector<unsigned int>& scaleWeightIDs = weightChoice->scaleWeightIDs;
449  const std::vector<unsigned int>& pdfWeightIDs = weightChoice->pdfWeightIDs;
450 
451  auto weights = genProd.weights();
452  double w0 = (weights.size() > 1) ? weights.at(1) : 1.;
453  double originalXWGTUP = (weights.size() > 1) ? weights.at(1) : 1.;
454 
455  std::vector<double> wScale, wPDF, wPS;
456  for (auto id : scaleWeightIDs)
457  wScale.push_back(weights.at(id) / w0);
458  for (auto id : pdfWeightIDs) {
459  wPDF.push_back(weights.at(id) / w0);
460  }
461 
462  std::string psWeightsDocStr;
463  setPSWeightInfo(genProd.weights(), weightChoice, wPS, psWeightsDocStr);
464 
465  outScale = std::make_unique<nanoaod::FlatTable>(wScale.size(), "LHEScaleWeight", false);
466  outScale->addColumn<float>("", wScale, weightChoice->scaleWeightsDoc, lheWeightPrecision_);
467 
468  outPdf = std::make_unique<nanoaod::FlatTable>(wPDF.size(), "LHEPdfWeight", false);
469  outPdf->addColumn<float>("", wPDF, weightChoice->pdfWeightsDoc, lheWeightPrecision_);
470 
471  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
472  outPS->addColumn<float>("", wPS, psWeightsDocStr, lheWeightPrecision_);
473 
474  outNamed = std::make_unique<nanoaod::FlatTable>(1, "LHEWeight", true);
475  outNamed->addColumnValue<float>("originalXWGTUP", originalXWGTUP, "Nominal event weight in the LHE file");
476  /*for (unsigned int i = 0, n = wNamed.size(); i < n; ++i) {
477  outNamed->addColumnValue<float>(namedWeightLabels_[i], wNamed[i], "LHE weight for id "+namedWeightIDs_[i]+", relative to nominal", lheWeightPrecision_);
478  }*/
479 
480  counter->incLHE(genWeight, wScale, wPDF, std::vector<double>(), std::vector<double>(), wPS);
481  }
482 
484  const DynamicWeightChoiceGenInfo* genWeightChoice,
485  double genWeight,
486  const GenEventInfoProduct& genProd,
487  std::unique_ptr<nanoaod::FlatTable>& outPS) const {
488  std::vector<double> wPS;
489  std::string psWeightDocStr;
490  setPSWeightInfo(genProd.weights(), genWeightChoice, wPS, psWeightDocStr);
491  outPS = std::make_unique<nanoaod::FlatTable>(wPS.size(), "PSWeight", false);
492  outPS->addColumn<float>("", wPS, psWeightDocStr, lheWeightPrecision_);
493 
494  counter->incGenOnly(genWeight);
495  counter->incPSOnly(genWeight, wPS);
496  }
497 
498  void setPSWeightInfo(const std::vector<double>& genWeights,
499  const DynamicWeightChoiceGenInfo* genWeightChoice,
500  std::vector<double>& wPS,
501  std::string& psWeightDocStr) const {
502  wPS.clear();
503  // isRegularPSSet = keeping all weights and the weights are a usual size, ie
504  // all weights are PS weights (don't use header incase missing names)
505  bool isRegularPSSet = keepAllPSWeights_ && (genWeights.size() == 14 || genWeights.size() == 46);
506  if (!genWeightChoice->psWeightIDs.empty() && !isRegularPSSet) {
507  psWeightDocStr = genWeightChoice->psWeightsDoc;
508  double psNom = genWeights.at(genWeightChoice->psBaselineID);
509  for (auto wgtidx : genWeightChoice->psWeightIDs) {
510  wPS.push_back(genWeights.at(wgtidx) / psNom);
511  }
512  } else {
513  int vectorSize =
514  keepAllPSWeights_ ? (genWeights.size() - 2) : ((genWeights.size() == 14 || genWeights.size() == 46) ? 4 : 1);
515 
516  if (vectorSize > 1) {
517  double nominal = genWeights.at(1); // Called 'Baseline' in GenLumiInfoHeader
518  if (keepAllPSWeights_) {
519  for (int i = 0; i < vectorSize; i++) {
520  wPS.push_back(genWeights.at(i + 2) / nominal);
521  }
522  psWeightDocStr = "All PS weights (w_var / w_nominal)";
523  } else {
524  if (!psWeightWarning_.exchange(true))
525  edm::LogWarning("LHETablesProducer")
526  << "GenLumiInfoHeader not found: Central PartonShower weights will fill with the 6-10th entries \n"
527  << " This may incorrect for some mcs (madgraph 2.6.1 with its `isr:murfact=0.5` have a differnt "
528  "order )";
529  for (std::size_t i = 6; i < 10; i++) {
530  wPS.push_back(genWeights.at(i) / nominal);
531  }
532  psWeightDocStr =
533  "PS weights (w_var / w_nominal); [0] is ISR=2 FSR=1; [1] is ISR=1 FSR=2"
534  "[2] is ISR=0.5 FSR=1; [3] is ISR=1 FSR=0.5;";
535  }
536  } else {
537  wPS.push_back(1.0);
538  psWeightDocStr = "dummy PS weight (1.0) ";
539  }
540  }
541  }
542 
543  // create an empty counter
544  std::shared_ptr<DynamicWeightChoice> globalBeginRun(edm::Run const& iRun, edm::EventSetup const&) const override {
546 
547  bool lheDebug = debugRun_.exchange(
548  false); // make sure only the first thread dumps out this (even if may still be mixed up with other output, but nevermind)
549  auto weightChoice = std::make_shared<DynamicWeightChoice>();
550 
551  // getByToken throws since we're not in the endRun (see https://github.com/cms-sw/cmssw/pull/18499)
552  //if (iRun.getByToken(lheRunTag_, lheInfo)) {
553  for (const auto& lheLabel : lheLabel_) {
554  iRun.getByLabel(lheLabel, lheInfo);
555  if (lheInfo.isValid()) {
556  break;
557  }
558  }
559  if (lheInfo.isValid()) {
560  std::vector<ScaleVarWeight> scaleVariationIDs;
561  std::vector<PDFSetWeights> pdfSetWeightIDs;
562  std::vector<std::string> lheReweighingIDs;
563  bool isFirstGroup = true;
564 
565  std::regex weightgroupmg26x("<weightgroup\\s+(?:name|type)=\"(.*)\"\\s+combine=\"(.*)\"\\s*>");
566  std::regex weightgroup("<weightgroup\\s+combine=\"(.*)\"\\s+(?:name|type)=\"(.*)\"\\s*>");
567  std::regex weightgroupRwgt("<weightgroup\\s+(?:name|type)=\"(.*)\"\\s*>");
568  std::regex endweightgroup("</weightgroup>");
569  std::regex scalewmg26x(
570  "<weight\\s+(?:.*\\s+)?id=\"(\\d+)\"\\s*(?:lhapdf=\\d+|dyn=\\s*-?\\d+)?\\s*((?:[mM][uU][rR]|renscfact)=\"("
571  "\\S+)\"\\s+(?:[mM][uU][Ff]|facscfact)=\"(\\S+)\")(\\s+.*)?</weight>");
572  std::regex scalewmg26xNew(
573  "<weight\\s*((?:[mM][uU][fF]|facscfact)=\"(\\S+)\"\\s+(?:[mM][uU][Rr]|renscfact)=\"(\\S+)\").+id=\"(\\d+)\"(."
574  "*)?</weight>");
575 
576  //<weight MUF="1.0" MUR="2.0" PDF="306000" id="1006"> MUR=2.0 </weight>
577  std::regex scalew(
578  "<weight\\s+(?:.*\\s+)?id=\"(\\d+|\\d+-NNLOPS)\">\\s*(?:lhapdf=\\d+|dyn=\\s*-?\\d+)?\\s*((?:mu[rR]|renscfact)"
579  "=(\\S+)\\s+(?:mu[Ff]|facscfact)=(\\S+)(\\s+.*)?)</weight>");
580  std::regex pdfw(
581  "<weight\\s+id=\"(\\d+)\">\\s*(?:PDF set|lhapdf|PDF|pdfset)\\s*=\\s*(\\d+)\\s*(?:\\s.*)?</weight>");
582  std::regex pdfwOld("<weight\\s+(?:.*\\s+)?id=\"(\\d+)\">\\s*Member \\s*(\\d+)\\s*(?:.*)</weight>");
583  std::regex pdfwmg26x(
584  "<weight\\s+id=\"(\\d+)\"\\s*MUR=\"(?:\\S+)\"\\s*MUF=\"(?:\\S+)\"\\s*(?:PDF "
585  "set|lhapdf|PDF|pdfset)\\s*=\\s*\"(\\d+)\"\\s*>\\s*(?:PDF=(\\d+)\\s*MemberID=(\\d+))?\\s*(?:\\s.*)?</"
586  "weight>");
587  //<weightgroup combine="symmhessian+as" name="NNPDF31_nnlo_as_0118_mc_hessian_pdfas">
588 
589  //<weight MUF="1.0" MUR="1.0" PDF="325300" id="1048"> PDF=325300 MemberID=0 </weight>
590  std::regex pdfwmg26xNew(
591  "<weight\\s+MUF=\"(?:\\S+)\"\\s*MUR=\"(?:\\S+)\"\\s*PDF=\"(?:\\S+)\"\\s*id=\"(\\S+)\"\\s*>"
592  "\\s*(?:PDF=(\\d+)\\s*MemberID=(\\d+))?\\s*(?:\\s.*)?</"
593  "weight>");
594 
595  std::regex rwgt("<weight\\s+id=\"(.+)\">(.+)?(</weight>)?");
596  std::smatch groups;
597  for (auto iter = lheInfo->headers_begin(), end = lheInfo->headers_end(); iter != end; ++iter) {
598  if (iter->tag() != "initrwgt") {
599  if (lheDebug)
600  std::cout << "Skipping LHE header with tag" << iter->tag() << std::endl;
601  continue;
602  }
603  if (lheDebug)
604  std::cout << "Found LHE header with tag" << iter->tag() << std::endl;
605  std::vector<std::string> lines = iter->lines();
606  bool missed_weightgroup =
607  false; //Needed because in some of the samples ( produced with MG26X ) a small part of the header info is ordered incorrectly
608  bool ismg26x = false;
609  bool ismg26xNew = false;
610  for (unsigned int iLine = 0, nLines = lines.size(); iLine < nLines;
611  ++iLine) { //First start looping through the lines to see which weightgroup pattern is matched
612  boost::replace_all(lines[iLine], "&lt;", "<");
613  boost::replace_all(lines[iLine], "&gt;", ">");
614  if (std::regex_search(lines[iLine], groups, weightgroupmg26x)) {
615  ismg26x = true;
616  } else if (std::regex_search(lines[iLine], groups, scalewmg26xNew) ||
617  std::regex_search(lines[iLine], groups, pdfwmg26xNew)) {
618  ismg26xNew = true;
619  }
620  }
621  for (unsigned int iLine = 0, nLines = lines.size(); iLine < nLines; ++iLine) {
622  if (lheDebug)
623  std::cout << lines[iLine];
624  if (std::regex_search(lines[iLine], groups, ismg26x ? weightgroupmg26x : weightgroup)) {
625  std::string groupname = groups.str(2);
626  if (ismg26x)
627  groupname = groups.str(1);
628  if (lheDebug)
629  std::cout << ">>> Looks like the beginning of a weight group for '" << groupname << "'" << std::endl;
630  if (groupname.find("scale_variation") == 0 || groupname == "Central scale variation" || isFirstGroup) {
631  if (lheDebug && groupname.find("scale_variation") != 0 && groupname != "Central scale variation")
632  std::cout << ">>> First weight is not scale variation, but assuming is the Central Weight" << std::endl;
633  else if (lheDebug)
634  std::cout << ">>> Looks like scale variation for theory uncertainties" << std::endl;
635  isFirstGroup = false;
636  for (++iLine; iLine < nLines; ++iLine) {
637  if (lheDebug) {
638  std::cout << " " << lines[iLine];
639  }
640  if (std::regex_search(
641  lines[iLine], groups, ismg26x ? scalewmg26x : (ismg26xNew ? scalewmg26xNew : scalew))) {
642  if (lheDebug)
643  std::cout << " >>> Scale weight " << groups[1].str() << " for " << groups[3].str() << " , "
644  << groups[4].str() << " , " << groups[5].str() << std::endl;
645  if (ismg26xNew) {
646  scaleVariationIDs.emplace_back(groups.str(4), groups.str(1), groups.str(3), groups.str(2));
647  } else {
648  scaleVariationIDs.emplace_back(groups.str(1), groups.str(2), groups.str(3), groups.str(4));
649  }
650  } else if (std::regex_search(lines[iLine], endweightgroup)) {
651  if (lheDebug)
652  std::cout << ">>> Looks like the end of a weight group" << std::endl;
653  if (!missed_weightgroup) {
654  break;
655  } else
656  missed_weightgroup = false;
657  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
658  if (lheDebug)
659  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
660  "of the group."
661  << std::endl;
662  if (ismg26x || ismg26xNew)
663  missed_weightgroup = true;
664  --iLine; // rewind by one, and go back to the outer loop
665  break;
666  }
667  }
668  } else if (groupname == "PDF_variation" || groupname.find("PDF_variation ") == 0) {
669  if (lheDebug)
670  std::cout << ">>> Looks like a new-style block of PDF weights for one or more pdfs" << std::endl;
671  for (++iLine; iLine < nLines; ++iLine) {
672  if (lheDebug)
673  std::cout << " " << lines[iLine];
674  if (std::regex_search(lines[iLine], groups, pdfw)) {
675  unsigned int lhaID = std::stoi(groups.str(2));
676  if (lheDebug)
677  std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID
678  << std::endl;
679  if (pdfSetWeightIDs.empty() || !pdfSetWeightIDs.back().maybe_add(groups.str(1), lhaID)) {
680  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
681  }
682  } else if (std::regex_search(lines[iLine], endweightgroup)) {
683  if (lheDebug)
684  std::cout << ">>> Looks like the end of a weight group" << std::endl;
685  if (!missed_weightgroup) {
686  break;
687  } else
688  missed_weightgroup = false;
689  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
690  if (lheDebug)
691  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
692  "of the group."
693  << std::endl;
694  if (ismg26x || ismg26xNew)
695  missed_weightgroup = true;
696  --iLine; // rewind by one, and go back to the outer loop
697  break;
698  }
699  }
700  } else if (groupname == "PDF_variation1" || groupname == "PDF_variation2") {
701  if (lheDebug)
702  std::cout << ">>> Looks like a new-style block of PDF weights for multiple pdfs" << std::endl;
703  unsigned int lastid = 0;
704  for (++iLine; iLine < nLines; ++iLine) {
705  if (lheDebug)
706  std::cout << " " << lines[iLine];
707  if (std::regex_search(lines[iLine], groups, pdfw)) {
708  unsigned int id = std::stoi(groups.str(1));
709  unsigned int lhaID = std::stoi(groups.str(2));
710  if (lheDebug)
711  std::cout << " >>> PDF weight " << groups.str(1) << " for " << groups.str(2) << " = " << lhaID
712  << std::endl;
713  if (id != (lastid + 1) || pdfSetWeightIDs.empty()) {
714  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
715  } else {
716  pdfSetWeightIDs.back().add(groups.str(1), lhaID);
717  }
718  lastid = id;
719  } else if (std::regex_search(lines[iLine], endweightgroup)) {
720  if (lheDebug)
721  std::cout << ">>> Looks like the end of a weight group" << std::endl;
722  if (!missed_weightgroup) {
723  break;
724  } else
725  missed_weightgroup = false;
726  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
727  if (lheDebug)
728  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
729  "of the group."
730  << std::endl;
731  if (ismg26x || ismg26xNew)
732  missed_weightgroup = true;
733  --iLine; // rewind by one, and go back to the outer loop
734  break;
735  }
736  }
737  } else if (lhaNameToID_.find(groupname) != lhaNameToID_.end()) {
738  if (lheDebug)
739  std::cout << ">>> Looks like an old-style PDF weight for an individual pdf" << std::endl;
740  unsigned int firstLhaID = lhaNameToID_.find(groupname)->second;
741  bool first = true;
742  for (++iLine; iLine < nLines; ++iLine) {
743  if (lheDebug)
744  std::cout << " " << lines[iLine];
745  if (std::regex_search(
746  lines[iLine], groups, ismg26x ? pdfwmg26x : (ismg26xNew ? pdfwmg26xNew : pdfwOld))) {
747  unsigned int member = 0;
748  if (!ismg26x && !ismg26xNew) {
749  member = std::stoi(groups.str(2));
750  } else if (ismg26xNew) {
751  if (!groups.str(3).empty()) {
752  member = std::stoi(groups.str(3));
753  }
754  } else {
755  if (!groups.str(4).empty()) {
756  member = std::stoi(groups.str(4));
757  }
758  }
759  unsigned int lhaID = member + firstLhaID;
760  if (lheDebug)
761  std::cout << " >>> PDF weight " << groups.str(1) << " for " << member << " = " << lhaID
762  << std::endl;
763  //if (member == 0) continue; // let's keep also the central value for now
764  if (first) {
765  pdfSetWeightIDs.emplace_back(groups.str(1), lhaID);
766  first = false;
767  } else {
768  pdfSetWeightIDs.back().add(groups.str(1), lhaID);
769  }
770  } else if (std::regex_search(lines[iLine], endweightgroup)) {
771  if (lheDebug)
772  std::cout << ">>> Looks like the end of a weight group" << std::endl;
773  if (!missed_weightgroup) {
774  break;
775  } else
776  missed_weightgroup = false;
777  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
778  if (lheDebug)
779  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
780  "of the group."
781  << std::endl;
782  if (ismg26x || ismg26xNew)
783  missed_weightgroup = true;
784  --iLine; // rewind by one, and go back to the outer loop
785  break;
786  }
787  }
788  } else if (groupname == "mass_variation" || groupname == "sthw2_variation" ||
789  groupname == "width_variation") {
790  if (lheDebug)
791  std::cout << ">>> Looks like an EW parameter weight" << std::endl;
792  for (++iLine; iLine < nLines; ++iLine) {
793  if (lheDebug)
794  std::cout << " " << lines[iLine];
795  if (std::regex_search(lines[iLine], groups, rwgt)) {
796  std::string rwgtID = groups.str(1);
797  if (lheDebug)
798  std::cout << " >>> LHE reweighting weight: " << rwgtID << std::endl;
799  if (std::find(lheReweighingIDs.begin(), lheReweighingIDs.end(), rwgtID) == lheReweighingIDs.end()) {
800  // we're only interested in the beggining of the block
801  lheReweighingIDs.emplace_back(rwgtID);
802  }
803  } else if (std::regex_search(lines[iLine], endweightgroup)) {
804  if (lheDebug)
805  std::cout << ">>> Looks like the end of a weight group" << std::endl;
806  }
807  }
808  } else {
809  for (++iLine; iLine < nLines; ++iLine) {
810  if (lheDebug)
811  std::cout << " " << lines[iLine];
812  if (std::regex_search(lines[iLine], groups, endweightgroup)) {
813  if (lheDebug)
814  std::cout << ">>> Looks like the end of a weight group" << std::endl;
815  if (!missed_weightgroup) {
816  break;
817  } else
818  missed_weightgroup = false;
819  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
820  if (lheDebug)
821  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
822  "of the group."
823  << std::endl;
824  if (ismg26x || ismg26xNew)
825  missed_weightgroup = true;
826  --iLine; // rewind by one, and go back to the outer loop
827  break;
828  }
829  }
830  }
831  } else if (std::regex_search(lines[iLine], groups, weightgroupRwgt)) {
832  std::string groupname = groups.str(1);
833  if (groupname.find("mg_reweighting") != std::string::npos) {
834  if (lheDebug)
835  std::cout << ">>> Looks like a LHE weights for reweighting" << std::endl;
836  for (++iLine; iLine < nLines; ++iLine) {
837  if (lheDebug)
838  std::cout << " " << lines[iLine];
839  if (std::regex_search(lines[iLine], groups, rwgt)) {
840  std::string rwgtID = groups.str(1);
841  if (lheDebug)
842  std::cout << " >>> LHE reweighting weight: " << rwgtID << std::endl;
843  if (std::find(lheReweighingIDs.begin(), lheReweighingIDs.end(), rwgtID) == lheReweighingIDs.end()) {
844  // we're only interested in the beggining of the block
845  lheReweighingIDs.emplace_back(rwgtID);
846  }
847  } else if (std::regex_search(lines[iLine], endweightgroup)) {
848  if (lheDebug)
849  std::cout << ">>> Looks like the end of a weight group" << std::endl;
850  if (!missed_weightgroup) {
851  break;
852  } else
853  missed_weightgroup = false;
854  } else if (std::regex_search(lines[iLine], ismg26x ? weightgroupmg26x : weightgroup)) {
855  if (lheDebug)
856  std::cout << ">>> Looks like the beginning of a new weight group, I will assume I missed the end "
857  "of the group."
858  << std::endl;
859  if (ismg26x)
860  missed_weightgroup = true;
861  --iLine; // rewind by one, and go back to the outer loop
862  break;
863  }
864  }
865  }
866  }
867  }
868  //std::cout << "============= END [ " << iter->tag() << " ] ============ \n\n" << std::endl;
869 
870  // ----- SCALE VARIATIONS -----
871  std::sort(scaleVariationIDs.begin(), scaleVariationIDs.end());
872  if (lheDebug)
873  std::cout << "Found " << scaleVariationIDs.size() << " scale variations: " << std::endl;
874  std::stringstream scaleDoc;
875  scaleDoc << "LHE scale variation weights (w_var / w_nominal); ";
876  for (unsigned int isw = 0, nsw = scaleVariationIDs.size(); isw < nsw; ++isw) {
877  const auto& sw = scaleVariationIDs[isw];
878  if (isw)
879  scaleDoc << "; ";
880  scaleDoc << "[" << isw << "] is " << sw.label;
881  weightChoice->scaleWeightIDs.push_back(sw.wid);
882  if (lheDebug)
883  printf(" id %s: scales ren = % .2f fact = % .2f text = %s\n",
884  sw.wid.c_str(),
885  sw.scales.first,
886  sw.scales.second,
887  sw.label.c_str());
888  }
889  if (!scaleVariationIDs.empty())
890  weightChoice->scaleWeightsDoc = scaleDoc.str();
891 
892  // ------ PDF VARIATIONS (take the preferred one) -----
893  if (lheDebug) {
894  std::cout << "Found " << pdfSetWeightIDs.size() << " PDF set errors: " << std::endl;
895  for (const auto& pw : pdfSetWeightIDs)
896  printf("lhaIDs %6d - %6d (%3lu weights: %s, ... )\n",
897  pw.lhaIDs.first,
898  pw.lhaIDs.second,
899  pw.wids.size(),
900  pw.wids.front().c_str());
901  }
902 
903  // ------ LHE REWEIGHTING -------
904  if (lheDebug) {
905  std::cout << "Found " << lheReweighingIDs.size() << " reweighting weights" << std::endl;
906  }
907  std::copy(lheReweighingIDs.begin(), lheReweighingIDs.end(), std::back_inserter(weightChoice->rwgtIDs));
908 
909  std::stringstream pdfDoc;
910  pdfDoc << "LHE pdf variation weights (w_var / w_nominal) for LHA IDs ";
911  bool found = false;
912  for (const auto& pw : pdfSetWeightIDs) {
913  for (uint32_t lhaid : preferredPDFLHAIDs_) {
914  if (pw.lhaIDs.first != lhaid && pw.lhaIDs.first != (lhaid + 1))
915  continue; // sometimes the first weight is not saved if that PDF is the nominal one for the sample
916  if (pw.wids.size() == 1)
917  continue; // only consider error sets
918  pdfDoc << pw.lhaIDs.first << " - " << pw.lhaIDs.second;
919  weightChoice->pdfWeightIDs = pw.wids;
920  if (maxPdfWeights_ < pw.wids.size()) {
921  weightChoice->pdfWeightIDs.resize(maxPdfWeights_); // drop some replicas
922  pdfDoc << ", truncated to the first " << maxPdfWeights_ << " replicas";
923  }
924  weightChoice->pdfWeightsDoc = pdfDoc.str();
925  found = true;
926  break;
927  }
928  if (found)
929  break;
930  }
931  }
932  }
933  return weightChoice;
934  }
935 
936  // create an empty counter
937  std::unique_ptr<LumiCacheInfoHolder> beginStream(edm::StreamID) const override {
938  return std::make_unique<LumiCacheInfoHolder>();
939  }
940  // inizialize to zero at begin run
941  void streamBeginRun(edm::StreamID id, edm::Run const&, edm::EventSetup const&) const override {
942  streamCache(id)->clear();
943  }
945  edm::LuminosityBlock const& lumiBlock,
946  edm::EventSetup const& eventSetup) const override {
947  auto counterMap = &(streamCache(id)->countermap);
948  edm::Handle<GenLumiInfoHeader> genLumiInfoHead;
949  lumiBlock.getByToken(genLumiInfoHeadTag_, genLumiInfoHead);
950  if (!genLumiInfoHead.isValid())
951  edm::LogWarning("LHETablesProducer")
952  << "No GenLumiInfoHeader product found, will not fill generator model string.\n";
953 
955  if (genLumiInfoHead.isValid()) {
956  label = genLumiInfoHead->configDescription();
957  boost::replace_all(label, "-", "_");
958  boost::replace_all(label, "/", "_");
959  }
960  counterMap->setLabel(label);
961 
962  if (genLumiInfoHead.isValid()) {
963  auto weightChoice = &(streamCache(id)->weightChoice);
964 
965  std::vector<ScaleVarWeight> scaleVariationIDs;
966  std::vector<PDFSetWeights> pdfSetWeightIDs;
967  weightChoice->psWeightIDs.clear();
968 
969  std::regex scalew("LHE,\\s+id\\s+=\\s+(\\d+),\\s+(.+)\\,\\s+mur=(\\S+)\\smuf=(\\S+)");
970  std::regex pdfw("LHE,\\s+id\\s+=\\s+(\\d+),\\s+(.+),\\s+Member\\s+(\\d+)\\s+of\\ssets\\s+(\\w+\\b)");
971  std::regex mainPSw("sr(Def|:murfac=)(Hi|Lo|_dn|_up|0.5|2.0)");
972  std::smatch groups;
973  auto weightNames = genLumiInfoHead->weightNames();
974  std::unordered_map<std::string, uint32_t> knownPDFSetsFromGenInfo_;
975  unsigned int weightIter = 0;
976  for (const auto& line : weightNames) {
977  if (std::regex_search(line, groups, scalew)) { // scale variation
978  auto id = groups.str(1);
979  auto group = groups.str(2);
980  auto mur = groups.str(3);
981  auto muf = groups.str(4);
982  if (group.find("Central scale variation") != std::string::npos)
983  scaleVariationIDs.emplace_back(groups.str(1), groups.str(2), groups.str(3), groups.str(4));
984  } else if (std::regex_search(line, groups, pdfw)) { // PDF variation
985  auto id = groups.str(1);
986  auto group = groups.str(2);
987  auto memberid = groups.str(3);
988  auto pdfset = groups.str(4);
989  if (group.find(pdfset) != std::string::npos) {
990  if (knownPDFSetsFromGenInfo_.find(pdfset) == knownPDFSetsFromGenInfo_.end()) {
991  knownPDFSetsFromGenInfo_[pdfset] = std::atoi(id.c_str());
992  pdfSetWeightIDs.emplace_back(id, std::atoi(id.c_str()));
993  } else
994  pdfSetWeightIDs.back().add(id, std::atoi(id.c_str()));
995  }
996  } else if (line == "Baseline") {
997  weightChoice->psBaselineID = weightIter;
998  } else if (line.find("isr") != std::string::npos || line.find("fsr") != std::string::npos) {
999  weightChoice->matchPS_alt = line.find("sr:") != std::string::npos; // (f/i)sr: for new weights
1000  if (keepAllPSWeights_) {
1001  weightChoice->psWeightIDs.push_back(weightIter); // PS variations
1002  } else if (std::regex_search(line, groups, mainPSw)) {
1003  if (weightChoice->psWeightIDs.empty())
1004  weightChoice->psWeightIDs = std::vector<unsigned int>(4, -1);
1005  int psIdx = (line.find("fsr") != std::string::npos) ? 1 : 0;
1006  psIdx += (groups.str(2) == "Hi" || groups.str(2) == "_up" || groups.str(2) == "2.0") ? 0 : 2;
1007  weightChoice->psWeightIDs[psIdx] = weightIter;
1008  }
1009  }
1010  weightIter++;
1011  }
1012  if (keepAllPSWeights_) {
1013  weightChoice->psWeightsDoc = "All PS weights (w_var / w_nominal)";
1014  } else if (weightChoice->psWeightIDs.size() == 4) {
1015  weightChoice->psWeightsDoc =
1016  "PS weights (w_var / w_nominal); [0] is ISR=2 FSR=1; [1] is ISR=1 FSR=2"
1017  "[2] is ISR=0.5 FSR=1; [3] is ISR=1 FSR=0.5;";
1018  for (int i = 0; i < 4; i++) {
1019  if (static_cast<int>(weightChoice->psWeightIDs[i]) == -1)
1020  weightChoice->setMissingWeight(i);
1021  }
1022  } else {
1023  weightChoice->psWeightsDoc = "dummy PS weight (1.0) ";
1024  }
1025 
1026  weightChoice->scaleWeightIDs.clear();
1027  weightChoice->pdfWeightIDs.clear();
1028 
1029  std::sort(scaleVariationIDs.begin(), scaleVariationIDs.end());
1030  std::stringstream scaleDoc;
1031  scaleDoc << "LHE scale variation weights (w_var / w_nominal); ";
1032  for (unsigned int isw = 0, nsw = scaleVariationIDs.size(); isw < nsw; ++isw) {
1033  const auto& sw = scaleVariationIDs[isw];
1034  if (isw)
1035  scaleDoc << "; ";
1036  scaleDoc << "[" << isw << "] is " << sw.label;
1037  weightChoice->scaleWeightIDs.push_back(std::atoi(sw.wid.c_str()));
1038  }
1039  if (!scaleVariationIDs.empty())
1040  weightChoice->scaleWeightsDoc = scaleDoc.str();
1041  std::stringstream pdfDoc;
1042  pdfDoc << "LHE pdf variation weights (w_var / w_nominal) for LHA names ";
1043  bool found = false;
1044  for (const auto& pw : pdfSetWeightIDs) {
1045  if (pw.wids.size() == 1)
1046  continue; // only consider error sets
1047  for (const auto& wantedpdf : lhaNameToID_) {
1048  auto pdfname = wantedpdf.first;
1049  if (knownPDFSetsFromGenInfo_.find(pdfname) == knownPDFSetsFromGenInfo_.end())
1050  continue;
1051  uint32_t lhaid = knownPDFSetsFromGenInfo_.at(pdfname);
1052  if (pw.lhaIDs.first != lhaid)
1053  continue;
1054  pdfDoc << pdfname;
1055  for (const auto& x : pw.wids)
1056  weightChoice->pdfWeightIDs.push_back(std::atoi(x.c_str()));
1057  if (maxPdfWeights_ < pw.wids.size()) {
1058  weightChoice->pdfWeightIDs.resize(maxPdfWeights_); // drop some replicas
1059  pdfDoc << ", truncated to the first " << maxPdfWeights_ << " replicas";
1060  }
1061  weightChoice->pdfWeightsDoc = pdfDoc.str();
1062  found = true;
1063  break;
1064  }
1065  if (found)
1066  break;
1067  }
1068  }
1069  }
1070  // create an empty counter
1071  std::shared_ptr<CounterMap> globalBeginRunSummary(edm::Run const&, edm::EventSetup const&) const override {
1072  return std::make_shared<CounterMap>();
1073  }
1074  // add this stream to the summary
1076  edm::Run const&,
1077  edm::EventSetup const&,
1078  CounterMap* runCounterMap) const override {
1079  runCounterMap->merge(streamCache(id)->countermap);
1080  }
1081  // nothing to do per se
1082  void globalEndRunSummary(edm::Run const&, edm::EventSetup const&, CounterMap* runCounterMap) const override {}
1083  // write the total to the run
1084  void globalEndRunProduce(edm::Run& iRun, edm::EventSetup const&, CounterMap const* runCounterMap) const override {
1085  auto out = std::make_unique<nanoaod::MergeableCounterTable>();
1086 
1087  for (const auto& x : runCounterMap->countermap) {
1088  auto runCounter = &(x.second);
1089  std::string label = (!x.first.empty()) ? (std::string("_") + x.first) : "";
1090  std::string doclabel = (!x.first.empty()) ? (std::string(", for model label ") + x.first) : "";
1091 
1092  out->addInt("genEventCount" + label, "event count" + doclabel, runCounter->num);
1093  out->addFloat("genEventSumw" + label, "sum of gen weights" + doclabel, runCounter->sumw);
1094  out->addFloat("genEventSumw2" + label, "sum of gen (weight^2)" + doclabel, runCounter->sumw2);
1095 
1096  double norm = runCounter->sumw ? 1.0 / runCounter->sumw : 1;
1097  auto sumScales = runCounter->sumScale;
1098  for (auto& val : sumScales)
1099  val *= norm;
1100  out->addVFloatWithNorm("LHEScaleSumw" + label,
1101  "Sum of genEventWeight * LHEScaleWeight[i], divided by genEventSumw" + doclabel,
1102  sumScales,
1103  runCounter->sumw);
1104  auto sumPDFs = runCounter->sumPDF;
1105  for (auto& val : sumPDFs)
1106  val *= norm;
1107  out->addVFloatWithNorm("LHEPdfSumw" + label,
1108  "Sum of genEventWeight * LHEPdfWeight[i], divided by genEventSumw" + doclabel,
1109  sumPDFs,
1110  runCounter->sumw);
1111  auto sumPS = runCounter->sumPS;
1112  for (auto& val : sumPS)
1113  val *= norm;
1114  out->addVFloatWithNorm("PSSumw" + label,
1115  "Sum of genEventWeight * PSWeight[i], divided by genEventSumw" + doclabel,
1116  sumPS,
1117  runCounter->sumw);
1118  if (!runCounter->sumRwgt.empty()) {
1119  auto sumRwgts = runCounter->sumRwgt;
1120  for (auto& val : sumRwgts)
1121  val *= norm;
1122  out->addVFloatWithNorm("LHEReweightingSumw" + label,
1123  "Sum of genEventWeight * LHEReweightingWeight[i], divided by genEventSumw" + doclabel,
1124  sumRwgts,
1125  runCounter->sumw);
1126  }
1127  if (!runCounter->sumNamed.empty()) { // it could be empty if there's no LHE info in the sample
1128  for (unsigned int i = 0, n = namedWeightLabels_.size(); i < n; ++i) {
1129  out->addFloatWithNorm(
1130  "LHESumw_" + namedWeightLabels_[i] + label,
1131  "Sum of genEventWeight * LHEWeight_" + namedWeightLabels_[i] + ", divided by genEventSumw" + doclabel,
1132  runCounter->sumNamed[i] * norm,
1133  runCounter->sumw);
1134  }
1135  }
1136  }
1137  iRun.put(std::move(out));
1138  }
1139  // nothing to do here
1140  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override {}
1141 
1144  desc.add<edm::InputTag>("genEvent", edm::InputTag("generator"))
1145  ->setComment("tag for the GenEventInfoProduct, to get the main weight");
1146  desc.add<edm::InputTag>("genLumiInfoHeader", edm::InputTag("generator"))
1147  ->setComment("tag for the GenLumiInfoProduct, to get the model string");
1148  desc.add<std::vector<edm::InputTag>>("lheInfo", std::vector<edm::InputTag>{{"externalLHEProducer"}, {"source"}})
1149  ->setComment("tag(s) for the LHE information (LHEEventProduct and LHERunInfoProduct)");
1150 
1152  prefpdf.add<std::string>("name");
1153  prefpdf.add<uint32_t>("lhaid");
1154  desc.addVPSet("preferredPDFs", prefpdf, std::vector<edm::ParameterSet>())
1155  ->setComment(
1156  "LHA PDF Ids of the preferred PDF sets, in order of preference (the first matching one will be used)");
1157  desc.add<std::vector<std::string>>("namedWeightIDs")->setComment("set of LHA weight IDs for named LHE weights");
1158  desc.add<std::vector<std::string>>("namedWeightLabels")
1159  ->setComment("output names for the namedWeightIDs (in the same order)");
1160  desc.add<int32_t>("lheWeightPrecision")->setComment("Number of bits in the mantissa for LHE weights");
1161  desc.add<uint32_t>("maxPdfWeights")->setComment("Maximum number of PDF weights to save (to crop NN replicas)");
1162  desc.add<bool>("keepAllPSWeights")->setComment("Store all PS weights found");
1163  desc.addOptionalUntracked<bool>("debug")->setComment("dump out all LHE information for one event");
1164  descriptions.add("genWeightsTable", desc);
1165  }
1166 
1167 protected:
1169  const std::vector<edm::InputTag> lheLabel_;
1170  const std::vector<edm::EDGetTokenT<LHEEventProduct>> lheTag_;
1171  const std::vector<edm::EDGetTokenT<LHERunInfoProduct>> lheRunTag_;
1173 
1174  std::vector<uint32_t> preferredPDFLHAIDs_;
1175  std::unordered_map<std::string, uint32_t> lhaNameToID_;
1176  std::vector<std::string> namedWeightIDs_;
1177  std::vector<std::string> namedWeightLabels_;
1179  unsigned int maxPdfWeights_;
1181 
1182  mutable std::atomic<bool> debug_, debugRun_, hasIssuedWarning_, psWeightWarning_;
1183 };
1184 
const std::vector< edm::EDGetTokenT< LHERunInfoProduct > > lheRunTag_
std::shared_ptr< CounterMap > globalBeginRunSummary(edm::Run const &, edm::EventSetup const &) const override
double originalXWGTUP() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T w() const
const edm::EDGetTokenT< GenEventInfoProduct > genTag_
void setPSWeightInfo(const std::vector< double > &genWeights, const DynamicWeightChoiceGenInfo *genWeightChoice, std::vector< double > &wPS, std::string &psWeightDocStr) const
void fillOnlyPSWeightTable(Counter *counter, const DynamicWeightChoiceGenInfo *genWeightChoice, double genWeight, const GenEventInfoProduct &genProd, std::unique_ptr< nanoaod::FlatTable > &outPS) const
Definition: weight.py:1
constexpr int pow(int x)
Definition: conifer.h:24
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:19
uint16_t size_type
const std::string & configDescription() const
const std::vector< std::string > & weightNames() const
const std::vector< edm::InputTag > lheLabel_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void fillLHEWeightTables(Counter *counter, const DynamicWeightChoice *weightChoice, const DynamicWeightChoiceGenInfo *genWeightChoice, 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
char const * label
std::function< unsigned int(align::ID)> Counter
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< LumiCacheInfoHolder > beginStream(edm::StreamID) const override
void globalEndRunProduce(edm::Run &iRun, edm::EventSetup const &, CounterMap const *runCounterMap) const override
int merge(int argc, char *argv[])
Definition: DMRmerge.cc:37
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:283
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::atomic< bool > psWeightWarning_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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
void globalEndRun(edm::Run const &, edm::EventSetup const &) const override
void streamEndRunSummary(edm::StreamID id, edm::Run const &, edm::EventSetup const &, CounterMap *runCounterMap) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
std::vector< std::string > namedWeightLabels_
GenWeightsTableProducer(edm::ParameterSet const &params)
std::shared_ptr< DynamicWeightChoice > globalBeginRun(edm::Run const &iRun, edm::EventSetup const &) const override
void put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Run.h:109
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
void globalEndRunSummary(edm::Run const &, edm::EventSetup const &, CounterMap *runCounterMap) const override
const std::vector< edm::EDGetTokenT< LHEEventProduct > > lheTag_
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
Definition: DTCELinkId.h:70
std::vector< uint32_t > preferredPDFLHAIDs_
void add(std::map< std::string, TH1 *> &h, TH1 *hist)
HLT enums.
std::vector< double > & weights()
void streamBeginLuminosityBlock(edm::StreamID id, edm::LuminosityBlock const &lumiBlock, edm::EventSetup const &eventSetup) const override
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 clear(EGIsoObj &c)
Definition: egamma.h:82
Log< level::Warning, false > LogWarning
void streamBeginRun(edm::StreamID id, edm::Run const &, edm::EventSetup const &) const override
#define str(s)
const edm::EDGetTokenT< GenLumiInfoHeader > genLumiInfoHeadTag_
def move(src, dest)
Definition: eostools.py:511
const std::vector< WGT > & weights() const
Definition: Run.h:45