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 
16 
20 
23 
30 
31 using namespace lhef;
32 
33 class LHEProducer : public edm::one::EDFilter<edm::EndRunProducer,
34  edm::one::WatchRuns> {
35  public:
36  explicit LHEProducer(const edm::ParameterSet &params);
37  virtual ~LHEProducer();
38 
39  protected:
40  virtual void beginJob() override;
41  virtual void endJob() override;
42  virtual void beginRun(edm::Run const& run, const edm::EventSetup &es) override;
43  virtual void endRun(edm::Run const&run, const edm::EventSetup &es) override;
44  virtual void endRunProduce(edm::Run &run, const edm::EventSetup &es) override;
45  virtual bool filter(edm::Event &event, const edm::EventSetup &es) override;
46 
47  private:
48  double matching(const HepMC::GenEvent *event, bool shower) const;
49 
50  bool showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event);
51  void onInit();
52  void onBeforeHadronisation();
53 
54  unsigned int eventsToPrint;
55  std::vector<int> removeResonances;
56  std::auto_ptr<Hadronisation> hadronisation;
57  std::auto_ptr<JetMatching> jetMatching;
58 
59  double extCrossSect;
60  double extFilterEff;
62 
63  boost::shared_ptr<LHEEvent> partonLevel;
64  boost::shared_ptr<LHERunInfo> runInfo;
65  unsigned int index;
67  double weight;
69 };
70 
72  eventsToPrint(params.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
73  extCrossSect(params.getUntrackedParameter<double>("crossSection", -1.0)),
74  extFilterEff(params.getUntrackedParameter<double>("filterEfficiency", -1.0)),
75  matchSummary(false)
76 {
78  params.getParameter<edm::ParameterSet>("hadronisation"));
79 
80  if (params.exists("removeResonances"))
82  params.getParameter<std::vector<int> >(
83  "removeResonances");
84 
85  std::set<std::string> matchingCapabilities;
86  if (params.exists("jetMatching")) {
87  edm::ParameterSet jetParams =
89  "jetMatching");
90  jetMatching = JetMatching::create(jetParams);
91 
92  matchingCapabilities = jetMatching->capabilities();
93  hadronisation->matchingCapabilities(matchingCapabilities);
94  }
95 
96  produces<edm::HepMCProduct>();
97  produces<GenEventInfoProduct>();
98  produces<GenRunInfoProduct, edm::InRun>();
99 
100  if (jetMatching.get()) {
101  if (params.getUntrackedParameter<bool>(
102  "preferShowerVetoCallback", true))
103  hadronisation->onShoweredEvent().connect(
104  sigc::mem_fun(*this,
106  hadronisation->onInit().connect(
107  sigc::mem_fun(*this, &LHEProducer::onInit));
108  hadronisation->onBeforeHadronisation().connect(
109  sigc::mem_fun(*this,
111 
112  matchSummary = matchingCapabilities.count("matchSummary");
113  if (matchSummary) {
114  produces< std::vector<double> >("matchDeltaR");
115  produces< std::vector<double> >("matchDeltaPRel");
116  }
117  }
118 
119  // force total branching ratio for QCD/QED to 1
120  for(int i = 0; i < 6; i++)
122  for(int i = 9; i < 23; i++)
124 }
125 
127 {
128 }
129 
131 {
132  hadronisation->init();
133 }
134 
136 {
137  hadronisation.reset();
138  jetMatching.reset();
139 }
140 
142 {
144  run.getByLabel("source", product);
145 
146  runInfo.reset(new LHERunInfo(*product));
147  index = 0;
148 }
150 {
151 }
152 
154 {
155  hadronisation->statistics();
156 
157  LHERunInfo::XSec crossSection;
158  if (runInfo) {
159  crossSection = runInfo->xsec();
160  runInfo->statistics();
161  }
162 
163  std::auto_ptr<GenRunInfoProduct> runInfo(new GenRunInfoProduct);
164 
165  runInfo->setInternalXSec(
166  GenRunInfoProduct::XSec(crossSection.value,
167  crossSection.error));
168  runInfo->setExternalXSecLO(extCrossSect);
169  runInfo->setFilterEfficiency(extFilterEff);
170 
171  run.put(runInfo);
172 
173  runInfo.reset();
174 }
175 
177 {
178  std::auto_ptr<edm::HepMCProduct> result(new edm::HepMCProduct);
179 
181  event.getByLabel("source", product);
182 
183  partonLevel.reset(new LHEEvent(runInfo, *product));
184  if (!removeResonances.empty())
185  partonLevel->removeResonances(removeResonances);
186 
187  if (product->pdf())
188  partonLevel->setPDF(
189  std::auto_ptr<LHEEvent::PDF>(
190  new LHEEvent::PDF(*product->pdf())));
191 
192  hadronisation->setEvent(partonLevel);
193 
194  double br = branchingRatios.getFactor(hadronisation.get(),
195  partonLevel);
196 
197  matchingDone = false;
198  weight = 1.0;
199  std::auto_ptr<HepMC::GenEvent> hadronLevel(hadronisation->hadronize());
200 
201  if (!hadronLevel.get()) {
202  if (matchingDone) {
203  if (weight == 0.0)
204  partonLevel->count(LHERunInfo::kSelected, br);
205  else
206  partonLevel->count(LHERunInfo::kKilled,
207  br, weight);
208  } else
209  partonLevel->count(LHERunInfo::kTried, br);
210  }
211 
212  if (!matchingDone && jetMatching.get() && hadronLevel.get())
213  weight = matching(hadronLevel.get(), false);
214 
215  if (weight == 0.0) {
216  edm::LogInfo("Generator|LHEInterface")
217  << "Event got rejected by the "
218  "jet matching." << std::endl;
219 
220  if (hadronLevel.get()) {
221  partonLevel->count(LHERunInfo::kSelected);
222  hadronLevel.reset();
223  }
224  }
225 
226  if (!hadronLevel.get()) {
227  event.put(result);
228  std::auto_ptr<GenEventInfoProduct> info(
229  new GenEventInfoProduct);
230  event.put(info);
231  return false;
232  }
233 
234  partonLevel->count(LHERunInfo::kAccepted, br, weight);
235 
236  hadronLevel->set_event_number(++index);
237 
238  if (eventsToPrint) {
239  eventsToPrint--;
240  hadronLevel->print();
241  }
242 
243  std::auto_ptr<GenEventInfoProduct> info(
244  new GenEventInfoProduct(hadronLevel.get()));
245  result->addHepMCData(hadronLevel.release());
246  event.put(result);
247  event.put(info);
248 
249  if (jetMatching.get() && matchSummary) {
250  std::auto_ptr< std::vector<double> > matchDeltaR(
251  new std::vector<double>);
252  std::auto_ptr< std::vector<double> > matchDeltaPRel(
253  new std::vector<double>);
254 
255  typedef std::vector<JetMatching::JetPartonMatch> Matches;
256  Matches matches = jetMatching->getMatchSummary();
257 
258  for(Matches::const_iterator iter = matches.begin();
259  iter != matches.end(); ++iter) {
260  if (!iter->isMatch())
261  continue;
262 
263  matchDeltaR->push_back(iter->delta);
264  matchDeltaPRel->push_back(iter->jet.rho() /
265  iter->parton.rho() - 1.0);
266  }
267 
268  event.put(matchDeltaR, "matchDeltaR");
269  event.put(matchDeltaPRel, "matchDeltaPRel");
270  }
271 
272  return true;
273 }
274 
275 double LHEProducer::matching(const HepMC::GenEvent *event, bool shower) const
276 {
277  if (!jetMatching.get())
278  return 1.0;
279 
280  return jetMatching->match(partonLevel->asHepMCEvent().get(),
281  event, shower);
282 }
283 
284 bool LHEProducer::showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event)
285 {
286  weight = matching(event.get(), true);
287  matchingDone = true;
288  return weight == 0.0;
289 }
290 
292 {
293  jetMatching->init(runInfo);
294 }
295 
297 {
298  jetMatching->beforeHadronisation(partonLevel);
299 }
300 
302 
virtual void endRun(edm::Run const &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:149
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:176
double extFilterEff
Definition: LHEProducer.cc:60
double extCrossSect
Definition: LHEProducer.cc:59
#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:126
void beginJob()
Definition: Breakpoints.cc:15
std::vector< int > removeResonances
Definition: LHEProducer.cc:55
double weight
Definition: LHEProducer.cc:67
BranchingRatios branchingRatios
Definition: LHEProducer.cc:68
void onBeforeHadronisation()
Definition: LHEProducer.cc:296
tuple result
Definition: query.py:137
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:141
LHEProducer(const edm::ParameterSet &params)
Definition: LHEProducer.cc:71
double matching(const HepMC::GenEvent *event, bool shower) const
Definition: LHEProducer.cc:275
virtual void endRunProduce(edm::Run &run, const edm::EventSetup &es) override
Definition: LHEProducer.cc:153
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:66
boost::shared_ptr< LHEEvent > partonLevel
Definition: LHEProducer.cc:63
bool matchSummary
Definition: LHEProducer.cc:61
void onInit()
Definition: LHEProducer.cc:291
virtual void endJob() override
Definition: LHEProducer.cc:135
virtual void beginJob() override
Definition: LHEProducer.cc:130
std::auto_ptr< JetMatching > jetMatching
Definition: LHEProducer.cc:57
unsigned int index
Definition: LHEProducer.cc:65
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:64
unsigned int eventsToPrint
Definition: LHEProducer.cc:54
std::auto_ptr< Hadronisation > hadronisation
Definition: LHEProducer.cc:56
SurfaceDeformation * create(int type, const std::vector< double > &params)
bool showeredEvent(const boost::shared_ptr< HepMC::GenEvent > &event)
Definition: LHEProducer.cc:284
Definition: Run.h:41