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 
23  public:
24  ThePEGHadronizer(const edm::ParameterSet &params);
25  virtual ~ThePEGHadronizer();
26 
27  bool readSettings( int ) { return true; }
29 // bool initializeForExternalPartons();
30  bool declareStableParticles(const std::vector<int> &pdgIds);
31  bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
32 
33  void statistics();
34 
36 // bool hadronize();
37  bool decay();
38  bool residualDecay();
39  void finalizeEvent();
40 
41  const char *classname() const { return "ThePEGHadronizer"; }
42 
43  private:
44  unsigned int eventsToPrint;
45 
46  ThePEG::EventPtr thepegEvent;
47 };
48 
50  ThePEGInterface(pset),
51  BaseHadronizer(pset),
52  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0))
53 {
54  initRepository(pset);
55 
56 }
57 
59 {
60 }
61 
63 {
64  initGenerator();
65  return true;
66 }
67 
68 bool ThePEGHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
69 {
70  return false;
71 }
72 
74 {
76  eg_->integratedXSec() / ThePEG::picobarn,
77  eg_->integratedXSecErr() / ThePEG::picobarn));
78 }
79 
81 {
82  edm::LogInfo("Generator|ThePEGHadronizer") << "Start production";
83 
85 
86  thepegEvent = eg_->shoot();
87  if (!thepegEvent) {
88  edm::LogWarning("Generator|ThePEGHadronizer") << "thepegEvent not initialized";
89  return false;
90  }
91 
93  if (!event().get()) {
94  edm::LogWarning("Generator|ThePEGHadronizer") << "genEvent not initialized";
95  return false;
96  }
97 
98  return true;
99 }
100 
102 {
103  HepMC::PdfInfo pdf;
104  clearAuxiliary(event().get(), &pdf);
105  fillAuxiliary(event().get(), &pdf, thepegEvent);
106  event()->set_pdf_info(pdf);
107 
108  eventInfo().reset(new GenEventInfoProduct(event().get()));
109  eventInfo()->setBinningValues(
110  std::vector<double>(1, pthat(thepegEvent)));
111 
112  if (eventsToPrint) {
113  eventsToPrint--;
114  event()->print();
115  }
116 
117  if (iobc_.get())
118  iobc_->write_event(event().get());
119 
120  edm::LogInfo("Generator|ThePEGHadronizer") << "Event produced";
121 }
122 
124 {
125  return true;
126 }
127 
129 {
130  return true;
131 }
132 
134 
137 
138 #if 0
140 DEFINE_FWK_MODULE(ThePEGHadronizerFilter);
141 #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()
const char * classname() const
bool declareSpecialSettings(const std::vector< std::string >)
static void clearAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf)
static void fillAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf, const ThePEG::EventPtr &event)