CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ThePEGHadronizer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <sstream>
3 
4 #include <HepMC/GenEvent.h>
5 #include <HepMC/IO_BaseClass.h>
6 
7 #include <ThePEG/Repository/Repository.h>
8 #include <ThePEG/EventRecord/Event.h>
9 #include <ThePEG/Config/ThePEG.h>
10 #include <ThePEG/LesHouches/LesHouchesReader.h>
11 
13 
17 
21 
23 
25 
26 namespace CLHEP {
27  class HepRandomEngine;
28 }
29 
31  public:
32  ThePEGHadronizer(const edm::ParameterSet &params);
33  virtual ~ThePEGHadronizer();
34 
35  bool readSettings( int ) { return true; }
38  bool declareStableParticles(const std::vector<int> &pdgIds);
39  bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
40 
41  void statistics();
42 
44  bool hadronize();
45  bool decay();
46  bool residualDecay();
47  void finalizeEvent();
48 
49  const char *classname() const { return "ThePEGHadronizer"; }
50 
51  private:
52 
53  virtual void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }
54 
55  unsigned int eventsToPrint;
56 
57  ThePEG::EventPtr thepegEvent;
58 
59  boost::shared_ptr<lhef::LHEProxy> proxy_;
61 };
62 
64  ThePEGInterface(pset),
65  BaseHadronizer(pset),
66  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
67  handlerDirectory_(pset.getParameter<std::string>("eventHandlers"))
68 {
69  initRepository(pset);
70 
71 }
72 
74 {
75 }
76 
78 {
79  initGenerator();
80  return true;
81 }
82 
84 {
86 
87  std::ostringstream ss;
88  ss << proxy_->getID();
89 
90  std::ostringstream logstream;
91  ThePEG::Repository::exec("set " + handlerDirectory_ +
92  "/LHEReader:ProxyID " + ss.str(), logstream);
93  edm::LogInfo("Generator|LHEInterface") << logstream.str();
94 
95  proxy_->loadRunInfo(getLHERunInfo());
96 
97  initGenerator();
98  return true;
99 }
100 
101 bool ThePEGHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
102 {
103  return false;
104 }
105 
107 {
109  eg_->integratedXSec() / ThePEG::picobarn,
110  eg_->integratedXSecErr() / ThePEG::picobarn));
111 }
112 
114 {
115  edm::LogInfo("Generator|ThePEGHadronizer") << "Start production";
116 
118 
119  try {
120  thepegEvent = eg_->shoot();
121  } catch (std::exception& exc) {
122  edm::LogWarning("Generator|ThePEGHadronizer") << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
123  return false;
124  } catch (...) {
125  edm::LogWarning("Generator|ThePEGHadronizer") << "EGPtr::shoot() thrown an unknown exception, event skipped";
126  return false;
127  }
128 
129  if (!thepegEvent) {
130  edm::LogWarning("Generator|ThePEGHadronizer") << "thepegEvent not initialized";
131  return false;
132  }
133 
135  if (!event().get()) {
136  edm::LogWarning("Generator|ThePEGHadronizer") << "genEvent not initialized";
137  return false;
138  }
139 
140  return true;
141 }
142 
144 {
145  edm::LogInfo("Generator|ThePEGHadronizer") << "Start production";
146 
148 
149  //need to copy lhe event here unfortunately because of interface mismatch
150  proxy_->loadEvent(boost::shared_ptr<lhef::LHEEvent>(new lhef::LHEEvent(*lheEvent())));
151 
152  //dummy for now
153  double mergeweight = 1.0;
154 
155  try {
156  thepegEvent = eg_->shoot();
157  } catch (std::exception& exc) {
158  edm::LogWarning("Generator|ThePEGHadronizer") << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
159  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
160  return false;
161  } catch (...) {
162  edm::LogWarning("Generator|ThePEGHadronizer") << "EGPtr::shoot() thrown an unknown exception, event skipped";
163  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
164  return false;
165  }
166 
167  if (!thepegEvent) {
168  edm::LogWarning("Generator|ThePEGHadronizer") << "thepegEvent not initialized";
169  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
170  return false;
171  }
172 
174  if (!event().get()) {
175  edm::LogWarning("Generator|ThePEGHadronizer") << "genEvent not initialized";
176  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
177  return false;
178  }
179 
180  //Fill LHE weight (since it's not otherwise propagated)
181  event()->weights()[0] *= lheEvent()->getHEPEUP()->XWGTUP;
182 
183  HepMC::PdfInfo pdf;
184  clearAuxiliary(event().get(), &pdf);
185  lheEvent()->fillPdfInfo(&pdf);
186  fillAuxiliary(event().get(), &pdf, thepegEvent);
187  event()->set_pdf_info(pdf);
188 
189  // update LHE matching statistics
190  //
191  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
192 
193  return true;
194 
195 }
196 
198 {
199  HepMC::PdfInfo pdf;
200  clearAuxiliary(event().get(), &pdf);
201  fillAuxiliary(event().get(), &pdf, thepegEvent);
202  event()->set_pdf_info(pdf);
203 
204  eventInfo().reset(new GenEventInfoProduct(event().get()));
205  eventInfo()->setBinningValues(
206  std::vector<double>(1, pthat(thepegEvent)));
207 
208  if (eventsToPrint) {
209  eventsToPrint--;
210  event()->print();
211  }
212 
213  if (iobc_.get())
214  iobc_->write_event(event().get());
215 
216  edm::LogInfo("Generator|ThePEGHadronizer") << "Event produced";
217 }
218 
220 {
221  return true;
222 }
223 
225 {
226  return true;
227 }
228 
230 
233 
unsigned int eventsToPrint
void initRepository(const edm::ParameterSet &params) const
edm::HadronizerFilter< ThePEGHadronizer, gen::ExternalDecayDriver > ThePEGHadronizerFilter
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:207
const std::string handlerDirectory_
static double pthat(const ThePEG::EventPtr &event)
void setInternalXSec(const XSec &xsec)
const HEPEUP * getHEPEUP() const
Definition: LHEEvent.h:43
std::auto_ptr< HepMC::GenEvent > & event()
bool declareStableParticles(const std::vector< int > &pdgIds)
GenRunInfoProduct & runInfo()
bool initializeForInternalPartons()
lhef::LHEEvent * lheEvent()
edm::GeneratorFilter< ThePEGHadronizer, gen::ExternalDecayDriver > ThePEGGeneratorFilter
void fillPdfInfo(HepMC::PdfInfo *info) const
Definition: LHEEvent.cc:221
void flushRandomNumberGenerator()
ThePEG::EventPtr thepegEvent
std::auto_ptr< HepMC::IO_BaseClass > iobc_
ThePEG::EGPtr eg_
static std::auto_ptr< HepMC::GenEvent > convert(const ThePEG::EventPtr &event)
bool generatePartonsAndHadronize()
std::auto_ptr< GenEventInfoProduct > & eventInfo()
ThePEGHadronizer(const edm::ParameterSet &params)
virtual ~ThePEGHadronizer()
bool initializeForExternalPartons()
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
boost::shared_ptr< lhef::LHEProxy > proxy_
const char * classname() const
bool declareSpecialSettings(const std::vector< std::string >)
void setPEGRandomEngine(CLHEP::HepRandomEngine *)
double XWGTUP
Definition: LesHouches.h:196
const boost::shared_ptr< lhef::LHERunInfo > & getLHERunInfo() const
static void clearAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf)
static boost::shared_ptr< LHEProxy > create()
Definition: LHEProxy.cc:41
static void fillAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf, const ThePEG::EventPtr &event)