CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Herwig6Hadronizer.cc
Go to the documentation of this file.
1 #include <cstring>
2 #include <sstream>
3 #include <string>
4 #include <vector>
5 #include <memory>
6 #include <map>
7 #include <set>
8 
9 #include <boost/shared_ptr.hpp>
10 #include <boost/algorithm/string/classification.hpp>
11 #include <boost/algorithm/string/split.hpp>
12 
13 #include <HepMC/GenEvent.h>
14 #include <HepMC/GenParticle.h>
15 #include <HepMC/GenVertex.h>
16 #include <HepMC/PdfInfo.h>
17 #include <HepMC/HerwigWrapper6_4.h>
18 #include <HepMC/HEPEVT_Wrapper.h>
19 #include <HepMC/IO_HERWIG.h>
20 #include "HepPID/ParticleIDTranslations.hh"
21 
23 
26 
29 
34 
36 
39 
41 
42 extern "C" {
43  void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8]);
44  double hwualf_(int *mode, double* scale);
45  double hwuaem_(double* scale);
46 }
47 
48 // helpers
49 namespace {
50  static inline bool call_hwmsct()
51  {
52  int result;
53  hwmsct(&result);
54  return result;
55  }
56 
57  static int pdgToHerwig(int ipdg, char *nwig)
58  {
59  int iopt = 1;
60  int iwig = 0;
61  hwuidt_(&iopt, &ipdg, &iwig, nwig);
62  return ipdg ? iwig : 0;
63  }
64 
65  static bool markStable(int pdgId)
66  {
67  char nwig[9] = " ";
68  if (!pdgToHerwig(pdgId, nwig))
69  return false;
70  hwusta(nwig, 1);
71  return true;
72  }
73 }
74 
76  public gen::Herwig6Instance {
77  public:
78  Herwig6Hadronizer(const edm::ParameterSet &params);
80 
81  void setSLHAFromHeader(const std::vector<std::string> &lines);
82  bool initialize(const lhef::HEPRUP *heprup);
83 
85  bool initializeForExternalPartons() { return initialize(lheRunInfo()->getHEPRUP()); }
86 
87  bool declareStableParticles( const std::vector<int> &pdgIds);
88  bool declareSpecialSettings( const std::vector<std::string> );
89 
90  void statistics();
91 
92 
94  bool hadronize();
95  bool decay();
96  bool residualDecay();
97  void finalizeEvent();
98 
99  const char *classname() const { return "Herwig6Hadronizer"; }
100 
101  private:
102  void clear();
103 
104  int pythiaStatusCode(const HepMC::GenParticle *p) const;
105  void pythiaStatusCodes();
106 
107  virtual void upInit();
108  virtual void upEvnt();
109 
110  HepMC::IO_HERWIG conv;
111  bool needClear;
113 
120  double comEnergy;
121  bool useJimmy;
124  bool fConvertToPDG; // convert PIDs
125 
127 
128 };
129 
130 extern "C" {
131  void hwwarn_(const char *method, int *id);
132  void setherwpdf_(void);
133  void mysetpdfpath_(const char *path);
134 } // extern "C"
135 
137  BaseHadronizer(params),
138  needClear(false),
139  parameters(params.getParameter<edm::ParameterSet>("HerwigParameters")),
140  herwigVerbosity(params.getUntrackedParameter<int>("herwigVerbosity", 0)),
141  hepmcVerbosity(params.getUntrackedParameter<int>("hepmcVerbosity", 0)),
142  maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
143  printCards(params.getUntrackedParameter<bool>("printCards", false)),
144  emulatePythiaStatusCodes(params.getUntrackedParameter<bool>("emulatePythiaStatusCodes", false)),
145  comEnergy(params.getParameter<double>("comEnergy")),
146  useJimmy(params.getParameter<bool>("useJimmy")),
147  doMPInteraction(params.getParameter<bool>("doMPInteraction")),
148  numTrials(params.getUntrackedParameter<int>("numTrialsMPI", 100)),
149  readMCatNLOfile(false)
150 {
151 
152  fConvertToPDG = false;
153  if ( params.exists( "doPDGConvert" ) )
154  fConvertToPDG = params.getParameter<bool>("doPDGConvert");
155 }
156 
158 {
159  clear();
160 }
161 
163 {
164  if (!needClear)
165  return;
166 
167  // teminate elementary process
168  call(hwefin);
169  if (useJimmy)
170  call(jmefin);
171 
172  needClear = false;
173 }
174 
176  const std::vector<std::string> &lines)
177 {
178  std::set<std::string> blocks;
179  std::string block;
180  for(std::vector<std::string>::const_iterator iter = lines.begin();
181  iter != lines.end(); ++iter) {
182  std::string line = *iter;
183  std::transform(line.begin(), line.end(),
184  line.begin(), (int(*)(int))std::toupper);
185  std::string::size_type pos = line.find('#');
186  if (pos != std::string::npos)
187  line.resize(pos);
188 
189  if (line.empty())
190  continue;
191 
192  if (!boost::algorithm::is_space()(line[0])) {
193  std::vector<std::string> tokens;
194  boost::split(tokens, line,
195  boost::algorithm::is_space(),
196  boost::token_compress_on);
197  if (!tokens.size())
198  continue;
199  block.clear();
200  if (tokens.size() < 2)
201  continue;
202  if (tokens[0] == "BLOCK")
203  block = tokens[1];
204  else if (tokens[0] == "DECAY")
205  block = "DECAY";
206 
207  if (block.empty())
208  continue;
209 
210  if (!blocks.count(block)) {
211  blocks.insert(block);
212  edm::LogWarning("Generator|Herwig6Hadronzier")
213  << "Unsupported SLHA block \"" << block
214  << "\". It will be ignored."
215  << std::endl;
216  }
217  }
218  }
219 }
220 
222 {
223  clear();
224 
225  externalPartons = (heprup != 0);
226 
227  std::ostringstream info;
228  info << "---------------------------------------------------\n";
229  info << "Initializing Herwig6Hadronizer for "
230  << (externalPartons ? "external" : "internal") << " partons\n";
231  info << "---------------------------------------------------\n";
232 
233  info << " Herwig verbosity level = " << herwigVerbosity << "\n";
234  info << " HepMC verbosity = " << hepmcVerbosity << "\n";
235  info << " Number of events to be printed = " << maxEventsToPrint << "\n";
236 
237  // Call hwudat to set up HERWIG block data
238  hwudat();
239 
240  // Setting basic parameters
241  if (externalPartons) {
242  hwproc.PBEAM1 = heprup->EBMUP.first;
243  hwproc.PBEAM2 = heprup->EBMUP.second;
244  pdgToHerwig(heprup->IDBMUP.first, hwbmch.PART1);
245  pdgToHerwig(heprup->IDBMUP.second, hwbmch.PART2);
246  } else {
247  hwproc.PBEAM1 = 0.5 * comEnergy;
248  hwproc.PBEAM2 = 0.5 * comEnergy;
249  pdgToHerwig(2212, hwbmch.PART1);
250  pdgToHerwig(2212, hwbmch.PART2);
251  }
252 
253  if (useJimmy) {
254  info << " HERWIG will be using JIMMY for UE/MI.\n";
255  jmparm.MSFLAG = 1;
256  if (doMPInteraction)
257  info << " JIMMY trying to generate multiple interactions.\n";
258  }
259 
260 
261  // set the IPROC already here... needed for VB pairs
262 
263  bool iprocFound=false;
264 
266  line != parameters.end(); ++line) {
267  if(!strcmp((line->substr(0,5)).c_str(),"IPROC")) {
268  if (!give(*line))
270  << "Herwig 6 did not accept the following: \""
271  << *line << "\"." << std::endl;
272  else iprocFound=true;
273  }
274  }
275 
276  if (!iprocFound && !externalPartons)
278  << "You have to define the process with IPROC." << std::endl;
279 
280  // initialize other common blocks ...
281  call(hwigin);
282  hwevnt.MAXER = 100000000; // O(inf)
283  hwpram.LWSUD = 0; // don't write Sudakov form factors
284  hwdspn.LWDEC = 0; // don't write three/four body decays
285  // (no fort.77 and fort.88 ...)
286 
287  // init LHAPDF glue
288 
289  std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
290  for(unsigned int i = 0; i < 2; i++) {
291  hwpram.MODPDF[i] = -111;
292  std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
293  }
294 
295  if (useJimmy)
296  call(jimmin);
297 
298  hwevnt.MAXPR = maxEventsToPrint;
299  hwpram.IPRINT = herwigVerbosity;
300 // hwprop.RMASS[6] = 175.0; //FIXME
301 
302  if (printCards) {
303  info << "\n";
304  info << "------------------------------------\n";
305  info << "Reading HERWIG parameters\n";
306  info << "------------------------------------\n";
307 
308  }
310  line != parameters.end(); ++line) {
311  if (printCards)
312  info << " " << *line << "\n";
313  if (!give(*line))
315  << "Herwig 6 did not accept the following: \""
316  << *line << "\"." << std::endl;
317  }
318 
319  if (printCards)
320  info << "\n";
321 
322  if (externalPartons) {
323  std::vector<std::string> slha =
324  lheRunInfo()->findHeader("slha");
325  if (!slha.empty())
326  setSLHAFromHeader(slha);
327  }
328 
329  needClear = true;
330 
331  std::pair<int, int> pdfs(-1, -1);
332  if (externalPartons)
333  pdfs = lheRunInfo()->pdfSetTranslation();
334 
335  if (hwpram.MODPDF[0] != -111 || hwpram.MODPDF[1] != -111) {
336  for(unsigned int i = 0; i < 2; i++)
337  if (hwpram.MODPDF[i] == -111)
338  hwpram.MODPDF[i] = -1;
339 
340  if (pdfs.first != -1 || pdfs.second != -1)
341  edm::LogError("Generator|Herwig6Hadronzier")
342  << "Both external Les Houches event and "
343  "config file specify a PDF set. "
344  "User PDF will override external one."
345  << std::endl;
346 
347  pdfs.first = hwpram.MODPDF[0] != -111 ? hwpram.MODPDF[0] : -1;
348  pdfs.second = hwpram.MODPDF[1] != -111 ? hwpram.MODPDF[1] : -1;
349  }
350 
351  hwpram.MODPDF[0] = pdfs.first;
352  hwpram.MODPDF[1] = pdfs.second;
353 
354  if (externalPartons)
355  hwproc.IPROC = -1;
356 
357  // HERWIG preparations ...
358  call(hwuinc);
359  markStable(13); // MU+
360  markStable(-13); // MU-
361  markStable(3112); // SIGMA+
362  markStable(-3112); // SIGMABAR+
363  markStable(3222); // SIGMA-
364  markStable(-3222); // SIGMABAR-
365  markStable(3122); // LAMBDA0
366  markStable(-3122); // LAMBDABAR0
367  markStable(3312); // XI-
368  markStable(-3312); // XIBAR+
369  markStable(3322); // XI0
370  markStable(-3322); // XI0BAR
371  markStable(3334); // OMEGA-
372  markStable(-3334); // OMEGABAR+
373  markStable(211); // PI+
374  markStable(-211); // PI-
375  markStable(321); // K+
376  markStable(-321); // K-
377  markStable(310); // K_S0
378  markStable(130); // K_L0
379 
380  // better: merge with declareStableParticles
381  // and get the list from configuration / Geant4 / Core somewhere
382 
383  // initialize HERWIG event generation
384  call(hweini);
385 
386  if (useJimmy)
387  call(jminit);
388 
389  edm::LogInfo(info.str());
390 
391  return true;
392 }
393 
394 bool Herwig6Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
395 {
396  for(std::vector<int>::const_iterator iter = pdgIds.begin();
397  iter != pdgIds.end(); ++iter)
398  if (!markStable(*iter))
399  return false;
400  return true;
401 }
402 
403 bool Herwig6Hadronizer::declareSpecialSettings( const std::vector<std::string> )
404 {
405  return true;
406 }
407 
409 {
410  if (!runInfo().internalXSec()) {
411  // not set via LHE, so get it from HERWIG
412  // the reason is that HERWIG doesn't compute the xsec
413  // in all LHE modes
414 
415  double RNWGT = 1. / hwevnt.NWGTS;
416  double AVWGT = hwevnt.WGTSUM * RNWGT;
417 
418  double xsec = 1.0e3 * AVWGT;
419 
420  runInfo().setInternalXSec(xsec);
421  }
422 }
423 
425 {
426  // hard process generation, parton shower, hadron formation
427 
428  InstanceWrapper wrapper(this); // safe guard
429 
430  event().reset();
431 
432  // call herwig routines to create HEPEVT
433 
434  hwuine(); // initialize event
435 
436  if (callWithTimeout(10, hwepro)) { // process event and PS
437  // We hung for more than 10 seconds
438  int error = 199;
439  hwwarn_("HWHGUP", &error);
440  }
441 
442  hwbgen(); // parton cascades
443 
444  // call jimmy ... only if event is not killed yet by HERWIG
445  if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct()) {
447  return false;
448  }
449 
450  hwdhob(); // heavy quark decays
451  hwcfor(); // cluster formation
452  hwcdec(); // cluster decays
453 
454  // if event *not* killed by HERWIG, return true
455  if (hwevnt.IERROR) {
456  hwufne(); // finalize event, to keep system clean
458  return false;
459  }
460 
462  return true;
463 }
464 
466 {
468 
469  HepMC::PdfInfo pdfInfo;
470  if(externalPartons) {
471  lheEvent()->fillEventInfo( event().get() );
472  lheEvent()->fillPdfInfo( &pdfInfo );
473 
474  // for MC@NLO: IDWRUP is not filled...
475  if(event()->signal_process_id()==0)
476  event()->set_signal_process_id( abs(hwproc.IPROC) );
477 
478  }
479 
480  HepMC::GenParticle* incomingParton = NULL;
481  HepMC::GenParticle* targetParton = NULL;
482 
483  HepMC::GenParticle* incomingProton = NULL;
484  HepMC::GenParticle* targetProton = NULL;
485 
486  // find incoming parton (first entry with IST=121)
487  for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin();
488  (it != event()->particles_end() && incomingParton==NULL); it++)
489  if((*it)->status()==121) incomingParton = (*it);
490 
491  // find target parton (first entry with IST=122)
492  for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin();
493  (it != event()->particles_end() && targetParton==NULL); it++)
494  if((*it)->status()==122) targetParton = (*it);
495 
496  // find incoming Proton (first entry ID=2212, IST=101)
497  for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin();
498  (it != event()->particles_end() && incomingProton==NULL); it++)
499  if((*it)->status()==101 && (*it)->pdg_id()==2212) incomingProton = (*it);
500 
501  // find target Proton (first entry ID=2212, IST=102)
502  for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin();
503  (it != event()->particles_end() && targetProton==NULL); it++)
504  if((*it)->status()==102 && (*it)->pdg_id()==2212) targetProton = (*it);
505 
506  // find hard scale Q (computed from colliding partons)
507  if( incomingParton && targetParton ) {
508  math::XYZTLorentzVector totMomentum(0,0,0,0);
509  totMomentum+=incomingParton->momentum();
510  totMomentum+=targetParton->momentum();
511  double evScale = totMomentum.mass();
512  double evScale2 = evScale*evScale;
513 
514  // find alpha_QED & alpha_QCD
515  int one=1;
516  double alphaQCD=hwualf_(&one,&evScale);
517  double alphaQED=hwuaem_(&evScale2);
518 
519  if(!externalPartons || event()->event_scale() < 0) event()->set_event_scale(evScale);
520  if(!externalPartons || event()->alphaQCD() < 0) event()->set_alphaQCD(alphaQCD);
521  if(!externalPartons || event()->alphaQED() < 0) event()->set_alphaQED(alphaQED);
522 
523  if(!externalPartons || pdfInfo.x1() < 0) {
524  // get the PDF information
525  pdfInfo.set_id1( incomingParton->pdg_id()==21 ? 0 : incomingParton->pdg_id());
526  pdfInfo.set_id2( targetParton->pdg_id()==21 ? 0 : targetParton->pdg_id());
527  if( incomingProton && targetProton ) {
528  double x1 = incomingParton->momentum().pz()/incomingProton->momentum().pz();
529  double x2 = targetParton->momentum().pz()/targetProton->momentum().pz();
530  pdfInfo.set_x1(x1);
531  pdfInfo.set_x2(x2);
532  }
533  // we do not fill pdf1 & pdf2, since they are not easily available (what are they needed for anyways???)
534  pdfInfo.set_scalePDF(evScale); // the same as Q above... does this make sense?
535  }
536 
537  if(!externalPartons || event()->signal_process_id() < 0) event()->set_signal_process_id( abs(hwproc.IPROC) );
538  event()->set_pdf_info( pdfInfo );
539  }
540 
541  // add event weight & PDF information
542  if (lheRunInfo() != 0 && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
543  // in LHE weighting mode 4 the weight is an xsec, so convert form HERWIG
544  // to standard CMS unit "picobarn"
545  event()->weights().push_back( 1.0e3 * hwevnt.EVWGT );
546  else
547  event()->weights().push_back( hwevnt.EVWGT );
548 
549 
550  // find final parton (first entry with IST=123)
551  HepMC::GenParticle* finalParton = NULL;
552  for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin();
553  (it != event()->particles_end() && finalParton==NULL); it++)
554  if((*it)->status()==123) finalParton = (*it);
555 
556 
557  // add GenEventInfo & binning Values
558  eventInfo().reset(new GenEventInfoProduct(event().get()));
559  if(finalParton) {
560  double thisPt=finalParton->momentum().perp();
561  eventInfo()->setBinningValues(std::vector<double>(1, thisPt));
562  }
563 
564  // emulate PY6 status codes, if switched on...
567 
568 }
569 
570 
572 {
573  // hadron decays
574 
575  InstanceWrapper wrapper(this); // safe guard
576 
577  hwdhad(); // unstable particle decays
578  hwdhvy(); // heavy flavour decays
579  hwmevt(); // soft underlying event
580 
581  hwufne(); // finalize event
582 
583  if (hwevnt.IERROR)
584  return false;
585 
586  event().reset(new HepMC::GenEvent);
587  if (!conv.fill_next_event(event().get()))
588  throw cms::Exception("Herwig6Error")
589  << "HepMC Conversion problems in event." << std::endl;
590 
591  // do particle ID conversion Herwig->PDG, if requested
592  if ( fConvertToPDG ) {
593  for ( HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end(); ++part) {
594  if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
595  (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
596  }
597  }
598 
599  return true;
600 }
601 
603 {
604  return true;
605 }
606 
608 {
610  heprup_.pdfgup[0] = heprup_.pdfgup[1] = -1;
611  heprup_.pdfsup[0] = heprup_.pdfsup[1] = -1;
612  // we set up the PDFs ourselves
613 
614  // pass HERWIG paramaters fomr header (if present)
615  std::string mcnloHeader="herwig6header";
616  std::vector<lhef::LHERunInfo::Header> headers=lheRunInfo()->getHeaders();
617  for(std::vector<lhef::LHERunInfo::Header>::const_iterator hIter=headers.begin();hIter!=headers.end(); ++hIter) {
618  if(hIter->tag()==mcnloHeader){
619  readMCatNLOfile=true;
620  for(lhef::LHERunInfo::Header::const_iterator lIter=hIter->begin(); lIter != hIter->end(); ++lIter) {
621  if((lIter->c_str())[0]!='#' && (lIter->c_str())[0]!='\n') { // it's not a comment)
622  if (!give(*lIter))
624  << "Herwig 6 did not accept the following: \""
625  << *lIter << "\"." << std::endl;
626  }
627  }
628  }
629  }
630 }
631 
633 {
635 
636  // if MCatNLO external file is read, read comment & pass IHPRO to HERWIG
637  if(readMCatNLOfile) {
638  for(std::vector<std::string>::const_iterator iter=lheEvent()->getComments().begin();
639  iter!=lheEvent()->getComments().end(); ++iter) {
640  std::string toParse(iter->substr(1));
641  if (!give(toParse))
643  << "Herwig 6 did not accept the following: \""
644  << toParse << "\"." << std::endl;
645  }
646  }
647 
648 }
649 
651 {
652  int status = p->status();
653 
654  // weird 9922212 particles...
655  if (status == 3 && !p->end_vertex())
656  return 2;
657 
658  if (status >= 1 && status <= 3)
659  return status;
660 
661  if (!p->end_vertex())
662  return 1;
663 
664  // let's prevent particles having status 3, if the identical
665  // particle downstream is a better status 3 candidate
666  int currentId = p->pdg_id();
667  int orig = status;
668  if (status == 123 || status == 124 ||
669  status == 155 || status == 156 || status == 160 ||
670  (status >= 195 && status <= 197)) {
671  for(const HepMC::GenParticle *q = p;;) {
672  const HepMC::GenVertex *vtx = q->end_vertex();
673  if (!vtx)
674  break;
675 
676  HepMC::GenVertex::particles_out_const_iterator iter;
677  for(iter = vtx->particles_out_const_begin();
678  iter != vtx->particles_out_const_end(); ++iter)
679  if ((*iter)->pdg_id() == currentId)
680  break;
681 
682  if (iter == vtx->particles_out_const_end())
683  break;
684 
685  q = *iter;
686  if (q->status() == 3 ||
687  ((status == 120 || status == 123 ||
688  status == 124) && orig > 124))
689  return 4;
690  }
691  }
692 
693  int nesting = 0;
694  for(;;) {
695  if ((status >= 120 && status <= 122) || status == 3) {
696  // avoid flagging status 3 if there is a
697  // better status 3 candidate upstream
698  if (externalPartons)
699  return ((orig >= 121 && orig <= 124) ||
700  orig == 3) ? 3 : 4;
701  else
702  return (nesting ||
703  (status != 3 && orig <= 124)) ? 3 : 4;
704  }
705 
706  // check whether we are leaving the hard process
707  // including heavy resonance decays
708  if (!(status == 4 || status == 123 || status == 124 ||
709  status == 155 || status == 156 || status == 160 ||
710  (status >= 195 && status <= 197)))
711  break;
712 
713  const HepMC::GenVertex *vtx = p->production_vertex();
714  if (!vtx || !vtx->particles_in_size())
715  break;
716 
717  p = *vtx->particles_in_const_begin();
718  status = p->status();
719 
720  int newId = p->pdg_id();
721 
722  if (!newId)
723  break;
724 
725  // nesting increases if we move to the next-best mother
726  if (newId != currentId) {
727  if (++nesting > 1 && externalPartons)
728  break;
729  currentId = newId;
730  }
731  }
732 
733  return 2;
734 }
735 
737 {
738  for(HepMC::GenEvent::particle_iterator iter =
739  event()->particles_begin();
740  iter != event()->particles_end(); iter++)
741  (*iter)->set_status(pythiaStatusCode(*iter));
742 
743  for(HepMC::GenEvent::particle_iterator iter =
744  event()->particles_begin();
745  iter != event()->particles_end(); iter++)
746  if ((*iter)->status() == 4)
747  (*iter)->set_status(2);
748 }
749 
751 
753 DEFINE_FWK_MODULE(Herwig6GeneratorFilter);
754 
756 DEFINE_FWK_MODULE(Herwig6HadronizerFilter);
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void call(void(&fn)())
void setSLHAFromHeader(const std::vector< std::string > &lines)
bool declareStableParticles(const std::vector< int > &pdgIds)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void fillEventInfo(HepMC::GenEvent *hepmc) const
Definition: LHEEvent.cc:240
std::pair< double, double > EBMUP
Definition: LesHouches.h:78
const char * classname() const
#define hwdspn
Definition: herwig.h:223
#define jminit
Definition: herwig.h:320
int pdfgup[2]
edm::HadronizerFilter< Herwig6Hadronizer, gen::ExternalDecayDriver > Herwig6HadronizerFilter
bool exists(std::string const &parameterName) const
checks if a parameter exists
HepMC::IO_HERWIG conv
#define abs(x)
Definition: mlp_lapack.h:159
static void fixHepMCEventTimeOrdering(HepMC::GenEvent *event)
Definition: LHEEvent.cc:503
#define NULL
Definition: scimark2.h:8
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:198
std::pair< int, int > IDBMUP
Definition: LesHouches.h:73
struct HEPRUP_ heprup_
bool generatePartonsAndHadronize()
std::vector< std::string >::const_iterator const_iterator
edm::GeneratorFilter< Herwig6Hadronizer, gen::ExternalDecayDriver > Herwig6GeneratorFilter
virtual void upInit()
bool initialize(const lhef::HEPRUP *heprup)
void setInternalXSec(const XSec &xsec)
#define hwprch
Definition: herwig.h:63
uint16_t size_type
std::auto_ptr< HepMC::GenEvent > & event()
int pythiaStatusCode(const HepMC::GenParticle *p) const
#define hwmsct
Definition: herwig.h:321
double hwualf_(int *mode, double *scale)
int path() const
Definition: HLTadd.h:3
double hwuaem_(double *scale)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
GenRunInfoProduct & runInfo()
bool initializeForInternalPartons()
#define jmparm
Definition: herwig.h:301
void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8])
const std::vector< std::string > & getComments() const
Definition: LHEEvent.h:45
void mysetpdfpath_(const char *path)
lhef::LHEEvent * lheEvent()
bool callWithTimeout(unsigned int secs, void(*fn)())
void hwwarn_(const char *method, int *id)
tuple result
Definition: query.py:137
static void fillHEPEUP(const HEPEUP *hepeup)
void fillPdfInfo(HepMC::PdfInfo *info) const
Definition: LHEEvent.cc:212
const std::vector< Header > & getHeaders() const
Definition: LHERunInfo.h:58
block
Formating index page&#39;s pieces.
Definition: Association.py:187
std::auto_ptr< GenEventInfoProduct > & eventInfo()
bool give(const std::string &line)
lhef::LHERunInfo * lheRunInfo()
#define jmefin
Definition: herwig.h:322
virtual void upEvnt()
part
Definition: HCALResponse.h:21
std::pair< int, int > pdfSetTranslation() const
Definition: LHERunInfo.cc:502
int mode
Definition: AMPTWrapper.h:139
bool declareSpecialSettings(const std::vector< std::string >)
#define begin
Definition: vmac.h:31
const_iterator end() const
Herwig6Hadronizer(const edm::ParameterSet &params)
const_iterator begin() const
void setherwpdf_(void)
std::vector< std::string > findHeader(const std::string &tag) const
Definition: LHERunInfo.cc:390
#define jimmin
Definition: herwig.h:319
int pdfsup[2]
tuple status
Definition: ntuplemaker.py:245
gen::ParameterCollector parameters
double split
Definition: MVATrainer.cc:139
bool initializeForExternalPartons()
static void fillHEPRUP(const HEPRUP *heprup)
static HepMC::HEPEVT_Wrapper wrapper