CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PomwigHadronizer.cc
Go to the documentation of this file.
2 
3 #include <cstring>
4 #include <memory>
5 
6 #include <map>
7 #include <memory>
8 #include <set>
9 #include <sstream>
10 #include <string>
11 #include <vector>
12 
13 #include <boost/algorithm/string/classification.hpp>
14 #include <boost/algorithm/string/split.hpp>
15 
16 #include <HepMC/GenEvent.h>
17 #include <HepMC/GenParticle.h>
18 #include <HepMC/GenVertex.h>
19 #include <HepMC/PdfInfo.h>
20 #include <HepMC/HerwigWrapper.h>
21 #include <HepMC/HEPEVT_Wrapper.h>
22 #include <HepMC/IO_HERWIG.h>
23 #include "HepPID/ParticleIDTranslations.hh"
24 
27 
30 
33 
37 
39 
42 
43 namespace gen {
44  extern "C" {
45  void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8]);
46  }
47 
48  // helpers
49  namespace {
50  inline bool call_hwmsct() {
51  int result;
52  hwmsct(&result);
53  return result;
54  }
55 
56  int pdgToHerwig(int ipdg, char *nwig) {
57  int iopt = 1;
58  int iwig = 0;
59  hwuidt_(&iopt, &ipdg, &iwig, nwig);
60  return ipdg ? iwig : 0;
61  }
62 
63  bool markStable(int pdgId) {
64  char nwig[9] = " ";
65  if (!pdgToHerwig(pdgId, nwig))
66  return false;
67  hwusta(nwig, 1);
68  return true;
69  }
70  } // namespace
71 
72 #define qcd_1994 qcd_1994_
73  extern "C" {
74  void qcd_1994(double &, double &, double *, int &);
75  }
76 // For H1 2006 fits
77 #define qcd_2006 qcd_2006_
78  extern "C" {
79  void qcd_2006(double &, double &, int &, double *, double *, double *, double *, double *);
80  }
81 
82  extern "C" {
83  void hwwarn_(const char *method, int *id);
84  void setherwpdf_(void);
85  void mysetpdfpath_(const char *path);
86  }
87 
90 
92  : BaseHadronizer(params),
93  needClear(false),
94  parameters(params.getParameter<edm::ParameterSet>("HerwigParameters")),
95  herwigVerbosity(params.getUntrackedParameter<int>("herwigVerbosity", 0)),
96  hepmcVerbosity(params.getUntrackedParameter<int>("hepmcVerbosity", 0)),
97  maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
98  printCards(params.getUntrackedParameter<bool>("printCards", false)),
99  comEnergy(params.getParameter<double>("comEnergy")),
100  survivalProbability(params.getParameter<double>("survivalProbability")),
101  diffTopology(params.getParameter<int>("diffTopology")),
102  h1fit(params.getParameter<int>("h1fit")),
103  useJimmy(params.getParameter<bool>("useJimmy")),
104  doMPInteraction(params.getParameter<bool>("doMPInteraction")),
105  numTrials(params.getUntrackedParameter<int>("numTrialsMPI", 100)),
106  doPDGConvert(false) {
107  if (params.exists("doPDGConvert"))
108  doPDGConvert = params.getParameter<bool>("doPDGConvert");
109  }
110 
112 
113  void PomwigHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine *v) { setHerwigRandomEngine(v); }
114 
116  if (!needClear)
117  return;
118 
119  // teminate elementary process
120  call(hwefin);
121  if (useJimmy)
122  call(jmefin);
123 
124  needClear = false;
125  }
126 
128 
130  clear();
131 
132  edm::LogVerbatim("") << "----------------------------------------------\n"
133  << "Initializing PomwigHadronizer\n"
134  << "----------------------------------------------\n";
135 
136  // Call hwudat to set up HERWIG block data
137  hwudat();
138 
139  // Setting basic parameters ...
140  hwproc.PBEAM1 = comEnergy / 2.;
141  hwproc.PBEAM2 = comEnergy / 2.;
142  // Choose beam particles for POMWIG depending on topology
143  switch (diffTopology) {
144  case 0: //DPE
145  hwbmch.PART1[0] = 'E';
146  hwbmch.PART1[1] = '-';
147  hwbmch.PART2[0] = 'E';
148  hwbmch.PART2[1] = '-';
149  break;
150  case 1: //SD survive PART1
151  hwbmch.PART1[0] = 'E';
152  hwbmch.PART1[1] = '-';
153  hwbmch.PART2[0] = 'P';
154  hwbmch.PART2[1] = ' ';
155  break;
156  case 2: //SD survive PART2
157  hwbmch.PART1[0] = 'P';
158  hwbmch.PART1[1] = ' ';
159  hwbmch.PART2[0] = 'E';
160  hwbmch.PART2[1] = '-';
161  break;
162  case 3: //Non diffractive
163  hwbmch.PART1[0] = 'P';
164  hwbmch.PART1[1] = ' ';
165  hwbmch.PART2[0] = 'P';
166  hwbmch.PART2[1] = ' ';
167  break;
168  default:
169  throw edm::Exception(edm::errors::Configuration, "PomwigError")
170  << " Invalid Diff. Topology. Must be DPE(diffTopology = 0), SD particle 1 (diffTopology = 1), SD particle "
171  "2 (diffTopology = 2) and Non diffractive (diffTopology = 3)";
172  break;
173  }
174  for (int i = 2; i < 8; ++i) {
175  hwbmch.PART1[i] = ' ';
176  hwbmch.PART2[i] = ' ';
177  }
178 
179  // initialize other common blocks ...
180  call(hwigin);
181 
182  hwevnt.MAXER = 100000000; // O(inf)
183  hwpram.LWSUD = 0; // don't write Sudakov form factors
184  hwdspn.LWDEC = 0; // don't write three/four body decays
185  // (no fort.77 and fort.88 ...)a
186 
187  std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
188  for (unsigned int i = 0; i < 2; i++) {
189  hwpram.MODPDF[i] = -111;
190  std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
191  }
192 
193  hwevnt.MAXPR = maxEventsToPrint;
194  hwpram.IPRINT = herwigVerbosity;
195 
196  edm::LogVerbatim("") << "------------------------------------\n"
197  << "Reading HERWIG parameters\n"
198  << "------------------------------------\n";
199 
201  edm::LogVerbatim("") << " " << *line;
202  if (!give(*line))
204  << "Herwig 6 did not accept the following: \"" << *line << "\"." << std::endl;
205  }
206 
207  needClear = true;
208 
209  return true;
210  }
211 
213  call(hwuinc);
214 
215  hwusta("PI0 ", 1);
216 
217  if (!initializeDPDF())
218  return false;
219 
220  call(hweini);
221 
222  return true;
223  }
224 
226  // Initialize H1 pomeron/reggeon
227 
228  if (diffTopology == 3)
229  return true;
230 
231  if ((diffTopology != 0) && (diffTopology != 1) && (diffTopology != 2))
232  return false;
233 
234  int nstru = hwpram.NSTRU;
235  int ifit = h1fit;
236  if ((nstru == 9) || (nstru == 10)) {
237  if ((ifit <= 0) || (ifit >= 7)) {
238  throw edm::Exception(edm::errors::Configuration, "PomwigError")
239  << " Attempted to set non existant H1 1997 fit index. Has to be 1...6";
240  }
241  std::string aux((nstru == 9) ? "Pomeron" : "Reggeon");
242  edm::LogVerbatim("") << " H1 1997 pdf's: " << aux << "\n"
243  << " IFIT = " << ifit;
244  double xp = 0.1;
245  double Q2 = 75.0;
246  double xpq[13];
247  qcd_1994(xp, Q2, xpq, ifit);
248  } else if ((nstru >= 12) && (nstru <= 15)) {
249  bool isPom = (nstru == 12) || (nstru == 14);
250  bool isFitA = (nstru == 12) || (nstru == 13);
251  ifit = (isFitA) ? 1 : 2;
252  std::string aux_0((isFitA) ? "A" : "B");
253  std::string aux_1((isPom) ? "Pomeron" : "Reggeon");
254  edm::LogVerbatim("") << " H1 2006 Fit " << aux_0 << " " << aux_1 << "\n"
255  << " IFIT = " << ifit;
256  double xp = 0.1;
257  double Q2 = 75.0;
258  double xpq[13];
259  double f2[2];
260  double fl[2];
261  double c2[2];
262  double cl[2];
263  qcd_2006(xp, Q2, ifit, xpq, f2, fl, c2, cl);
264  } else {
265  throw edm::Exception(edm::errors::Configuration, "PomwigError")
266  << " Only running Pomeron H1 1997 (NSTRU=9), H1 2006 fit A (NSTRU=12) and H1 2006 fit B (NSTRU=14) or "
267  "Reggeon H1 1997 (NSTRU=10), H1 2006 fit A (NSTRU=13) and H1 2006 fit B (NSTRU=15)";
268  }
269 
270  return true;
271  }
272 
273  bool PomwigHadronizer::declareStableParticles(const std::vector<int> &pdgIds) {
274  for (std::vector<int>::const_iterator iter = pdgIds.begin(); iter != pdgIds.end(); ++iter)
275  if (!markStable(*iter))
276  return false;
277  return true;
278  }
279 
281  double RNWGT = 1. / hwevnt.NWGTS;
282  double AVWGT = hwevnt.WGTSUM * RNWGT;
283 
284  double xsec = 1.0e3 * AVWGT;
285  xsec = survivalProbability * xsec;
286 
287  runInfo().setInternalXSec(xsec);
288  }
289 
290  bool PomwigHadronizer::hadronize() { return false; }
291 
293  // hard process generation, parton shower, hadron formation
294 
295  InstanceWrapper wrapper(this); // safe guard
296 
297  event().reset();
298 
299  // call herwig routines to create HEPEVT
300 
301  hwuine(); // initialize event
302 
303  if (callWithTimeout(10, hwepro)) { // process event and PS
304  // We hung for more than 10 seconds
305  int error = 199;
306  hwwarn_("HWHGUP", &error);
307  }
308 
309  hwbgen(); // parton cascades
310 
311  // call jimmy ... only if event is not killed yet by HERWIG
312  if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct())
313  return false;
314 
315  hwdhob(); // heavy quark decays
316  hwcfor(); // cluster formation
317  hwcdec(); // cluster decays
318 
319  // if event *not* killed by HERWIG, return true
320  if (!hwevnt.IERROR)
321  return true;
322 
323  hwufne(); // finalize event
324  return false;
325  }
326 
329 
330  event()->set_signal_process_id(hwproc.IPROC);
331 
332  event()->weights().push_back(hwevnt.EVWGT);
333  }
334 
336  // hadron decays
337 
338  InstanceWrapper wrapper(this); // safe guard
339 
340  hwdhad(); // unstable particle decays
341  hwdhvy(); // heavy flavour decays
342  hwmevt(); // soft underlying event
343 
344  hwufne(); // finalize event
345 
346  if (hwevnt.IERROR)
347  return false;
348 
349  event() = std::make_unique<HepMC::GenEvent>();
350  if (!conv.fill_next_event(event().get()))
351  throw cms::Exception("PomwigError") << "HepMC Conversion problems in event." << std::endl;
352 
353  // do particle ID conversion Herwig->PDG, if requested
354  if (doPDGConvert) {
355  for (HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end();
356  ++part) {
357  if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
358  (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
359  }
360  }
361 
362  return true;
363  }
364 
365  bool PomwigHadronizer::residualDecay() { return true; }
366 
367 } //namespace gen
bool declareStableParticles(const std::vector< int > &pdgIds)
Log< level::Info, true > LogVerbatim
void call(void(&fn)())
void mysetpdfpath_(const char *path)
bool exists(std::string const &parameterName) const
checks if a parameter exists
static void fixHepMCEventTimeOrdering(HepMC::GenEvent *event)
Definition: LHEEvent.cc:453
gen::ParameterCollector parameters
void setInternalXSec(const XSec &xsec)
double v[5][pyjets_maxn]
tuple cl
Definition: haddnano.py:49
tuple result
Definition: mps_fire.py:311
GenRunInfoProduct & runInfo()
static const std::vector< std::string > theSharedResources
#define hwprch
Definition: herwig.h:60
bool callWithTimeout(unsigned int secs, void(*fn)())
void hwwarn_(const char *method, int *id)
void setherwpdf_(void)
#define qcd_2006
#define jmefin
Definition: herwig.h:318
std::unique_ptr< HepMC::GenEvent > & event()
HepMC::IO_HERWIG conv
bool give(const std::string &line)
static const std::string kFortranInstance
#define qcd_1994
void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8])
part
Definition: HCALResponse.h:20
#define hwdspn
Definition: herwig.h:220
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
#define hwmsct
Definition: herwig.h:317
const_iterator end() const
void setHerwigRandomEngine(CLHEP::HepRandomEngine *v)
const_iterator begin() const
PomwigHadronizer(const edm::ParameterSet &params)
static const std::string kHerwig6
static HepMC::HEPEVT_Wrapper wrapper