CMS 3D CMS Logo

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