CMS 3D CMS Logo

GenXSecAnalyzer.cc
Go to the documentation of this file.
3 #include "TMath.h"
4 #include <iostream>
5 #include <iomanip>
6 
8  : nMCs_(0),
9  hepidwtup_(-9999),
10  totalWeightPre_(0),
11  thisRunWeightPre_(0),
12  totalWeight_(0),
13  thisRunWeight_(0),
14  xsecPreFilter_(-1, -1),
15  xsec_(-1, -1),
16  product_(GenLumiInfoProduct(-9999)),
17  filterOnlyEffRun_(0, 0, 0, 0, 0., 0., 0., 0.),
18  hepMCFilterEffRun_(0, 0, 0, 0, 0., 0., 0., 0.),
19  filterOnlyEffStat_(0, 0, 0, 0, 0., 0., 0., 0.),
20  hepMCFilterEffStat_(0, 0, 0, 0, 0., 0., 0., 0.) {
21  xsecBeforeMatching_.clear();
22  xsecAfterMatching_.clear();
23  jetMatchEffStat_.clear();
26 
27  genFilterInfoToken_ = consumes<GenFilterInfo, edm::InLumi>(edm::InputTag("genFilterEfficiencyProducer", ""));
28  hepMCFilterInfoToken_ = consumes<GenFilterInfo, edm::InLumi>(edm::InputTag("generator", ""));
29  genLumiInfoToken_ = consumes<GenLumiInfoProduct, edm::InLumi>(edm::InputTag("generator", ""));
30  lheRunInfoToken_ = consumes<LHERunInfoProduct, edm::InRun>(edm::InputTag("externalLHEProducer", ""));
31 }
32 
34 
36  xsecBeforeMatching_.clear();
37  xsecAfterMatching_.clear();
38  jetMatchEffStat_.clear();
41 }
42 
44  // initialization for every different physics MC
45 
46  nMCs_++;
47 
49  thisRunWeight_ = 0;
50 
52 
53  filterOnlyEffRun_ = GenFilterInfo(0, 0, 0, 0, 0., 0., 0., 0.);
54  hepMCFilterEffRun_ = GenFilterInfo(0, 0, 0, 0, 0., 0., 0., 0.);
55 
56  xsecBeforeMatching_.clear();
57  xsecAfterMatching_.clear();
58  jetMatchEffStat_.clear();
61 
62  return;
63 }
64 
66 
68 
71  iLumi.getByToken(genLumiInfoToken_, genLumiInfo);
72  if (!genLumiInfo.isValid())
73  return;
74  hepidwtup_ = genLumiInfo->getHEPIDWTUP();
75 
76  std::vector<GenLumiInfoProduct::ProcessInfo> theProcesses = genLumiInfo->getProcessInfos();
77 
78  unsigned int theProcesses_size = theProcesses.size();
79 
80  // if it's a pure parton-shower generator, check there should be only one element in thisProcessInfos
81  // the error of lheXSec is -1
82  if (hepidwtup_ == -1) {
83  if (theProcesses_size != 1) {
84  edm::LogError("GenXSecAnalyzer::endLuminosityBlock") << "Pure parton shower has thisProcessInfos size!=1";
85  return;
86  }
87  }
88 
89  for (unsigned int ip = 0; ip < theProcesses_size; ip++) {
90  if (theProcesses[ip].lheXSec().value() < 0) {
91  edm::LogError("GenXSecAnalyzer::endLuminosityBlock")
92  << "cross section of process " << ip << " value = " << theProcesses[ip].lheXSec().value();
93  return;
94  }
95  }
96 
97  product_.mergeProduct(*genLumiInfo);
98 
100  iLumi.getByToken(genFilterInfoToken_, genFilter);
101  if (genFilter.isValid()) {
102  filterOnlyEffStat_.mergeProduct(*genFilter);
103  filterOnlyEffRun_.mergeProduct(*genFilter);
104  thisRunWeight_ += genFilter->sumPassWeights();
105  }
106 
107  edm::Handle<GenFilterInfo> hepMCFilter;
108  iLumi.getByToken(hepMCFilterInfoToken_, hepMCFilter);
109 
110  if (hepMCFilter.isValid()) {
111  hepMCFilterEffStat_.mergeProduct(*hepMCFilter);
112  hepMCFilterEffRun_.mergeProduct(*hepMCFilter);
113  }
114 
115  // doing generic summing for jet matching statistics
116  // and computation of combined LHE information
117  for (unsigned int ip = 0; ip < theProcesses_size; ip++) {
118  int id = theProcesses[ip].process();
121  GenLumiInfoProduct::FinalStat temp_killed = theProcesses[ip].killed();
122  GenLumiInfoProduct::FinalStat temp_selected = theProcesses[ip].selected();
123  double passw = temp_killed.sum();
124  double passw2 = temp_killed.sum2();
125  double totalw = temp_selected.sum();
126  double totalw2 = temp_selected.sum2();
127  GenFilterInfo tempInfo(theProcesses[ip].nPassPos(),
128  theProcesses[ip].nPassNeg(),
129  theProcesses[ip].nTotalPos(),
130  theProcesses[ip].nTotalNeg(),
131  passw,
132  passw2,
133  totalw,
134  totalw2);
135 
136  // matching statistics for all processes
137  jetMatchEffStat_[10000].mergeProduct(tempInfo);
138  double currentValue = theProcesses[ip].lheXSec().value();
139  double currentError = theProcesses[ip].lheXSec().error();
140 
141  // this process ID has occurred before
142  if (y.value() > 0) {
143  x.mergeProduct(tempInfo);
144  double previousValue = previousLumiBlockLHEXSec_[id].value();
145 
146  if (currentValue != previousValue) // transition of cross section
147  {
148  double xsec = y.value();
149  double err = y.error();
150  combine(xsec, err, thisRunWeightPre_, currentValue, currentError, totalw);
151  y = GenLumiInfoProduct::XSec(xsec, err);
152  } else // LHE cross section is the same as previous lumiblock
153  thisRunWeightPre_ += totalw;
154 
155  }
156  // this process ID has never occurred before
157  else {
158  x = tempInfo;
159  y = theProcesses[ip].lheXSec();
160  thisRunWeightPre_ += totalw;
161  }
162 
163  previousLumiBlockLHEXSec_[id] = theProcesses[ip].lheXSec();
164  } // end
165 
166  return;
167 }
168 
170  //xsection before matching
172 
173  if (iRun.getByToken(lheRunInfoToken_, run)) {
174  const lhef::HEPRUP thisHeprup_ = run->heprup();
175 
176  for (unsigned int iSize = 0; iSize < thisHeprup_.XSECUP.size(); iSize++) {
177  std::cout << std::setw(14) << std::fixed << thisHeprup_.XSECUP[iSize] << std::setw(14) << std::fixed
178  << thisHeprup_.XERRUP[iSize] << std::setw(14) << std::fixed << thisHeprup_.XMAXUP[iSize]
179  << std::setw(14) << std::fixed << thisHeprup_.LPRUP[iSize] << std::endl;
180  }
181  std::cout << " " << std::endl;
182  }
183 
184  // compute cross section for this run first
185  // set the correct combined LHE+filter cross sections
186  unsigned int i = 0;
187  std::vector<GenLumiInfoProduct::ProcessInfo> newInfos;
188  for (std::map<int, GenLumiInfoProduct::XSec>::const_iterator iter = currentLumiBlockLHEXSec_.begin();
189  iter != currentLumiBlockLHEXSec_.end();
190  ++iter, i++) {
192  temp.setLheXSec(iter->second.value(), iter->second.error());
193  newInfos.push_back(temp);
194  }
195  product_.setProcessInfo(newInfos);
196 
197  const GenLumiInfoProduct::XSec thisRunXSecPre = compute(product_);
198  // xsection after matching before filters
200 
201  double thisHepFilterEff = 1;
202  double thisHepFilterErr = 0;
203 
204  if (hepMCFilterEffRun_.sumWeights2() > 0) {
205  thisHepFilterEff = hepMCFilterEffRun_.filterEfficiency(hepidwtup_);
207  if (thisHepFilterEff < 0) {
208  thisHepFilterEff = 1;
209  thisHepFilterErr = 0;
210  }
211  }
212 
213  double thisGenFilterEff = 1;
214  double thisGenFilterErr = 0;
215 
216  if (filterOnlyEffRun_.sumWeights2() > 0) {
217  thisGenFilterEff = filterOnlyEffRun_.filterEfficiency(hepidwtup_);
219  if (thisGenFilterEff < 0) {
220  thisGenFilterEff = 1;
221  thisGenFilterErr = 0;
222  }
223  }
224  double thisXsec = thisRunXSecPre.value() > 0 ? thisHepFilterEff * thisGenFilterEff * thisRunXSecPre.value() : 0;
225  double thisErr =
226  thisRunXSecPre.value() > 0
227  ? thisXsec * sqrt(pow(TMath::Max(thisRunXSecPre.error(), (double)0) / thisRunXSecPre.value(), 2) +
228  pow(thisHepFilterErr / thisHepFilterEff, 2) + pow(thisGenFilterErr / thisGenFilterEff, 2))
229  : 0;
230  const GenLumiInfoProduct::XSec thisRunXSec = GenLumiInfoProduct::XSec(thisXsec, thisErr);
231  combine(xsec_, totalWeight_, thisRunXSec, thisRunWeight_);
232 }
233 
234 void GenXSecAnalyzer::combine(double& finalValue,
235  double& finalError,
236  double& finalWeight,
237  const double& currentValue,
238  const double& currentError,
239  const double& currentWeight) {
240  if (finalValue <= 0) {
241  finalValue = currentValue;
242  finalError = currentError;
243  finalWeight += currentWeight;
244  } else {
245  double wgt1 = (finalError <= 0 || currentError <= 0) ? finalWeight : 1 / (finalError * finalError);
246  double wgt2 = (finalError <= 0 || currentError <= 0) ? currentWeight : 1 / (currentError * currentError);
247  double xsec = (wgt1 * finalValue + wgt2 * currentValue) / (wgt1 + wgt2);
248  double err = (finalError <= 0 || currentError <= 0) ? 0 : 1.0 / std::sqrt(wgt1 + wgt2);
249  finalValue = xsec;
250  finalError = err;
251  finalWeight += currentWeight;
252  }
253  return;
254 }
255 
257  double& totalw,
258  const GenLumiInfoProduct::XSec& thisRunXSec,
259  const double& thisw) {
260  double value = finalXSec.value();
261  double error = finalXSec.error();
262  double thisValue = thisRunXSec.value();
263  double thisError = thisRunXSec.error();
264  combine(value, error, totalw, thisValue, thisError, thisw);
265  finalXSec = GenLumiInfoProduct::XSec(value, error);
266  return;
267 }
268 
270  // sum of cross sections and errors over different processes
271  double sigSelSum = 0.0;
272  double err2SelSum = 0.0;
273  double sigSum = 0.0;
274  double err2Sum = 0.0;
275 
276  std::vector<GenLumiInfoProduct::XSec> tempVector_before;
277  std::vector<GenLumiInfoProduct::XSec> tempVector_after;
278 
279  // loop over different processes for each sample
280  unsigned int vectorSize = iLumiInfo.getProcessInfos().size();
281  for (unsigned int ip = 0; ip < vectorSize; ip++) {
283  double hepxsec_value = proc.lheXSec().value();
284  double hepxsec_error = proc.lheXSec().error() <= 0 ? 0 : proc.lheXSec().error();
285  tempVector_before.push_back(GenLumiInfoProduct::XSec(hepxsec_value, hepxsec_error));
286 
287  sigSelSum += hepxsec_value;
288  err2SelSum += hepxsec_error * hepxsec_error;
289 
290  // skips computation if jet matching efficiency=0
291  if (proc.killed().n() < 1) {
292  tempVector_after.push_back(GenLumiInfoProduct::XSec(0.0, 0.0));
293  continue;
294  }
295 
296  // computing jet matching efficiency for this process
297  double fracAcc = 0;
298  double ntotal = proc.nTotalPos() - proc.nTotalNeg();
299  double npass = proc.nPassPos() - proc.nPassNeg();
300  switch (hepidwtup_) {
301  case 3:
302  case -3:
303  fracAcc = ntotal > 0 ? npass / ntotal : -1;
304  break;
305  default:
306  fracAcc = proc.selected().sum() > 0 ? proc.killed().sum() / proc.selected().sum() : -1;
307  break;
308  }
309 
310  if (fracAcc <= 0) {
311  tempVector_after.push_back(GenLumiInfoProduct::XSec(0.0, 0.0));
312  continue;
313  }
314 
315  // cross section after matching for this particular process
316  double sigmaFin = hepxsec_value * fracAcc;
317 
318  // computing error on jet matching efficiency
319  double relErr = 1.0;
320  double efferr2 = 0;
321  switch (hepidwtup_) {
322  case 3:
323  case -3: {
324  double ntotal_pos = proc.nTotalPos();
325  double effp = ntotal_pos > 0 ? (double)proc.nPassPos() / ntotal_pos : 0;
326  double effp_err2 = ntotal_pos > 0 ? (1 - effp) * effp / ntotal_pos : 0;
327 
328  double ntotal_neg = proc.nTotalNeg();
329  double effn = ntotal_neg > 0 ? (double)proc.nPassNeg() / ntotal_neg : 0;
330  double effn_err2 = ntotal_neg > 0 ? (1 - effn) * effn / ntotal_neg : 0;
331 
332  efferr2 = ntotal > 0
333  ? (ntotal_pos * ntotal_pos * effp_err2 + ntotal_neg * ntotal_neg * effn_err2) / ntotal / ntotal
334  : 0;
335  break;
336  }
337  default: {
338  double denominator = pow(proc.selected().sum(), 4);
339  double passw = proc.killed().sum();
340  double passw2 = proc.killed().sum2();
341  double failw = proc.selected().sum() - passw;
342  double failw2 = proc.selected().sum2() - passw2;
343  double numerator = (passw2 * failw * failw + failw2 * passw * passw);
344 
345  efferr2 = denominator > 0 ? numerator / denominator : 0;
346  break;
347  }
348  }
349  double delta2Veto = efferr2 / fracAcc / fracAcc;
350 
351  // computing total error on cross section after matching efficiency
352 
353  double sigma2Sum, sigma2Err;
354  sigma2Sum = hepxsec_value * hepxsec_value;
355  sigma2Err = hepxsec_error * hepxsec_error;
356 
357  double delta2Sum = delta2Veto + sigma2Err / sigma2Sum;
358  relErr = (delta2Sum > 0.0 ? std::sqrt(delta2Sum) : 0.0);
359  double deltaFin = sigmaFin * relErr;
360 
361  tempVector_after.push_back(GenLumiInfoProduct::XSec(sigmaFin, deltaFin));
362 
363  // sum of cross sections and errors over different processes
364  sigSum += sigmaFin;
365  err2Sum += deltaFin * deltaFin;
366 
367  } // end of loop over different processes
368  tempVector_before.push_back(GenLumiInfoProduct::XSec(sigSelSum, sqrt(err2SelSum)));
369 
370  double total_matcheff = jetMatchEffStat_[10000].filterEfficiency(hepidwtup_);
371  double total_matcherr = jetMatchEffStat_[10000].filterEfficiencyError(hepidwtup_);
372 
373  double xsec_after = sigSelSum * total_matcheff;
374  double xsecerr_after = (total_matcheff > 0 && sigSelSum > 0)
375  ? xsec_after * sqrt(err2SelSum / sigSelSum / sigSelSum +
376  total_matcherr * total_matcherr / total_matcheff / total_matcheff)
377  : 0;
378 
379  GenLumiInfoProduct::XSec result(xsec_after, xsecerr_after);
380  tempVector_after.push_back(result);
381 
382  xsecBeforeMatching_ = tempVector_before;
383  xsecAfterMatching_ = tempVector_after;
384 
385  return result;
386 }
387 
389  edm::LogPrint("GenXSecAnalyzer") << "\n"
390  << "------------------------------------"
391  << "\n"
392  << "GenXsecAnalyzer:"
393  << "\n"
394  << "------------------------------------";
395 
396  if (jetMatchEffStat_.empty()) {
397  edm::LogPrint("GenXSecAnalyzer") << "------------------------------------"
398  << "\n"
399  << "Cross-section summary not available"
400  << "\n"
401  << "------------------------------------";
402  return;
403  }
404 
405  // fraction of negative weights
406  double final_fract_neg_w = 0;
407  double final_fract_neg_w_unc = 0;
408 
409  // below print out is only for combination of same physics MC samples and ME+Pythia MCs
410 
411  if (nMCs_ == 1 && hepidwtup_ != -1) {
412  edm::LogPrint("GenXSecAnalyzer")
413  << "-----------------------------------------------------------------------------------------------------------"
414  "--------------------------------------------------------------- \n"
415  << "Overall cross-section summary \n"
416  << "-----------------------------------------------------------------------------------------------------------"
417  "---------------------------------------------------------------";
418  edm::LogPrint("GenXSecAnalyzer") << "Process\t\txsec_before [pb]\t\tpassed\tnposw\tnnegw\ttried\tnposw\tnnegw "
419  "\txsec_match [pb]\t\t\taccepted [%]\t event_eff [%]";
420 
421  const unsigned sizeOfInfos = jetMatchEffStat_.size();
422  const unsigned last = sizeOfInfos - 1;
423  std::string* title = new std::string[sizeOfInfos];
424  unsigned int i = 0;
425  double jetmatch_eff = 0;
426  double jetmatch_err = 0;
427  double matching_eff = 1;
428  double matching_efferr = 1;
429 
430  for (std::map<int, GenFilterInfo>::const_iterator iter = jetMatchEffStat_.begin(); iter != jetMatchEffStat_.end();
431  ++iter, i++) {
432  GenFilterInfo thisJetMatchStat = iter->second;
433  GenFilterInfo thisEventEffStat =
434  GenFilterInfo(thisJetMatchStat.numPassPositiveEvents() + thisJetMatchStat.numPassNegativeEvents(),
435  0,
436  thisJetMatchStat.numTotalPositiveEvents() + thisJetMatchStat.numTotalNegativeEvents(),
437  0,
438  thisJetMatchStat.numPassPositiveEvents() + thisJetMatchStat.numPassNegativeEvents(),
439  thisJetMatchStat.numPassPositiveEvents() + thisJetMatchStat.numPassNegativeEvents(),
440  thisJetMatchStat.numTotalPositiveEvents() + thisJetMatchStat.numTotalNegativeEvents(),
441  thisJetMatchStat.numTotalPositiveEvents() + thisJetMatchStat.numTotalNegativeEvents());
442 
443  jetmatch_eff = thisJetMatchStat.filterEfficiency(hepidwtup_);
444  jetmatch_err = thisJetMatchStat.filterEfficiencyError(hepidwtup_);
445 
446  if (i == last) {
447  title[i] = "Total";
448 
449  edm::LogPrint("GenXSecAnalyzer")
450  << "-------------------------------------------------------------------------------------------------------"
451  "------------------------------------------------------------------- ";
452 
453  // fill negative fraction of negative weights and uncertainty after matching
454  final_fract_neg_w = thisEventEffStat.numEventsPassed() > 0
455  ? thisJetMatchStat.numPassNegativeEvents() / thisEventEffStat.numEventsPassed()
456  : 0;
457  final_fract_neg_w_unc =
458  thisJetMatchStat.numPassNegativeEvents() > 0
459  ? final_fract_neg_w * final_fract_neg_w / thisEventEffStat.numEventsPassed() *
460  sqrt(thisJetMatchStat.numPassPositiveEvents() * thisJetMatchStat.numPassPositiveEvents() /
461  thisJetMatchStat.numPassNegativeEvents() +
462  thisJetMatchStat.numPassPositiveEvents())
463  : 0;
464  } else {
465  title[i] = Form("%d", i);
466  }
467 
468  edm::LogPrint("GenXSecAnalyzer") << title[i] << "\t\t" << std::scientific << std::setprecision(3)
469  << xsecBeforeMatching_[i].value() << " +/- " << xsecBeforeMatching_[i].error()
470  << "\t\t" << thisEventEffStat.numEventsPassed() << "\t"
471  << thisJetMatchStat.numPassPositiveEvents() << "\t"
472  << thisJetMatchStat.numPassNegativeEvents() << "\t"
473  << thisEventEffStat.numEventsTotal() << "\t"
474  << thisJetMatchStat.numTotalPositiveEvents() << "\t"
475  << thisJetMatchStat.numTotalNegativeEvents() << "\t" << std::scientific
476  << std::setprecision(3) << xsecAfterMatching_[i].value() << " +/- "
477  << xsecAfterMatching_[i].error() << "\t\t" << std::fixed << std::setprecision(1)
478  << (jetmatch_eff * 100) << " +/- " << (jetmatch_err * 100) << "\t" << std::fixed
479  << std::setprecision(1) << (thisEventEffStat.filterEfficiency(+3) * 100)
480  << " +/- " << (thisEventEffStat.filterEfficiencyError(+3) * 100);
481 
482  matching_eff = thisEventEffStat.filterEfficiency(+3);
483  matching_efferr = thisEventEffStat.filterEfficiencyError(+3);
484  }
485  delete[] title;
486 
487  edm::LogPrint("GenXSecAnalyzer")
488  << "-----------------------------------------------------------------------------------------------------------"
489  "---------------------------------------------------------------";
490 
491  edm::LogPrint("GenXSecAnalyzer") << "Before matching: total cross section = " << std::scientific
492  << std::setprecision(3) << xsecBeforeMatching_[last].value() << " +- "
493  << xsecBeforeMatching_[last].error() << " pb";
494 
495  edm::LogPrint("GenXSecAnalyzer") << "After matching: total cross section = " << std::scientific
496  << std::setprecision(3) << xsecAfterMatching_[last].value() << " +- "
497  << xsecAfterMatching_[last].error() << " pb";
498 
499  edm::LogPrint("GenXSecAnalyzer") << "Matching efficiency = " << std::fixed << std::setprecision(1) << matching_eff
500  << " +/- " << matching_efferr << " [TO BE USED IN MCM]";
501 
502  } else if (hepidwtup_ == -1)
503  edm::LogPrint("GenXSecAnalyzer") << "Before Filter: total cross section = " << std::scientific
504  << std::setprecision(3) << xsecPreFilter_.value() << " +- "
505  << xsecPreFilter_.error() << " pb";
506 
507  // hepMC filter efficiency
508  double hepMCFilter_eff = 1.0;
509  double hepMCFilter_err = 0.0;
510  if (hepMCFilterEffStat_.sumWeights2() > 0) {
511  hepMCFilter_eff = hepMCFilterEffStat_.filterEfficiency(-1);
512  hepMCFilter_err = hepMCFilterEffStat_.filterEfficiencyError(-1);
513  edm::LogPrint("GenXSecAnalyzer") << "HepMC filter efficiency (taking into account weights)= "
514  << "(" << hepMCFilterEffStat_.sumPassWeights() << ")"
515  << " / "
516  << "(" << hepMCFilterEffStat_.sumWeights() << ")"
517  << " = " << std::scientific << std::setprecision(3) << hepMCFilter_eff << " +- "
518  << hepMCFilter_err;
519 
520  double hepMCFilter_event_total =
522  double hepMCFilter_event_pass =
524  double hepMCFilter_event_eff = hepMCFilter_event_total > 0 ? hepMCFilter_event_pass / hepMCFilter_event_total : 0;
525  double hepMCFilter_event_err =
526  hepMCFilter_event_total > 0
527  ? sqrt((1 - hepMCFilter_event_eff) * hepMCFilter_event_eff / hepMCFilter_event_total)
528  : -1;
529  edm::LogPrint("GenXSecAnalyzer") << "HepMC filter efficiency (event-level)= "
530  << "(" << hepMCFilter_event_pass << ")"
531  << " / "
532  << "(" << hepMCFilter_event_total << ")"
533  << " = " << std::scientific << std::setprecision(3) << hepMCFilter_event_eff
534  << " +- " << hepMCFilter_event_err;
535  }
536 
537  // gen-particle filter efficiency
538  if (filterOnlyEffStat_.sumWeights2() > 0) {
539  double filterOnly_eff = filterOnlyEffStat_.filterEfficiency(-1);
540  double filterOnly_err = filterOnlyEffStat_.filterEfficiencyError(-1);
541 
542  edm::LogPrint("GenXSecAnalyzer") << "Filter efficiency (taking into account weights)= "
543  << "(" << filterOnlyEffStat_.sumPassWeights() << ")"
544  << " / "
545  << "(" << filterOnlyEffStat_.sumWeights() << ")"
546  << " = " << std::scientific << std::setprecision(3) << filterOnly_eff << " +- "
547  << filterOnly_err;
548 
549  double filterOnly_event_total =
551  double filterOnly_event_pass =
553  double filterOnly_event_eff = filterOnly_event_total > 0 ? filterOnly_event_pass / filterOnly_event_total : 0;
554  double filterOnly_event_err = filterOnly_event_total > 0
555  ? sqrt((1 - filterOnly_event_eff) * filterOnly_event_eff / filterOnly_event_total)
556  : -1;
557  edm::LogPrint("GenXSecAnalyzer") << "Filter efficiency (event-level)= "
558  << "(" << filterOnly_event_pass << ")"
559  << " / "
560  << "(" << filterOnly_event_total << ")"
561  << " = " << std::scientific << std::setprecision(3) << filterOnly_event_eff
562  << " +- " << filterOnly_event_err << " [TO BE USED IN MCM]";
563 
564  // fill negative fraction of negative weights and uncertainty after filter
565  final_fract_neg_w =
566  filterOnly_event_pass > 0 ? filterOnlyEffStat_.numPassNegativeEvents() / (filterOnly_event_pass) : 0;
567  final_fract_neg_w_unc =
569  ? final_fract_neg_w * final_fract_neg_w / filterOnly_event_pass *
573  : 0;
574  }
575 
576  edm::LogPrint("GenXSecAnalyzer") << "\nAfter filter: final cross section = " << std::scientific
577  << std::setprecision(3) << xsec_.value() << " +- " << xsec_.error() << " pb";
578 
579  edm::LogPrint("GenXSecAnalyzer") << "After filter: final fraction of events with negative weights = "
580  << std::scientific << std::setprecision(3) << final_fract_neg_w << " +- "
581  << final_fract_neg_w_unc;
582 
583  // L=[N*(1-2f)^2]/s
584  double lumi_1M_evts =
585  xsec_.value() > 0 ? 1e6 * (1 - 2 * final_fract_neg_w) * (1 - 2 * final_fract_neg_w) / xsec_.value() / 1e3 : 0;
586  double lumi_1M_evts_unc =
587  xsec_.value() > 0 ? (1 - 2 * final_fract_neg_w) * lumi_1M_evts *
588  sqrt(1e-6 + 16 * pow(final_fract_neg_w_unc, 2) / pow(1 - 2 * final_fract_neg_w, 2) +
589  pow(xsec_.error() / xsec_.value(), 2))
590  : 0;
591  edm::LogPrint("GenXSecAnalyzer") << "After filter: final equivalent lumi for 1M events (1/fb) = " << std::scientific
592  << std::setprecision(3) << lumi_1M_evts << " +- " << lumi_1M_evts_unc;
593 }
std::vector< GenLumiInfoProduct::XSec > xsecBeforeMatching_
const std::vector< ProcessInfo > & getProcessInfos() const
std::map< int, GenFilterInfo > jetMatchEffStat_
void endJob() override
void combine(GenLumiInfoProduct::XSec &, double &, const GenLumiInfoProduct::XSec &, const double &)
edm::EDGetTokenT< GenFilterInfo > hepMCFilterInfoToken_
void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
unsigned int numTotalPositiveEvents() const
Definition: GenFilterInfo.h:26
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
double filterEfficiency(int idwtup=+3) const
GenFilterInfo filterOnlyEffStat_
edm::EDGetTokenT< GenLumiInfoProduct > genLumiInfoToken_
std::map< int, GenLumiInfoProduct::XSec > currentLumiBlockLHEXSec_
const lhef::HEPRUP & heprup() const
GenLumiInfoProduct::XSec xsec_
GenFilterInfo hepMCFilterEffRun_
GenFilterInfo hepMCFilterEffStat_
bool mergeProduct(GenFilterInfo const &other)
double sumPassWeights() const
Definition: GenFilterInfo.h:31
unsigned int numEventsPassed() const
Definition: GenFilterInfo.h:22
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Run.h:315
T sqrt(T t)
Definition: SSEVec.h:19
void analyze(const edm::Event &, const edm::EventSetup &) override
void beginJob() override
Definition: value.py:1
void beginRun(edm::Run const &, edm::EventSetup const &) override
bool isValid() const
Definition: HandleBase.h:70
GenLumiInfoProduct::XSec compute(const GenLumiInfoProduct &)
unsigned int numTotalNegativeEvents() const
Definition: GenFilterInfo.h:29
std::vector< double > XERRUP
Definition: LesHouches.h:118
unsigned int numEventsTotal() const
Definition: GenFilterInfo.h:23
T Max(T a, T b)
Definition: MathUtil.h:44
std::map< int, GenLumiInfoProduct::XSec > previousLumiBlockLHEXSec_
std::vector< double > XMAXUP
Definition: LesHouches.h:123
double filterEfficiencyError(int idwtup=+3) const
virtual bool mergeProduct(const GenLumiInfoProduct &other)
unsigned int numPassPositiveEvents() const
Definition: GenFilterInfo.h:25
void setLheXSec(double value, double err)
double sumWeights2() const
Definition: GenFilterInfo.h:38
GenFilterInfo filterOnlyEffRun_
edm::EDGetTokenT< GenFilterInfo > genFilterInfoToken_
double thisRunWeightPre_
GenXSecAnalyzer(const edm::ParameterSet &)
std::vector< GenLumiInfoProduct::XSec > xsecAfterMatching_
~GenXSecAnalyzer() override
void endRun(edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< LHERunInfoProduct > lheRunInfoToken_
std::vector< double > XSECUP
Definition: LesHouches.h:112
double sumWeights() const
Definition: GenFilterInfo.h:37
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
void setProcessInfo(const std::vector< ProcessInfo > &processes)
GenLumiInfoProduct::XSec xsecPreFilter_
Definition: Run.h:45
GenLumiInfoProduct product_
const int getHEPIDWTUP() const
unsigned int numPassNegativeEvents() const
Definition: GenFilterInfo.h:28
std::vector< int > LPRUP
Definition: LesHouches.h:128