CMS 3D CMS Logo

SherpaHadronizer.cc
Go to the documentation of this file.
1 #include <cassert>
2 #include <iostream>
3 #include <sstream>
4 #include <string>
5 #include <memory>
6 #include <cstdint>
7 #include <vector>
8 
9 #include "SHERPA/Main/Sherpa.H"
10 #include "ATOOLS/Math/Random.H"
11 
12 #include "ATOOLS/Org/Run_Parameter.H"
13 #include "ATOOLS/Org/MyStrStream.H"
14 #include "ATOOLS/Org/CXXFLAGS.H"
15 #include "ATOOLS/Org/CXXFLAGS_PACKAGES.H"
16 #include "ATOOLS/Org/My_MPI.H"
17 
23 
24 #include "CLHEP/Random/RandomEngine.h"
25 
27 
28 //This unnamed namespace is used (instead of static variables) to pass the
29 //randomEngine passed to doSetRandomEngine to the External Random
30 //Number Generator CMS_SHERPA_RNG of sherpa
31 //The advantage of the unnamed namespace over static variables is
32 //that it is only accessible in this file
33 
34 namespace {
35  CLHEP::HepRandomEngine *ExternalEngine = nullptr;
36  CLHEP::HepRandomEngine *GetExternalEngine() { return ExternalEngine; }
37  void SetExternalEngine(CLHEP::HepRandomEngine *v) { ExternalEngine = v; }
38 } // namespace
39 
41 public:
43  ~SherpaHadronizer() override;
44 
45  bool readSettings(int) { return true; }
47  bool declareStableParticles(const std::vector<int> &pdgIds);
48  bool declareSpecialSettings(const std::vector<std::string> &) { return true; }
49  void statistics();
51  bool decay();
53  bool residualDecay();
54  void finalizeEvent();
55  std::unique_ptr<GenLumiInfoHeader> getGenLumiInfoHeader() const override;
56  const char *classname() const { return "SherpaHadronizer"; }
57 
58 private:
59  void doSetRandomEngine(CLHEP::HepRandomEngine *v) override;
60 
68  unsigned int maxEventsToPrint;
69  std::vector<std::string> arguments;
73  std::vector<std::string> weightlist;
74  std::vector<std::string> variationweightlist;
75 };
76 
77 class CMS_SHERPA_RNG : public ATOOLS::External_RNG {
78 public:
80  edm::LogVerbatim("SherpaHadronizer") << "Use stored reference for the external RNG";
81  setRandomEngine(GetExternalEngine());
82  }
83  void setRandomEngine(CLHEP::HepRandomEngine *v) { randomEngine = v; }
84 
85 private:
86  double Get() override;
87  CLHEP::HepRandomEngine *randomEngine;
88 };
89 
90 void SherpaHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine *v) {
91  CMS_SHERPA_RNG *cmsSherpaRng = dynamic_cast<CMS_SHERPA_RNG *>(ATOOLS::ran->GetExternalRng());
92  if (cmsSherpaRng == nullptr) {
93  //First time call to this function makes the interface store the reference in the unnamed namespace
94  if (!isRNGinitialized) {
95  isRNGinitialized = true;
96  edm::LogVerbatim("SherpaHadronizer") << "Store assigned reference of the randomEngine";
97  SetExternalEngine(v);
98  // Throw exception if there is no reference to an external RNG and it is not the first call!
99  } else {
100  if (isInitialized and v != nullptr) {
102  << "The Sherpa interface got a randomEngine reference but there is "
103  "no reference to the external RNG to hand it over to\n";
104  }
105  }
106  } else {
107  cmsSherpaRng->setRandomEngine(v);
108  }
109 }
110 
112  : BaseHadronizer(params),
113  SherpaParameterSet(params.getParameter<edm::ParameterSet>("SherpaParameters")),
114  isRNGinitialized(false) {
115  if (!params.exists("SherpaProcess"))
116  SherpaProcess = "";
117  else
118  SherpaProcess = params.getParameter<std::string>("SherpaProcess");
119  if (!params.exists("SherpaPath"))
120  SherpaPath = "";
121  else
122  SherpaPath = params.getParameter<std::string>("SherpaPath");
123  if (!params.exists("SherpaPathPiece"))
124  SherpaPathPiece = "";
125  else
126  SherpaPathPiece = params.getParameter<std::string>("SherpaPathPiece");
127  if (!params.exists("SherpaResultDir"))
128  SherpaResultDir = "Result";
129  else
130  SherpaResultDir = params.getParameter<std::string>("SherpaResultDir");
131  if (!params.exists("SherpaDefaultWeight"))
132  SherpaDefaultWeight = 1.;
133  else
134  SherpaDefaultWeight = params.getParameter<double>("SherpaDefaultWeight");
135  if (!params.exists("maxEventsToPrint"))
136  maxEventsToPrint = 0;
137  else
138  maxEventsToPrint = params.getParameter<int>("maxEventsToPrint");
139  // if hepmcextendedweights is used the event weights have to be reordered ( unordered list can be accessed via event->weights().write() )
140  // two lists have to be provided:
141  // 1) SherpaWeights
142  // - containing nominal event weight, combined matrix element and phase space weight, event normalization, and possibly other sherpa weights
143  // 2) SherpaVariationsWeights
144  // - containing weights from scale and PDF variations ( have to be defined in the runcard )
145  // - in case of unweighted events these weights are also divided by the event normalization (see list 1 )
146  // Sherpa Documentation: http://sherpa.hepforge.org/doc/SHERPA-MC-2.2.0.html#Scale-and-PDF-variations
147  if (!params.exists("SherpaWeightsBlock")) {
148  rearrangeWeights = false;
149  } else {
150  rearrangeWeights = true;
151  edm::ParameterSet WeightsBlock = params.getParameter<edm::ParameterSet>("SherpaWeightsBlock");
152  if (WeightsBlock.exists("SherpaWeights"))
153  weightlist = WeightsBlock.getParameter<std::vector<std::string> >("SherpaWeights");
154  else
155  throw cms::Exception("SherpaInterface") << "SherpaWeights does not exists in SherpaWeightsBlock" << std::endl;
156  if (WeightsBlock.exists("SherpaVariationWeights"))
157  variationweightlist = WeightsBlock.getParameter<std::vector<std::string> >("SherpaVariationWeights");
158  else
159  throw cms::Exception("SherpaInterface")
160  << "SherpaVariationWeights does not exists in SherpaWeightsBlock" << std::endl;
161  edm::LogVerbatim("SherpaHadronizer") << "SherpaHadronizer will try rearrange the event weights according to "
162  "SherpaWeights and SherpaVariationWeights";
163  }
164 
165  spf::SherpackFetcher Fetcher(params);
166  int retval = Fetcher.Fetch();
167  if (retval != 0) {
168  throw cms::Exception("SherpaInterface") << "SherpaHadronizer: Preparation of Sherpack failed ... ";
169  }
170  // The ids (names) of parameter sets to be read (Analysis,Run) to create Analysis.dat, Run.dat
171  //They are given as a vstring.
172  std::vector<std::string> setNames = SherpaParameterSet.getParameter<std::vector<std::string> >("parameterSets");
173  //Loop all set names...
174  for (unsigned i = 0; i < setNames.size(); ++i) {
175  // ...and read the parameters for each set given in vstrings
176  std::vector<std::string> pars = SherpaParameterSet.getParameter<std::vector<std::string> >(setNames[i]);
177  edm::LogVerbatim("SherpaHadronizer") << "Write Sherpa parameter set " << setNames[i] << " to " << setNames[i]
178  << ".dat ";
179  std::string datfile = SherpaPath + "/" + setNames[i] + ".dat";
180  std::ofstream os(datfile.c_str());
181  // Loop over all strings and write the according *.dat
182  for (std::vector<std::string>::const_iterator itPar = pars.begin(); itPar != pars.end(); ++itPar) {
183  os << (*itPar) << std::endl;
184  }
185  }
186 
187  //To be conform to the default Sherpa usage create a command line:
188  //name of executable (only for demonstration, could also be empty)
189  std::string shRun = "./Sherpa";
190  //Path where the Sherpa libraries are stored
191  std::string shPath = "PATH=" + SherpaPath;
192  // new for Sherpa 1.3.0, suggested by authors
193  std::string shPathPiece = "PATH_PIECE=" + SherpaPathPiece;
194  //Path where results are stored
195  std::string shRes = "RESULT_DIRECTORY=" + SherpaResultDir; // from Sherpa 1.2.0 on
196  //Name of the external random number class
197  std::string shRng = "EXTERNAL_RNG=CMS_SHERPA_RNG";
198 
199  //create the command line
200  arguments.push_back(shRun);
201  arguments.push_back(shPath);
202  arguments.push_back(shPathPiece);
203  arguments.push_back(shRes);
204  arguments.push_back(shRng);
205  isInitialized = false;
206  //initialization of Sherpa moved to initializeForInternalPartons
207 #ifdef USING__MPI
208  MPI::Init();
209 #endif
210 }
211 
213  Generator->~Sherpa();
214 #ifdef USING__MPI
215  MPI::Finalize();
216 #endif
217 }
218 
220  //initialize Sherpa but only once
221  if (!isInitialized) {
222  int argc = arguments.size();
223  char *argv[argc];
224  for (int l = 0; l < argc; l++)
225  argv[l] = (char *)arguments[l].c_str();
226  Generator->InitializeTheRun(argc, argv);
227  Generator->InitializeTheEventHandler();
228  isInitialized = true;
229  }
230  return true;
231 }
232 
233 bool SherpaHadronizer::declareStableParticles(const std::vector<int> &pdgIds) { return false; }
234 
236  //calculate statistics
237  Generator->SummarizeRun();
238 
239  //get the xsec & err
240  double xsec_val = Generator->TotalXS();
241  double xsec_err = Generator->TotalErr();
242 
243  //set the internal cross section in pb in GenRunInfoProduct
244  runInfo().setInternalXSec(GenRunInfoProduct::XSec(xsec_val, xsec_err));
245 }
246 
248  //get the next event and check if it produced
249  bool rc = false;
250  int itry = 0;
251  bool gen_event = true;
252  while ((itry < 3) && gen_event) {
253  try {
254  rc = Generator->GenerateOneEvent();
255  gen_event = false;
256  } catch (...) {
257  ++itry;
258  std::cerr << "Exception from Generator->GenerateOneEvent() catch. Call # " << itry << " for this event\n";
259  }
260  }
261  if (rc) {
262  //convert it to HepMC2
263  auto evt = std::make_unique<HepMC::GenEvent>();
264  Generator->FillHepMCEvent(*evt);
265 
266  // in case of unweighted events sherpa puts the max weight as event weight.
267  // this is not optimal, we want 1 for unweighted events, so we check
268  // whether we are producing unweighted events ("EVENT_GENERATION_MODE" == "1")
269  // the information about the weights to the HepMC weight vector:
270  // [0] event weight
271  // [1] combined matrix element and phase space weight (missing only PDF information, thus directly suitable for PDF reweighting)
272  // [2] event weight normalisation (in case of unweighted events event weights of ~ +/-1 can be obtained by (event weight)/(event weight normalisation))
273  // [3] number of trials.
274  // see also: https://sherpa.hepforge.org/doc/SHERPA-MC-2.1.0.html#Event-output-formats
275  bool unweighted = false;
276  double weight_normalization = -1;
277  if (ATOOLS::ToType<int>(ATOOLS::rpa->gen.Variable("EVENT_GENERATION_MODE")) == 1) {
278  if (evt->weights().size() > 2) {
279  unweighted = true;
280  weight_normalization = evt->weights()[2];
281  } else {
282  throw cms::Exception("SherpaInterface")
283  << "Requested unweighted production. Missing normalization weight." << std::endl;
284  }
285  }
286 
287  // vector to fill new weights in correct order
288  std::vector<double> newWeights;
289  if (rearrangeWeights) {
290  for (auto &i : weightlist) {
291  if (evt->weights().has_key(i)) {
292  newWeights.push_back(evt->weights()[i]);
293  } else {
294  throw cms::Exception("SherpaInterface")
295  << "Missing weights! Key " << i << " not found, please check the weight definition!" << std::endl;
296  }
297  }
298  for (auto &i : variationweightlist) {
299  if (evt->weights().has_key(i)) {
300  if (unweighted) {
301  newWeights.push_back(evt->weights()[i] / weight_normalization);
302  } else {
303  newWeights.push_back(evt->weights()[i]);
304  }
305  } else {
306  throw cms::Exception("SherpaInterface")
307  << "Missing weights! Key " << i << " not found, please check the weight definition!" << std::endl;
308  }
309  }
310 
311  //Change original weights for reordered ones
312  evt->weights().clear();
313  for (auto &elem : newWeights) {
314  evt->weights().push_back(elem);
315  }
316  }
317 
318  if (unweighted) {
319  evt->weights()[0] /= weight_normalization;
320  }
321  resetEvent(std::move(evt));
322  return true;
323  } else {
324  return false;
325  }
326 }
327 
328 bool SherpaHadronizer::decay() { return true; }
329 
330 bool SherpaHadronizer::residualDecay() { return true; }
331 
333  //******** Verbosity *******
334  if (maxEventsToPrint > 0) {
336  event()->print();
337  }
338 }
339 
340 //GETTER for the external random numbers
341 DECLARE_GETTER(CMS_SHERPA_RNG, "CMS_SHERPA_RNG", ATOOLS::External_RNG, ATOOLS::RNG_Key);
342 
343 ATOOLS::External_RNG *ATOOLS::Getter<ATOOLS::External_RNG, ATOOLS::RNG_Key, CMS_SHERPA_RNG>::operator()(
344  const ATOOLS::RNG_Key &) const {
345  return new CMS_SHERPA_RNG();
346 }
347 
348 void ATOOLS::Getter<ATOOLS::External_RNG, ATOOLS::RNG_Key, CMS_SHERPA_RNG>::PrintInfo(std::ostream &str,
349  const size_t) const {
350  str << "CMS_SHERPA_RNG interface";
351 }
352 
354  if (randomEngine == nullptr) {
355  throw edm::Exception(edm::errors::LogicError) << "The Sherpa code attempted to a generate random number while\n"
356  << "the engine pointer was null. This might mean that the code\n"
357  << "was modified to generate a random number outside the event and\n"
358  << "beginLuminosityBlock methods, which is not allowed.\n";
359  }
360  return randomEngine->flat();
361 }
362 std::unique_ptr<GenLumiInfoHeader> SherpaHadronizer::getGenLumiInfoHeader() const {
363  auto genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
364 
365  if (rearrangeWeights) {
366  edm::LogPrint("SherpaHadronizer") << "The order of event weights was changed!";
367  for (auto &i : weightlist) {
368  genLumiInfoHeader->weightNames().push_back(i);
369  edm::LogVerbatim("SherpaHadronizer") << i;
370  }
371  for (auto &i : variationweightlist) {
372  genLumiInfoHeader->weightNames().push_back(i);
373  edm::LogVerbatim("SherpaHadronizer") << i;
374  }
375  }
376  return genLumiInfoHeader;
377 }
379 
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
CLHEP::HepRandomEngine * randomEngine
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< std::string > weightlist
std::string SherpaPathPiece
unsigned int maxEventsToPrint
edm::GeneratorFilter< SherpaHadronizer, gen::ExternalDecayDriver > SherpaGeneratorFilter
void setInternalXSec(const XSec &xsec)
DECLARE_GETTER(CMS_SHERPA_RNG, "CMS_SHERPA_RNG", ATOOLS::External_RNG, ATOOLS::RNG_Key)
bool declareStableParticles(const std::vector< int > &pdgIds)
GenRunInfoProduct & runInfo()
edm::ParameterSet SherpaParameterSet
bool initializeForInternalPartons()
void resetEvent(std::unique_ptr< HepMC::GenEvent > event)
std::unique_ptr< GenLumiInfoHeader > getGenLumiInfoHeader() const override
std::vector< std::string > variationweightlist
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::unique_ptr< HepMC::GenEvent > & event()
Log< level::Warning, true > LogPrint
bool generatePartonsAndHadronize()
std::string SherpaPath
~SherpaHadronizer() override
bool declareSpecialSettings(const std::vector< std::string > &)
void setRandomEngine(CLHEP::HepRandomEngine *v)
const char * classname() const
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::string SherpaProcess
std::string SherpaResultDir
HLT enums.
std::vector< std::string > arguments
SherpaHadronizer(const edm::ParameterSet &params)
#define str(s)
def move(src, dest)
Definition: eostools.py:511
std::string SherpaChecksum
double Get() override