test
CMS 3D CMS Logo

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