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 
3 #include <HepMC/GenEvent.h>
4 #include <HepMC/IO_BaseClass.h>
5 
6 #include <ThePEG/Repository/Repository.h>
7 #include <ThePEG/EventRecord/Event.h>
8 #include <ThePEG/Config/ThePEG.h>
9 
11 
15 
19 
21 
22 namespace CLHEP {
23  class HepRandomEngine;
24 }
25 
27  public:
28  ThePEGHadronizer(const edm::ParameterSet &params);
29  virtual ~ThePEGHadronizer();
30 
31  bool readSettings( int ) { return true; }
33 // bool initializeForExternalPartons();
34  bool declareStableParticles(const std::vector<int> &pdgIds);
35  bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
36 
37  void statistics();
38 
40 // bool hadronize();
41  bool decay();
42  bool residualDecay();
43  void finalizeEvent();
44 
45  const char *classname() const { return "ThePEGHadronizer"; }
46 
47  private:
48 
49  virtual void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }
50 
51  unsigned int eventsToPrint;
52 
53  ThePEG::EventPtr thepegEvent;
54 };
55 
57  ThePEGInterface(pset),
58  BaseHadronizer(pset),
59  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0))
60 {
61  initRepository(pset);
62 
63 }
64 
66 {
67 }
68 
70 {
71  initGenerator();
72  return true;
73 }
74 
75 bool ThePEGHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
76 {
77  return false;
78 }
79 
81 {
83  eg_->integratedXSec() / ThePEG::picobarn,
84  eg_->integratedXSecErr() / ThePEG::picobarn));
85 }
86 
88 {
89  edm::LogInfo("Generator|ThePEGHadronizer") << "Start production";
90 
92 
93  thepegEvent = eg_->shoot();
94  if (!thepegEvent) {
95  edm::LogWarning("Generator|ThePEGHadronizer") << "thepegEvent not initialized";
96  return false;
97  }
98 
100  if (!event().get()) {
101  edm::LogWarning("Generator|ThePEGHadronizer") << "genEvent not initialized";
102  return false;
103  }
104 
105  return true;
106 }
107 
109 {
110  HepMC::PdfInfo pdf;
111  clearAuxiliary(event().get(), &pdf);
112  fillAuxiliary(event().get(), &pdf, thepegEvent);
113  event()->set_pdf_info(pdf);
114 
115  eventInfo().reset(new GenEventInfoProduct(event().get()));
116  eventInfo()->setBinningValues(
117  std::vector<double>(1, pthat(thepegEvent)));
118 
119  if (eventsToPrint) {
120  eventsToPrint--;
121  event()->print();
122  }
123 
124  if (iobc_.get())
125  iobc_->write_event(event().get());
126 
127  edm::LogInfo("Generator|ThePEGHadronizer") << "Event produced";
128 }
129 
131 {
132  return true;
133 }
134 
136 {
137  return true;
138 }
139 
141 
144 
145 #if 0
147 DEFINE_FWK_MODULE(ThePEGHadronizerFilter);
148 #endif
unsigned int eventsToPrint
void initRepository(const edm::ParameterSet &params) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static double pthat(const ThePEG::EventPtr &event)
void setInternalXSec(const XSec &xsec)
std::auto_ptr< HepMC::GenEvent > & event()
bool declareStableParticles(const std::vector< int > &pdgIds)
GenRunInfoProduct & runInfo()
bool initializeForInternalPartons()
edm::GeneratorFilter< ThePEGHadronizer, gen::ExternalDecayDriver > ThePEGGeneratorFilter
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()
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
const char * classname() const
bool declareSpecialSettings(const std::vector< std::string >)
void setPEGRandomEngine(CLHEP::HepRandomEngine *)
static void clearAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf)
static void fillAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf, const ThePEG::EventPtr &event)