CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LHEProducer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <memory>
3 
4 #include <boost/shared_ptr.hpp>
5 #include <sigc++/signal.h>
6 
7 #include "HepMC/GenEvent.h"
8 #include "HepMC/SimpleVector.h"
9 
17 
21 
24 
31 
32 using namespace lhef;
33 
34 class LHEProducer : public edm::one::EDFilter<edm::EndRunProducer,
35  edm::one::WatchRuns> {
36  public:
37  explicit LHEProducer(const edm::ParameterSet &params);
38  virtual ~LHEProducer();
39 
40  protected:
41  virtual void beginJob() override;
42  virtual void endJob() override;
43  virtual void beginRun(edm::Run const& run, const edm::EventSetup &es) override;
44  virtual void endRun(edm::Run const&run, const edm::EventSetup &es) override;
45  virtual void endRunProduce(edm::Run &run, const edm::EventSetup &es) override;
46  virtual bool filter(edm::Event &event, const edm::EventSetup &es) override;
47 
48  private:
49  double matching(const HepMC::GenEvent *event, bool shower) const;
50 
51  bool showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event);
52  void onInit();
53  void onBeforeHadronisation();
54 
55  unsigned int eventsToPrint;
56  std::vector<int> removeResonances;
57  std::auto_ptr<Hadronisation> hadronisation;
58  std::auto_ptr<JetMatching> jetMatching;
59 
60  double extCrossSect;
61  double extFilterEff;
63 
64  boost::shared_ptr<LHEEvent> partonLevel;
65  boost::shared_ptr<LHERunInfo> runInfo;
66  unsigned int index;
68  double weight;
70 };
71 
73  eventsToPrint(params.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
74  extCrossSect(params.getUntrackedParameter<double>("crossSection", -1.0)),
75  extFilterEff(params.getUntrackedParameter<double>("filterEfficiency", -1.0)),
76  matchSummary(false)
77 {
79  params.getParameter<edm::ParameterSet>("hadronisation"));
80 
81  if (params.exists("removeResonances"))
83  params.getParameter<std::vector<int> >(
84  "removeResonances");
85 
86  std::set<std::string> matchingCapabilities;
87  if (params.exists("jetMatching")) {
90  "jetMatching");
91  jetMatching = JetMatching::create(jetParams);
92 
93  matchingCapabilities = jetMatching->capabilities();
94  hadronisation->matchingCapabilities(matchingCapabilities);
95  }
96 
97  produces<edm::HepMCProduct>();
98  produces<GenEventInfoProduct>();
99  produces<GenRunInfoProduct, edm::InRun>();
100 
101  if (jetMatching.get()) {
102  if (params.getUntrackedParameter<bool>(
103  "preferShowerVetoCallback", true))
104  hadronisation->onShoweredEvent().connect(
105  sigc::mem_fun(*this,
107  hadronisation->onInit().connect(
108  sigc::mem_fun(*this, &LHEProducer::onInit));
109  hadronisation->onBeforeHadronisation().connect(
110  sigc::mem_fun(*this,
112 
113  matchSummary = matchingCapabilities.count("matchSummary");
114  if (matchSummary) {
115  produces< std::vector<double> >("matchDeltaR");
116  produces< std::vector<double> >("matchDeltaPRel");
117  }
118  }
119 
120  // force total branching ratio for QCD/QED to 1
121  for(int i = 0; i < 6; i++)
123  for(int i = 9; i < 23; i++)
125 }
126 
128 {
129 }
130 
132 {
133  hadronisation->init();
134 }
135 
137 {
138  hadronisation.reset();
139  jetMatching.reset();
140 }
141 
143 {
145  run.getByLabel("source", product);
146 
147  runInfo.reset(new LHERunInfo(*product));
148  index = 0;
149 }
151 {
152 }
153 
155 {
156  hadronisation->statistics();
157 
158  LHERunInfo::XSec crossSection;
159  if (runInfo) {
160  crossSection = runInfo->xsec();
161  runInfo->statistics();
162  }
163 
164  std::auto_ptr<GenRunInfoProduct> runInfo(new GenRunInfoProduct);
165 
166  runInfo->setInternalXSec(
167  GenRunInfoProduct::XSec(crossSection.value(),
168  crossSection.error()));
169  runInfo->setExternalXSecLO(extCrossSect);
170  runInfo->setFilterEfficiency(extFilterEff);
171 
172  run.put(runInfo);
173 
174  runInfo.reset();
175 }
176 
178 {
179  edm::RandomEngineSentry<Hadronisation> randomEngineSentry(hadronisation.get(), event.streamID());
180 
181  std::auto_ptr<edm::HepMCProduct> result(new edm::HepMCProduct);
182 
184  event.getByLabel("source", product);
185 
186  partonLevel.reset(new LHEEvent(runInfo, *product));
187  if (!removeResonances.empty())
188  partonLevel->removeResonances(removeResonances);
189 
190  if (product->pdf())
191  partonLevel->setPDF(
192  std::auto_ptr<LHEEvent::PDF>(
193  new LHEEvent::PDF(*product->pdf())));
194 
195  hadronisation->setEvent(partonLevel);
196 
197  double br = branchingRatios.getFactor(hadronisation.get(),
198  partonLevel);
199 
200  matchingDone = false;
201  weight = 1.0;
202  std::auto_ptr<HepMC::GenEvent> hadronLevel(hadronisation->hadronize());
203 
204  if (!hadronLevel.get()) {
205  if (matchingDone) {
206  if (weight == 0.0)
207  partonLevel->count(LHERunInfo::kSelected, br);
208  else
209  partonLevel->count(LHERunInfo::kKilled,
210  br, weight);
211  } else
212  partonLevel->count(LHERunInfo::kTried, br);
213  }
214 
215  if (!matchingDone && jetMatching.get() && hadronLevel.get())
216  weight = matching(hadronLevel.get(), false);
217 
218  if (weight == 0.0) {
219  edm::LogInfo("Generator|LHEInterface")
220  << "Event got rejected by the "
221  "jet matching." << std::endl;
222 
223  if (hadronLevel.get()) {
224  partonLevel->count(LHERunInfo::kSelected);
225  hadronLevel.reset();
226  }
227  }
228 
229  if (!hadronLevel.get()) {
230  event.put(result);
231  std::auto_ptr<GenEventInfoProduct> info(
232  new GenEventInfoProduct);
233  event.put(info);
234  return false;
235  }
236 
237  partonLevel->count(LHERunInfo::kAccepted, br, weight);
238 
239  hadronLevel->set_event_number(++index);
240 
241  if (eventsToPrint) {
242  eventsToPrint--;
243  hadronLevel->print();
244  }
245 
246  std::auto_ptr<GenEventInfoProduct> info(
247  new GenEventInfoProduct(hadronLevel.get()));
248  result->addHepMCData(hadronLevel.release());
249  event.put(result);
250  event.put(info);
251 
252  if (jetMatching.get() && matchSummary) {
253  std::auto_ptr< std::vector<double> > matchDeltaR(
254  new std::vector<double>);
255  std::auto_ptr< std::vector<double> > matchDeltaPRel(
256  new std::vector<double>);
257 
258  typedef std::vector<JetMatching::JetPartonMatch> Matches;
259  Matches matches = jetMatching->getMatchSummary();
260 
261  for(Matches::const_iterator iter = matches.begin();
262  iter != matches.end(); ++iter) {
263  if (!iter->isMatch())
264  continue;
265 
266  matchDeltaR->push_back(iter->delta);
267  matchDeltaPRel->push_back(iter->jet.rho() /
268  iter->parton.rho() - 1.0);
269  }
270 
271  event.put(matchDeltaR, "matchDeltaR");
272  event.put(matchDeltaPRel, "matchDeltaPRel");
273  }
274 
275  return true;
276 }
277 
278 double LHEProducer::matching(const HepMC::GenEvent *event, bool shower) const
279 {
280  if (!jetMatching.get())
281  return 1.0;
282 
283  return jetMatching->match(partonLevel->asHepMCEvent().get(),
284  event, shower);
285 }
286 
287 bool LHEProducer::showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event)
288 {
289  weight = matching(event.get(), true);
290  matchingDone = true;
291  return weight == 0.0;
292 }
293 
295 {
296  jetMatching->init(runInfo);
297 }
298 
300 {
301  jetMatching->beforeHadronisation(partonLevel);
302 }
303 
305 
virtual void endRun(edm::Run const &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:150
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:200
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
static const TGPicture * info(bool iBackgroundIsBlack)
void set(int pdgId, bool both=true, double value=1.0)
virtual bool filter(edm::Event &event, const edm::EventSetup &es) override
Definition: LHEProducer.cc:177
double extFilterEff
Definition: LHEProducer.cc:61
double extCrossSect
Definition: LHEProducer.cc:60
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double getFactor(const Hadronisation *hadronisation, const boost::shared_ptr< LHEEvent > &event) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual ~LHEProducer()
Definition: LHEProducer.cc:127
void beginJob()
Definition: Breakpoints.cc:15
std::vector< int > removeResonances
Definition: LHEProducer.cc:56
double weight
Definition: LHEProducer.cc:68
BranchingRatios branchingRatios
Definition: LHEProducer.cc:69
void onBeforeHadronisation()
Definition: LHEProducer.cc:299
tuple result
Definition: query.py:137
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:142
LHEProducer(const edm::ParameterSet &params)
Definition: LHEProducer.cc:72
double matching(const HepMC::GenEvent *event, bool shower) const
Definition: LHEProducer.cc:278
virtual void endRunProduce(edm::Run &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:154
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool matchingDone
Definition: LHEProducer.cc:67
boost::shared_ptr< LHEEvent > partonLevel
Definition: LHEProducer.cc:64
bool matchSummary
Definition: LHEProducer.cc:62
void onInit()
Definition: LHEProducer.cc:294
virtual void endJob() override
Definition: LHEProducer.cc:136
virtual void beginJob() override
Definition: LHEProducer.cc:131
std::auto_ptr< JetMatching > jetMatching
Definition: LHEProducer.cc:58
unsigned int index
Definition: LHEProducer.cc:66
void put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Run.h:107
DEFINE_LHE_JETMATCHING_PLUGIN(JetMatchingMLM)
volatile std::atomic< bool > shutdown_flag false
boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEProducer.cc:65
unsigned int eventsToPrint
Definition: LHEProducer.cc:55
std::auto_ptr< Hadronisation > hadronisation
Definition: LHEProducer.cc:57
SurfaceDeformation * create(int type, const std::vector< double > &params)
bool showeredEvent(const boost::shared_ptr< HepMC::GenEvent > &event)
Definition: LHEProducer.cc:287
Definition: Run.h:41