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>
14 
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 
56  std::unique_ptr<GenLumiInfoHeader> getGenLumiInfoHeader() const override;
57  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&);
58  void randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine);
59 
60 
61  private:
62 
63  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }
64 
65  unsigned int eventsToPrint;
66 
67  ThePEG::EventPtr thepegEvent;
68 
69  std::shared_ptr<lhef::LHEProxy> proxy_;
73 
74  unsigned int firstLumiBlock=0;
75  unsigned int currentLumiBlock=0;
76 };
77 
79  Herwig7Interface(pset),
80  BaseHadronizer(pset),
81  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
82  handlerDirectory_(pset.getParameter<std::string>("eventHandlers")),
83  runFileName(pset.getParameter<std::string>("run"))
84 {
85  initRepository(pset);
87 }
88 
90 {
91 }
92 
94 {
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  std::ifstream runFile(runFileName+".run");
104  if (runFile.fail()) //required for showering of LHE files
105  {
107  }
108  if (!initGenerator())
109  {
110  edm::LogInfo("Generator|Herwig7Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
111  exit(0);
112  }
113  }
114  return true;
115 }
116 
118 {
119  edm::LogError("Herwig7 interface") << "Read in of LHE files is not supported in this way. You can read them manually if necessary.";
120  return false;
121 }
122 
123 bool Herwig7Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
124 {
125  return false;
126 }
127 
129 {
130  if(eg_){
132  eg_->integratedXSec() / ThePEG::picobarn,
133  eg_->integratedXSecErr() / ThePEG::picobarn));
134  }
135 }
136 
138 {
139  edm::LogInfo("Generator|Herwig7Hadronizer") << "Start production";
140 
141  try {
142  thepegEvent = eg_->shoot();
143  } catch (std::exception& exc) {
144  edm::LogWarning("Generator|Herwig7Hadronizer") << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
145  return false;
146  }
147 
148  if (!thepegEvent) {
149  edm::LogWarning("Generator|Herwig7Hadronizer") << "thepegEvent not initialized";
150  return false;
151  }
152 
154  if (!event().get()) {
155  edm::LogWarning("Generator|Herwig7Hadronizer") << "genEvent not initialized";
156  return false;
157  }
158 
159  return true;
160 }
161 
163 {
164 
165  edm::LogError("Herwig7 interface") << "Read in of LHE files is not supported in this way. You can read them manually if necessary.";
166  return false;
167 }
168 
169 std::unique_ptr<GenLumiInfoHeader> Herwig7Hadronizer::getGenLumiInfoHeader() const {
170  auto genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
171 
172  if (thepegEvent)
173  {
174  int weights_number = thepegEvent->optionalWeights().size();
175 
176  if(weights_number > 1){
177  genLumiInfoHeader->weightNames().reserve(weights_number + 1);
178  genLumiInfoHeader->weightNames().push_back("nominal");
179  std::map<std::string,double> weights_map = thepegEvent->optionalWeights();
180  for (std::map<std::string,double>::iterator it = weights_map.begin(); it != weights_map.end(); it++)
181  {
182  genLumiInfoHeader->weightNames().push_back(it->first);
183  }
184  }
185  }
186 
187  return genLumiInfoHeader;
188 }
189 
190 void Herwig7Hadronizer::randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine)
191 {
192  BaseHadronizer::randomizeIndex(lumi, rengine);
193 
194  if (firstLumiBlock==0)
195  {
197  }
199 
200 }
201 
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|Herwig7Hadronizer") << "Event produced";
217 }
218 
220 {
221  return true;
222 }
223 
225 {
226  return true;
227 }
228 
230 
233 
LuminosityBlockID id() const
~Herwig7Hadronizer() override
static double pthat(const ThePEG::EventPtr &event)
void initRepository(const edm::ParameterSet &params)
unsigned int firstLumiBlock
std::unique_ptr< HepMC::IO_BaseClass > iobc_
void setInternalXSec(const XSec &xsec)
GenRunInfoProduct & runInfo()
ThePEG::EGPtr eg_
bool declareSpecialSettings(const std::vector< std::string >)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::GeneratorFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7GeneratorFilter
unsigned int eventsToPrint
static std::unique_ptr< HepMC::GenEvent > convert(const ThePEG::EventPtr &event)
edm::ParameterSet paramSettings
edm::HadronizerFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7HadronizerFilter
bool initializeForInternalPartons()
std::unique_ptr< HepMC::GenEvent > & event()
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::shared_ptr< lhef::LHEProxy > proxy_
std::unique_ptr< GenEventInfoProduct > & eventInfo()
LuminosityBlockNumber_t luminosityBlock() const
unsigned int currentLumiBlock
Herwig7Hadronizer(const edm::ParameterSet &params)
ThePEG::EventPtr thepegEvent
const char * classname() const
std::unique_ptr< GenLumiInfoHeader > getGenLumiInfoHeader() const override
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)