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 
29 #include <Herwig/API/HerwigAPI.h>
30 #include "CLHEP/Random/RandomEngine.h"
31 
32 namespace CLHEP {
33  class HepRandomEngine;
34 }
35 
37 public:
39  ~Herwig7Hadronizer() override;
40 
41  bool readSettings(int) { return true; }
44  bool declareStableParticles(const std::vector<int>& pdgIds);
45  bool declareSpecialSettings(const std::vector<std::string>) { return true; }
46 
47  void statistics();
48 
50  bool hadronize();
51  bool decay();
52  bool residualDecay();
53  void finalizeEvent();
54 
55  const char* classname() const { return "Herwig7Hadronizer"; }
56  std::unique_ptr<GenLumiInfoHeader> getGenLumiInfoHeader() const override;
58  void randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine);
59 
60 private:
61  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }
62 
63  unsigned int eventsToPrint;
64 
65  ThePEG::EventPtr thepegEvent;
66  bool haveEvt = false;
67 
68  std::shared_ptr<lhef::LHEProxy> proxy_;
72 
73  unsigned int firstLumiBlock = 0;
74  unsigned int currentLumiBlock = 0;
75 };
76 
79  BaseHadronizer(pset),
80  eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
81  handlerDirectory_(pset.getParameter<std::string>("eventHandlers")),
82  runFileName(pset.getParameter<std::string>("run")) {
85 }
86 
88 
91  std::ifstream runFile(runFileName + ".run");
92  if (runFile.fail()) //required for showering of LHE files
93  {
95  }
96  if (!initGenerator()) {
97  edm::LogInfo("Generator|Herwig7Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
98  exit(0);
99  }
100  }
101  return true;
102 }
103 
106  std::ifstream runFile(runFileName + ".run");
107  if (runFile.fail()) //required for showering of LHE files
108  {
110  }
111  if (!initGenerator()) {
112  edm::LogInfo("Generator|Herwig7Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
113  exit(0);
114  }
115  }
116  return true;
117 }
118 
119 bool Herwig7Hadronizer::declareStableParticles(const std::vector<int>& pdgIds) { return false; }
120 
122  if (eg_) {
124  GenRunInfoProduct::XSec(eg_->integratedXSec() / ThePEG::picobarn, eg_->integratedXSecErr() / ThePEG::picobarn));
125  }
126 }
127 
129  edm::LogInfo("Generator|Herwig7Hadronizer") << "Start production";
130 
131  try {
132  thepegEvent = eg_->shoot();
133  } catch (std::exception& exc) {
134  edm::LogWarning("Generator|Herwig7Hadronizer")
135  << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
136  return false;
137  }
138 
139  if (!thepegEvent) {
140  edm::LogWarning("Generator|Herwig7Hadronizer") << "thepegEvent not initialized";
141  return false;
142  }
143 
145  if (!event().get()) {
146  edm::LogWarning("Generator|Herwig7Hadronizer") << "genEvent not initialized";
147  return false;
148  }
149 
150  return true;
151 }
152 
154  if (!haveEvt) {
155  try {
156  thepegEvent = eg_->shoot();
157  haveEvt = true;
158  } catch (std::exception& exc) {
159  edm::LogWarning("Generator|Herwig7Hadronizer")
160  << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
161  return false;
162  }
163  }
164  int evtnum = lheEvent()->evtnum();
165  if (evtnum == -1) {
166  edm::LogError("Generator|Herwig7Hadronizer")
167  << "Event number not set in lhe file, needed for correctly aligning Herwig and LHE events!";
168  return false;
169  }
170  if (thepegEvent->number() < evtnum) {
171  edm::LogError("Herwig7 interface") << "Herwig does not seem to be generating events in order, did you set "
172  "/Herwig/EventHandlers/FxFxLHReader:AllowedToReOpen Yes?";
173  return false;
174  } else if (thepegEvent->number() == evtnum) {
175  haveEvt = false;
176  if (!thepegEvent) {
177  edm::LogWarning("Generator|Herwig7Hadronizer") << "thepegEvent not initialized";
178  return false;
179  }
180 
182  if (!event().get()) {
183  edm::LogWarning("Generator|Herwig7Hadronizer") << "genEvent not initialized";
184  return false;
185  }
186  return true;
187  }
188  edm::LogWarning("Generator|Herwig7Hadronizer") << "Event " << evtnum << " not generated (likely skipped in merging)";
189  return false;
190 }
191 
193  eventInfo() = std::make_unique<GenEventInfoProduct>(event().get());
194  eventInfo()->setBinningValues(std::vector<double>(1, pthat(thepegEvent)));
195 
196  if (eventsToPrint) {
197  eventsToPrint--;
198  event()->print();
199  }
200 
201  if (iobc_.get())
202  iobc_->write_event(event().get());
203 
204  edm::LogInfo("Generator|Herwig7Hadronizer") << "Event produced";
205 }
206 
207 bool Herwig7Hadronizer::decay() { return true; }
208 
209 bool Herwig7Hadronizer::residualDecay() { return true; }
210 
211 std::unique_ptr<GenLumiInfoHeader> Herwig7Hadronizer::getGenLumiInfoHeader() const {
212  auto genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
213 
214  if (thepegEvent) {
215  int weights_number = thepegEvent->optionalWeights().size();
216 
217  if (weights_number > 1) {
218  genLumiInfoHeader->weightNames().reserve(weights_number + 1);
219  genLumiInfoHeader->weightNames().push_back("nominal");
220  std::map<std::string, double> weights_map = thepegEvent->optionalWeights();
221  for (std::map<std::string, double>::iterator it = weights_map.begin(); it != weights_map.end(); it++) {
222  genLumiInfoHeader->weightNames().push_back(it->first);
223  }
224  }
225  }
226 
227  return genLumiInfoHeader;
228 }
229 
230 void Herwig7Hadronizer::randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine) {
231  BaseHadronizer::randomizeIndex(lumi, rengine);
232 
233  if (firstLumiBlock == 0) {
234  firstLumiBlock = lumi.id().luminosityBlock();
235  }
236  currentLumiBlock = lumi.id().luminosityBlock();
237 }
238 
240 
243 
~Herwig7Hadronizer() override
static double pthat(const ThePEG::EventPtr &event)
void initRepository(const edm::ParameterSet &params)
unsigned int firstLumiBlock
void setPEGRandomEngine(CLHEP::HepRandomEngine *)
std::unique_ptr< HepMC::IO_BaseClass > iobc_
Log< level::Error, false > LogError
void setInternalXSec(const XSec &xsec)
GenRunInfoProduct & runInfo()
ThePEG::EGPtr eg_
bool declareSpecialSettings(const std::vector< std::string >)
const char * classname() const
lhef::LHEEvent * lheEvent()
edm::GeneratorFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7GeneratorFilter
unsigned int eventsToPrint
static std::unique_ptr< HepMC::GenEvent > convert(const ThePEG::EventPtr &event)
std::unique_ptr< GenLumiInfoHeader > getGenLumiInfoHeader() const override
edm::ParameterSet paramSettings
edm::HadronizerFilter< Herwig7Hadronizer, gen::ExternalDecayDriver > Herwig7HadronizerFilter
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool initializeForInternalPartons()
std::unique_ptr< HepMC::GenEvent > & event()
int evtnum() const
Definition: LHEEvent.h:55
Log< level::Info, false > LogInfo
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::shared_ptr< lhef::LHEProxy > proxy_
std::unique_ptr< GenEventInfoProduct > & eventInfo()
unsigned int currentLumiBlock
Herwig7Hadronizer(const edm::ParameterSet &params)
ThePEG::EventPtr thepegEvent
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
const std::string handlerDirectory_
const std::string runFileName
bool declareStableParticles(const std::vector< int > &pdgIds)
Log< level::Warning, false > LogWarning
Definition: event.py:1
void randomizeIndex(edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine)
def exit(msg="")