CMS 3D CMS Logo

HepMCFilterDriver.cc
Go to the documentation of this file.
7 
9 
11  filter_(nullptr),
12  numEventsPassPos_(0),
13  numEventsPassNeg_(0),
14  numEventsTotalPos_(0),
15  numEventsTotalNeg_(0),
16  sumpass_w_(0.),
17  sumpass_w2_(0.),
18  sumtotal_w_(0.),
19  sumtotal_w2_(0.)
20 {
21 
22  std::string filterName = pset.getParameter<std::string>("filterName");
24 
25  if (filterName=="GenericDauHepMCFilter") {
26  filter_ = new GenericDauHepMCFilter(filterParameters);
27  }
28  else if (filterName=="PartonShowerBsHepMCFilter") {
29  filter_ = new PartonShowerBsHepMCFilter(filterParameters);
30  }
31  else if (filterName=="PartonShowerCsHepMCFilter") {
32  filter_ = new PartonShowerCsHepMCFilter(filterParameters);
33  }
34 
35  else if (filterName=="EmbeddingHepMCFilter") {
36  filter_ = new EmbeddingHepMCFilter(filterParameters);
37  }
38  else if (filterName=="PythiaHepMCFilterGammaGamma") {
39  filter_ = new PythiaHepMCFilterGammaGamma(filterParameters);
40  }
41  else {
42  throw edm::Exception(edm::errors::Configuration,"HepMCFilterDriver") << "Invalid HepMCFilter name:" << filterName;
43  }
44 
45 }
46 
48 {
49  if (filter_) delete filter_;
50 
51 }
52 
54 {
55  if(weight>0)
57  else
59 
61  sumtotal_w2_ += weight*weight;
62 
63 
64  bool accepted = filter_->filter(evt);
65 
66  if (accepted) {
67 
68  if(weight>0)
70  else
72  sumpass_w_ += weight;
73  sumpass_w2_ += weight*weight;
74 
75  }
76 
77  return accepted;
78 }
79 
81 {
82 
83  unsigned int ntried_ = numEventsTotalPos_ + numEventsTotalNeg_;
84  unsigned int naccepted_ = numEventsPassPos_ + numEventsPassNeg_;
85  printf("ntried = %i, naccepted = %i, efficiency = %5f\n",ntried_,naccepted_,(double)naccepted_/(double)ntried_);
86  printf("weighttried = %5f, weightaccepted = %5f, efficiency = %5f\n",sumtotal_w_,sumpass_w_,sumpass_w_/sumtotal_w_);
87 
88 }
89 
90 
92 
97  sumpass_w_ = 0;
98  sumpass_w2_ = 0;
99  sumtotal_w_ = 0;
100  sumtotal_w2_ = 0;
101 
102 }
T getParameter(std::string const &) const
unsigned int numEventsPassPos_
unsigned int numEventsTotalPos_
#define nullptr
Definition: weight.py:1
HepMCFilterDriver(const edm::ParameterSet &)
unsigned int numEventsPassNeg_
BaseHepMCFilter * filter_
void statistics() const
bool accepted(std::vector< std::string_view > const &, std::string_view)
unsigned int numEventsTotalNeg_
virtual bool filter(const HepMC::GenEvent *evt)=0
bool filter(const HepMC::GenEvent *evt, double weight=1.)