test
CMS 3D CMS Logo

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