CMS 3D CMS Logo

Herwig7Hadronizer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <sstream>
3 #include <fstream>
4 
5 #include <HepMC/GenEvent.h>
6 #include <HepMC/IO_BaseClass.h>
7 
8 #include <ThePEG/Repository/Repository.h>
9 #include <ThePEG/EventRecord/Event.h>
10 #include <ThePEG/Config/ThePEG.h>
11 #include <ThePEG/LesHouches/LesHouchesReader.h>
12 
16 
20 
24 
26 
28 #include <Herwig/API/HerwigAPI.h>
29 #include "CLHEP/Random/RandomEngine.h"
30 
31 namespace CLHEP {
32  class HepRandomEngine;
33 }
34 
36  public:
37  Herwig7Hadronizer(const edm::ParameterSet &params);
38  ~Herwig7Hadronizer() override;
39 
40  bool readSettings( int ) { return true; }
41  bool initializeForInternalPartons();
42  bool initializeForExternalPartons();
43  bool declareStableParticles(const std::vector<int> &pdgIds);
44  bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
45 
46  void statistics();
47 
48  bool generatePartonsAndHadronize();
49  bool hadronize();
50  bool decay();
51  bool residualDecay();
52  void finalizeEvent();
53 
54  const char *classname() const { return "Herwig7Hadronizer"; }
55  GenLumiInfoHeader *getGenLumiInfoHeader() const override;
56  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&);
57  void randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine);
58 
59  private:
60 
61  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }
62 
63  unsigned int eventsToPrint;
64 
65  ThePEG::EventPtr thepegEvent;
66 
67  std::shared_ptr<lhef::LHEProxy> proxy_;
71  unsigned int firstLumiBlock=0;
72  unsigned int currentLumiBlock=0;
73 };
74 
76  Herwig7Interface(pset),
77  BaseHadronizer(pset),
78  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
79  handlerDirectory_(pset.getParameter<std::string>("eventHandlers")),
80  runFileName(pset.getParameter<std::string>("run"))
81 {
82  initRepository(pset);
84 }
85 
87 {
88 }
89 
91 {
93  {
94  std::ifstream runFile(runFileName+".run");
95  if (runFile.fail()) //required for showering of LHE files
96  {
98  }
99  if (!initGenerator())
100  {
101  edm::LogInfo("Generator|Herwig7Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
102  exit(0);
103  }
104  }
105  return true;
106 }
107 
109 {
110  edm::LogError("Herwig7 interface") << "Read in of LHE files is not supported in this way. You can read them manually if necessary.";
111  return false;
112 }
113 
114 bool Herwig7Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
115 {
116  return false;
117 }
118 
120 {
121  if(eg_){
123  eg_->integratedXSec() / ThePEG::picobarn,
124  eg_->integratedXSecErr() / ThePEG::picobarn));
125  }
126 }
127 
129 {
130  edm::LogInfo("Generator|Herwig7Hadronizer") << "Start production";
131  try {
132  thepegEvent = eg_->shoot();
133  } catch (std::exception& exc) {
134  edm::LogWarning("Generator|Herwig7Hadronizer") << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
135  return false;
136  }
137 
138  if (!thepegEvent) {
139  edm::LogWarning("Generator|Herwig7Hadronizer") << "thepegEvent not initialized";
140  return false;
141  }
142 
144  if (!event().get()) {
145  edm::LogWarning("Generator|Herwig7Hadronizer") << "genEvent not initialized";
146  return false;
147  }
148 
149  return true;
150 }
151 
153 {
154 
155  edm::LogError("Herwig7 interface") << "Read in of LHE files is not supported in this way. You can read them manually if necessary.";
156  return false;
157 }
158 
160 {
161  eventInfo().reset(new GenEventInfoProduct(event().get()));
162  eventInfo()->setBinningValues(
163  std::vector<double>(1, pthat(thepegEvent)));
164 
165  if (eventsToPrint) {
166  eventsToPrint--;
167  event()->print();
168  }
169 
170  if (iobc_.get())
171  iobc_->write_event(event().get());
172 
173  edm::LogInfo("Generator|Herwig7Hadronizer") << "Event produced";
174 }
175 
177 {
178  return true;
179 }
180 
182 {
183  return true;
184 }
185 
187  GenLumiInfoHeader *genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
188 
189  if (thepegEvent)
190  {
191  int weights_number = thepegEvent->optionalWeights().size();
192 
193  if(weights_number > 1){
194  genLumiInfoHeader->weightNames().reserve(weights_number + 1);
195  genLumiInfoHeader->weightNames().push_back("nominal");
196  std::map<std::string,double> weights_map = thepegEvent->optionalWeights();
197  for (std::map<std::string,double>::iterator it = weights_map.begin(); it != weights_map.end(); it++)
198  {
199  genLumiInfoHeader->weightNames().push_back(it->first);
200  }
201  }
202  }
203 
204  return genLumiInfoHeader;
205 }
206 
207 void Herwig7Hadronizer::randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine)
208 {
209  BaseHadronizer::randomizeIndex(lumi, rengine);
210 
211  if (firstLumiBlock==0)
212  {
214  }
216 
217 }
218 
219 
221 
224 
LuminosityBlockID id() const
~Herwig7Hadronizer() override
static double pthat(const ThePEG::EventPtr &event)
void initRepository(const edm::ParameterSet &params)
GenLumiInfoHeader * getGenLumiInfoHeader() const override
unsigned int firstLumiBlock
std::auto_ptr< HepMC::IO_BaseClass > iobc_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const std::vector< std::string > & weightNames() const
genLumiInfoHeader
Definition: nano_cff.py:97
void setInternalXSec(const XSec &xsec)
std::auto_ptr< HepMC::GenEvent > & event()
GenRunInfoProduct & runInfo()
ThePEG::EGPtr eg_
bool declareSpecialSettings(const std::vector< std::string >)
edm::GeneratorFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7GeneratorFilter
unsigned int eventsToPrint
edm::ParameterSet paramSettings
edm::HadronizerFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7HadronizerFilter
static std::auto_ptr< HepMC::GenEvent > convert(const ThePEG::EventPtr &event)
bool initializeForInternalPartons()
std::auto_ptr< GenEventInfoProduct > & eventInfo()
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::shared_ptr< lhef::LHEProxy > proxy_
LuminosityBlockNumber_t luminosityBlock() const
unsigned int currentLumiBlock
Herwig7Hadronizer(const edm::ParameterSet &params)
ThePEG::EventPtr thepegEvent
const char * classname() const
const std::string handlerDirectory_
const std::string runFileName
bool declareStableParticles(const std::vector< int > &pdgIds)
Definition: event.py:1
void randomizeIndex(edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine)