CMS 3D CMS Logo

HydjetHadronizer.cc
Go to the documentation of this file.
1 
7 #include <iostream>
8 #include <cmath>
9 
16 
22 
23 #include "HepMC/IO_HEPEVT.h"
24 #include "HepMC/PythiaWrapper6_4.h"
25 #include "HepMC/GenEvent.h"
26 #include "HepMC/HeavyIon.h"
27 #include "HepMC/SimpleVector.h"
28 
33 
34 using namespace edm;
35 using namespace std;
36 using namespace gen;
37 
38 int HydjetHadronizer::convertStatus(int st) {
39  if (st <= 0)
40  return 0;
41  if (st <= 10)
42  return 1;
43  if (st <= 20)
44  return 2;
45  if (st <= 30)
46  return 3;
47  else
48  return st;
49 }
50 
51 const std::vector<std::string> HydjetHadronizer::theSharedResources = {edm::SharedResourceNames::kPythia6,
53 
54 //_____________________________________________________________________
55 HydjetHadronizer::HydjetHadronizer(const ParameterSet& pset, edm::ConsumesCollector&& iC)
57  evt(nullptr),
58  pset_(pset),
59  abeamtarget_(pset.getParameter<double>("aBeamTarget")),
60  angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
61  bfixed_(pset.getParameter<double>("bFixed")),
62  bmax_(pset.getParameter<double>("bMax")),
63  bmin_(pset.getParameter<double>("bMin")),
64  cflag_(pset.getParameter<int>("cFlag")),
65  embedding_(pset.getParameter<int>("embeddingMode")),
66  comenergy(pset.getParameter<double>("comEnergy")),
67  doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
68  docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
69  fracsoftmult_(pset.getParameter<double>("fracSoftMultiplicity")),
70  hadfreeztemp_(pset.getParameter<double>("hadronFreezoutTemperature")),
71  hymode_(pset.getParameter<string>("hydjetMode")),
72  maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)),
73  maxlongy_(pset.getParameter<double>("maxLongitudinalRapidity")),
74  maxtrany_(pset.getParameter<double>("maxTransverseRapidity")),
75  nsub_(0),
76  nhard_(0),
77  nmultiplicity_(pset.getParameter<int>("nMultiplicity")),
78  nsoft_(0),
79  nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
80  pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
81  qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
82  qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
83  phi0_(0.),
84  sinphi0_(0.),
85  cosphi0_(1.),
86  rotate_(pset.getParameter<bool>("rotateEventPlane")),
87  shadowingswitch_(pset.getParameter<int>("shadowingSwitch")),
88  signn_(pset.getParameter<double>("sigmaInelNN")),
89  fVertex_(nullptr),
90  pythia6Service_(new Pythia6Service(pset)) {
91  // Default constructor
92 
93  if (pset.exists("signalVtx"))
94  signalVtx_ = pset.getUntrackedParameter<std::vector<double> >("signalVtx");
95 
96  if (signalVtx_.size() == 4) {
97  if (!fVertex_)
98  fVertex_ = new HepMC::FourVector();
99  LogDebug("EventSignalVertex") << "Setting event signal vertex "
100  << " x = " << signalVtx_.at(0) << " y = " << signalVtx_.at(1)
101  << " z= " << signalVtx_.at(2) << " t = " << signalVtx_.at(3) << endl;
102  fVertex_->set(signalVtx_.at(0), signalVtx_.at(1), signalVtx_.at(2), signalVtx_.at(3));
103  }
104 
105  // PYLIST Verbosity Level
106  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
107  pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
108  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_;
109 
110  //Max number of events printed on verbosity level
111  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint", 0);
112  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_;
113 
114  if (embedding_) {
115  cflag_ = 0;
116  src_ = iC.consumes<CrossingFrame<edm::HepMCProduct> >(
117  pset.getUntrackedParameter<edm::InputTag>("backgroundLabel", edm::InputTag("mix", "generatorSmeared")));
118  }
119 
120  int cm = 1, va, vb, vc;
121  HYJVER(cm, va, vb, vc);
122  HepMC::HEPEVT_Wrapper::set_max_number_entries(4000);
123 }
124 
125 //_____________________________________________________________________
127  // destructor
128  call_pystat(1);
129  delete pythia6Service_;
130 }
131 
132 //_____________________________________________________________________
134 
135 //_____________________________________________________________________
137  // heavy ion record in the final CMSSW Event
138  double npart = hyfpar.npart;
139  int nproj = static_cast<int>(npart / 2);
140  int ntarg = static_cast<int>(npart - nproj);
141 
142  HepMC::HeavyIon* hi = new HepMC::HeavyIon(nsub_, // Ncoll_hard/N of SubEvents
143  nproj, // Npart_proj
144  ntarg, // Npart_targ
145  static_cast<int>(hyfpar.nbcol), // Ncoll
146  0, // spectator_neutrons
147  0, // spectator_protons
148  0, // N_Nwounded_collisions
149  0, // Nwounded_N_collisions
150  0, // Nwounded_Nwounded_collisions
151  hyfpar.bgen * nuclear_radius(), // impact_parameter in [fm]
152  phi0_, // event_plane_angle
153  0, //hypsi3.psi3, // eccentricity
154  hyjpar.sigin // sigma_inel_NN
155  );
156 
157  evt->set_heavy_ion(*hi);
158  delete hi;
159 }
160 
161 //___________________________________________________________________
163  // Build particle object corresponding to index in hyjets (soft+hard)
164  double x0 = hyjets.phj[0][index];
165  double y0 = hyjets.phj[1][index];
166 
167  double x = x0 * cosphi0_ - y0 * sinphi0_;
168  double y = y0 * cosphi0_ + x0 * sinphi0_;
169 
170  HepMC::GenParticle* p = new HepMC::GenParticle(HepMC::FourVector(x, // px
171  y, // py
172  hyjets.phj[2][index], // pz
173  hyjets.phj[3][index]), // E
174  hyjets.khj[1][index], // id
175  convertStatus(hyjets.khj[0][index] // status
176  ));
177 
178  p->suggest_barcode(barcode);
179  return p;
180 }
181 
182 //___________________________________________________________________
183 HepMC::GenVertex* HydjetHadronizer::build_hyjet_vertex(int i, int id) {
184  // build verteces for the hyjets stored events
185  double x0 = hyjets.vhj[0][i];
186  double y0 = hyjets.vhj[1][i];
187  double x = x0 * cosphi0_ - y0 * sinphi0_;
188  double y = y0 * cosphi0_ + x0 * sinphi0_;
189  double z = hyjets.vhj[2][i];
190  double t = hyjets.vhj[4][i];
191 
192  HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x, y, z, t), id);
193  return vertex;
194 }
195 
196 //___________________________________________________________________
197 
200 
201  // generate single event
202  if (embedding_) {
203  const edm::Event& e = getEDMEvent();
204  HepMC::GenVertex* genvtx = nullptr;
205  const HepMC::GenEvent* inev = nullptr;
207  e.getByToken(src_, cf);
209  if (mix.size() < 1) {
210  throw cms::Exception("MatchVtx") << "Mixing has " << mix.size() << " sub-events, should have been at least 1"
211  << endl;
212  }
213  const HepMCProduct& bkg = mix.getObject(0);
214  if (!(bkg.isVtxGenApplied())) {
215  throw cms::Exception("MatchVtx") << "Input background does not have smeared vertex!" << endl;
216  } else {
217  inev = bkg.GetEvent();
218  }
219 
220  genvtx = inev->signal_process_vertex();
221 
222  if (!genvtx)
223  throw cms::Exception("MatchVtx") << "Input background does not have signal process vertex!" << endl;
224 
225  double aX, aY, aZ, aT;
226 
227  aX = genvtx->position().x();
228  aY = genvtx->position().y();
229  aZ = genvtx->position().z();
230  aT = genvtx->position().t();
231 
232  if (!fVertex_) {
233  fVertex_ = new HepMC::FourVector();
234  }
235  LogInfo("MatchVtx") << " setting vertex "
236  << " aX " << aX << " aY " << aY << " aZ " << aZ << " aT " << aT << endl;
237  fVertex_->set(aX, aY, aZ, aT);
238 
239  const HepMC::HeavyIon* hi = inev->heavy_ion();
240 
241  if (hi) {
242  bfixed_ = (hi->impact_parameter()) / nuclear_radius();
243  phi0_ = hi->event_plane_angle();
244  sinphi0_ = sin(phi0_);
245  cosphi0_ = cos(phi0_);
246  } else {
247  LogWarning("EventEmbedding") << "Background event does not have heavy ion record!";
248  }
249 
250  } else if (rotate_)
251  rotateEvtPlane();
252 
253  nsoft_ = 0;
254  nhard_ = 0;
255 
256  edm::LogInfo("HYDJETmode") << "##### HYDJET nhsel = " << hyjpar.nhsel;
257  edm::LogInfo("HYDJETfpart") << "##### HYDJET fpart = " << hyflow.fpart;
258  edm::LogInfo("HYDJETtf") << "##### HYDJET hadron freez-out temp, Tf = " << hyflow.Tf;
259  edm::LogInfo("HYDJETinTemp") << "##### HYDJET: QGP init temperature, T0 =" << pyqpar.T0u;
260  edm::LogInfo("HYDJETinTau") << "##### HYDJET: QGP formation time,tau0 =" << pyqpar.tau0u;
261 
262  int ntry = 0;
263  while (nsoft_ == 0 && nhard_ == 0) {
264  if (ntry > 100) {
265  edm::LogError("HydjetEmptyEvent") << "##### HYDJET: No Particles generated, Number of tries =" << ntry;
266 
267  // Throw an exception. Use the EventCorruption exception since it maps onto SkipEvent
268  // which is what we want to do here.
269 
270  std::ostringstream sstr;
271  sstr << "HydjetHadronizerProducer: No particles generated after " << ntry << " tries.\n";
272  edm::Exception except(edm::errors::EventCorruption, sstr.str());
273  throw except;
274  } else {
275  HYEVNT(bfixed_);
276  nsoft_ = hyfpar.nhyd;
277  nsub_ = hyjpar.njet;
278  nhard_ = hyfpar.npyt;
279  ++ntry;
280  }
281  }
282 
283  if (hyjpar.nhsel < 3)
284  nsub_++;
285 
286  // event information
287  std::unique_ptr<HepMC::GenEvent> evt = std::make_unique<HepMC::GenEvent>();
288  std::unique_ptr<edm::HepMCProduct> HepMCEvt = std::make_unique<edm::HepMCProduct>();
289 
290  if (nhard_ > 0 || nsoft_ > 0)
291  get_particles(evt.get());
292 
293  evt->set_signal_process_id(pypars.msti[0]); // type of the process
294  evt->set_event_scale(pypars.pari[16]); // Q^2
295  add_heavy_ion_rec(evt.get());
296 
297  if (fVertex_) {
298  // generate new vertex & apply the shift
299  // Copy the HepMC::GenEvent
300  HepMCEvt = std::make_unique<edm::HepMCProduct>(evt.get());
301  HepMCEvt->applyVtxGen(fVertex_);
302  evt = std::make_unique<HepMC::GenEvent>(*HepMCEvt->GetEvent());
303  }
304 
305  HepMC::HEPEVT_Wrapper::check_hepevt_consistency();
306  LogDebug("HEPEVT_info") << "Ev numb: " << HepMC::HEPEVT_Wrapper::event_number()
307  << " Entries number: " << HepMC::HEPEVT_Wrapper::number_entries() << " Max. entries "
308  << HepMC::HEPEVT_Wrapper::max_number_entries() << std::endl;
309 
310  event() = std::move(evt);
311  return true;
312 }
313 
314 //_____________________________________________________________________
316  // Hard particles. The first nhard_ lines from hyjets array.
317  // Pythia/Pyquen sub-events (sub-collisions) for a given event
318  // Return T/F if success/failure
319  // Create particles from lujet entries, assign them into vertices and
320  // put the vertices in the GenEvent, for each SubEvent
321  // The SubEvent information is kept by storing indeces of main vertices
322  // of subevents as a vector in GenHIEvent.
323 
324  LogDebug("Hydjet") << " Number of sub events " << nsub_;
325  LogDebug("Hydjet") << " Number of hard events " << hyjpar.njet;
326  LogDebug("Hydjet") << " Number of hard particles " << nhard_;
327  LogDebug("Hydjet") << " Number of soft particles " << nsoft_;
328  LogDebug("Hydjet") << " nhard_ + nsoft_ = " << nhard_ + nsoft_ << " hyjets.nhj = " << hyjets.nhj << endl;
329 
330  int ihy = 0;
331  int isub = -1;
332  int isub_l = -1;
333  int stab = 0;
334 
335  vector<HepMC::GenParticle*> primary_particle(hyjets.nhj);
336  vector<HepMC::GenParticle*> particle(hyjets.nhj);
337 
338  HepMC::GenVertex* sub_vertices = nullptr; // just initialization
339 
340  // contain the last index in for each subevent
341  vector<int> index(nsub_);
342 
343  while (ihy < hyjets.nhj) {
344  isub = std::floor((hyjets.khj[2][ihy] / 50000));
345  int hjoffset = isub * 50000;
346 
347  if (isub != isub_l) {
348  sub_vertices = new HepMC::GenVertex(HepMC::FourVector(0, 0, 0, 0), isub);
349  evt->add_vertex(sub_vertices);
350  if (!evt->signal_process_vertex())
351  evt->set_signal_process_vertex(sub_vertices);
352  index[isub] = ihy - 1;
353  isub_l = isub;
354  }
355 
356  if (convertStatus(hyjets.khj[0][ihy]) == 1)
357  stab++;
358  LogDebug("Hydjet_array") << ihy << " MULTin ev.:" << hyjets.nhj << " SubEv.#" << isub << " Part #" << ihy + 1
359  << ", PDG: " << hyjets.khj[1][ihy] << " (st. " << convertStatus(hyjets.khj[0][ihy])
360  << ") mother=" << hyjets.khj[2][ihy] - (isub * 50000) + index[isub] + 1 << " ("
361  << hyjets.khj[2][ihy] << "), childs ("
362  << hyjets.khj[3][ihy] - (isub * 50000) + index[isub] + 1 << "-"
363  << hyjets.khj[4][ihy] - (isub * 50000) + index[isub] + 1 << "), vtx ("
364  << hyjets.vhj[0][ihy] << "," << hyjets.vhj[1][ihy] << "," << hyjets.vhj[2][ihy] << ") "
365  << std::endl;
366 
367  if (hyjets.khj[2][ihy] == 0) {
368  primary_particle[ihy] = build_hyjet(ihy, ihy + 1);
369  if (!sub_vertices)
370  LogError("Hydjet_array") << "##### HYDJET2: Vertex not initialized!";
371  else
372  sub_vertices->add_particle_out(primary_particle[ihy]);
373  LogDebug("Hydjet_array") << " ---> " << ihy + 1 << std::endl;
374  } else {
375  particle[ihy] = build_hyjet(ihy, ihy + 1);
376  int mid = hyjets.khj[2][ihy] - hjoffset + index[isub];
377  int mid_t = mid;
378  while (((mid + 1) < ihy) && (std::abs(hyjets.khj[1][mid]) < 100) &&
379  (hyjets.khj[3][mid + 1] - hjoffset + index[isub] <= ihy))
380  mid++;
381  if (std::abs(hyjets.khj[1][mid]) < 100)
382  mid = mid_t;
383 
384  HepMC::GenParticle* mother = primary_particle.at(mid);
385  HepMC::GenVertex* prods = build_hyjet_vertex(ihy, isub);
386 
387  if (!mother) {
388  mother = particle[mid];
389  primary_particle[mid] = mother;
390  }
391 
392  HepMC::GenVertex* prod_vertex = mother->end_vertex();
393  if (!prod_vertex) {
394  prod_vertex = prods;
395  prod_vertex->add_particle_in(mother);
396  LogDebug("Hydjet_array") << " <--- " << mid + 1 << std::endl;
397  evt->add_vertex(prod_vertex);
398  prods = nullptr;
399  }
400 
401  prod_vertex->add_particle_out(particle[ihy]);
402  LogDebug("Hydjet_array") << " ---" << mid + 1 << "---> " << ihy + 1 << std::endl;
403  delete prods;
404  }
405  ihy++;
406  }
407  LogDebug("Hydjet_array") << " MULTin ev.:" << hyjets.nhj << ", last index: " << ihy - 1
408  << ", Sub events: " << isub + 1 << ", stable particles: " << stab << std::endl;
409 
410  return true;
411 }
412 
413 //______________________________________________________________
414 bool HydjetHadronizer::call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh) {
415  // initialize hydjet
416 
417  pydatr.mrpy[2] = 1;
418  HYINIT(energy, a, ifb, bmin, bmax, bfix, nh);
419  return true;
420 }
421 
422 //______________________________________________________________
424  // set hydjet options
425 
426  // hydjet running mode mode
427  // kHydroOnly --- nhsel=0 jet production off (pure HYDRO event), nhsel=0
428  // kHydroJets --- nhsle=1 jet production on, jet quenching off (HYDRO+njet*PYTHIA events)
429  // kHydroQJet --- nhsel=2 jet production & jet quenching on (HYDRO+njet*PYQUEN events)
430  // kJetsOnly --- nhsel=3 jet production on, jet quenching off, HYDRO off (njet*PYTHIA events)
431  // kQJetsOnly --- nhsel=4 jet production & jet quenching on, HYDRO off (njet*PYQUEN events)
432 
433  if (hymode_ == "kHydroOnly")
434  hyjpar.nhsel = 0;
435  else if (hymode_ == "kHydroJets")
436  hyjpar.nhsel = 1;
437  else if (hymode_ == "kHydroQJets")
438  hyjpar.nhsel = 2;
439  else if (hymode_ == "kJetsOnly")
440  hyjpar.nhsel = 3;
441  else if (hymode_ == "kQJetsOnly")
442  hyjpar.nhsel = 4;
443  else
444  hyjpar.nhsel = 2;
445 
446  // fraction of soft hydro induced multiplicity
447  hyflow.fpart = fracsoftmult_;
448 
449  // hadron freez-out temperature
450  hyflow.Tf = hadfreeztemp_;
451 
452  // maximum longitudinal collective rapidity
453  hyflow.ylfl = maxlongy_;
454 
455  // maximum transverse collective rapidity
456  hyflow.ytfl = maxtrany_;
457 
458  // shadowing on=1, off=0
459  hyjpar.ishad = shadowingswitch_;
460 
461  // set inelastic nucleon-nucleon cross section
462  hyjpar.sigin = signn_;
463 
464  // angular emitted gluon spectrum selection
465  pyqpar.ianglu = angularspecselector_;
466 
467  // number of active quark flavors in qgp
468  pyqpar.nfu = nquarkflavor_;
469 
470  // initial temperature of QGP
471  pyqpar.T0u = qgpt0_;
472 
473  // proper time of QGP formation
474  pyqpar.tau0u = qgptau0_;
475 
476  // type of medium induced partonic energy loss
478  edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
479  pyqpar.ienglu = 0;
480  } else if (doradiativeenloss_) {
481  edm::LogInfo("HydjetenLoss") << "##### Only RADIATIVE partonic energy loss ON ####";
482  pyqpar.ienglu = 1;
483  } else if (docollisionalenloss_) {
484  edm::LogInfo("HydjetEnLoss") << "##### Only COLLISIONAL partonic energy loss ON ####";
485  pyqpar.ienglu = 2;
486  } else {
487  edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
488  pyqpar.ienglu = 0;
489  }
490  return true;
491 }
492 
493 //_____________________________________________________________________
494 
498 
499  return true;
500 }
501 
502 //_____________________________________________________________________
503 
506  // pythia6Service_->setGeneralParams();
507 
508  // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage
509  const float ra = nuclear_radius();
510  LogInfo("RAScaling") << "Nuclear radius(RA) = " << ra;
511  bmin_ /= ra;
512  bmax_ /= ra;
513  bfixed_ /= ra;
514 
515  // hydjet running options
517  // initialize hydjet
518  LogInfo("HYDJETinAction") << "##### Calling HYINIT(" << comenergy << "," << abeamtarget_ << "," << cflag_ << ","
519  << bmin_ << "," << bmax_ << "," << bfixed_ << "," << nmultiplicity_ << ") ####";
521  return true;
522 }
523 
524 bool HydjetHadronizer::declareStableParticles(const std::vector<int>& _pdg) {
525  std::vector<int> pdg = _pdg;
526  for (size_t i = 0; i < pdg.size(); i++) {
527  int pyCode = pycomp_(pdg[i]);
528  std::ostringstream pyCard;
529  pyCard << "MDCY(" << pyCode << ",1)=0";
530  std::cout << pyCard.str() << std::endl;
531  call_pygive(pyCard.str());
532  }
533  return true;
534 }
535 
536 //________________________________________________________________
538  const double pi = 3.14159265358979;
539  phi0_ = 2. * pi * gen::pyr_(nullptr) - pi;
540  sinphi0_ = sin(phi0_);
541  cosphi0_ = cos(phi0_);
542 }
543 
544 //________________________________________________________________
545 bool HydjetHadronizer::hadronize() { return false; }
546 
547 bool HydjetHadronizer::decay() { return true; }
548 
549 bool HydjetHadronizer::residualDecay() { return true; }
550 
552 
554 
555 const char* HydjetHadronizer::classname() const { return "gen::HydjetHadronizer"; }
#define HYEVNT
Definition: HydjetWrapper.h:20
double bfixed_
fixed impact param (fm); valid only if cflag_=0
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
HepMC::GenEvent * evt
bool docollisionalenloss_
DEFAULT = true.
bool isVtxGenApplied() const
Definition: HepMCProduct.h:39
double pyr_(int *idummy)
bool hydjet_init(const edm::ParameterSet &pset)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T const * product() const
Definition: Handle.h:70
bool call_pygive(const std::string &line)
double npart
Definition: HydjetWrapper.h:46
void add_heavy_ion_rec(HepMC::GenEvent *evt)
#define HYJVER
Definition: HydjetWrapper.h:25
Log< level::Error, false > LogError
double comenergy
collision energy
#define hyfpar
Definition: HydjetWrapper.h:51
double v[5][pyjets_maxn]
const Double_t pi
double nuclear_radius() const
static const std::string kPythia6
#define HYINIT
Definition: HydjetWrapper.h:13
double p[5][pyjets_maxn]
bool embedding_
Switch for embedding mode.
edm::Event & getEDMEvent() const
Definition: EPCuts.h:4
std::string hymode_
Hydjet running mode.
int nsoft_
multiplicity of HYDRO-induced particles in event
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const char * classname() const
unsigned int maxEventsToPrint_
Events to print if verbosity.
uint32_t nh
std::unique_ptr< HepMC::GenEvent > & event()
HepMC::GenParticle * build_hyjet(int index, int barcode)
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
Pythia6Service * pythia6Service_
#define hyjpar
Definition: HydjetWrapper.h:93
int nsub_
number of sub-events
std::vector< double > signalVtx_
Pset double vector to set event signal vertex.
Log< level::Info, false > LogInfo
static const std::string kFortranInstance
int nhard_
multiplicity of PYTHIA(+PYQUEN)-induced particles in event
int pycomp_(int &)
#define pypars
HepMC::GenVertex * build_hyjet_vertex(int i, int id)
#define hyjets
Definition: HydjetWrapper.h:80
bool call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh)
unsigned int pythiaPylistVerbosity_
pythia verbosity; def=1
void setRandomEngine(CLHEP::HepRandomEngine *v)
double abeamtarget_
beam/target atomic mass number
#define pyqpar
HLT enums.
#define hyflow
Definition: HydjetWrapper.h:40
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > src_
double a
Definition: hdecay.h:119
edm::ParameterSet pset_
unsigned int shadowingswitch_
bool doradiativeenloss_
DEFAULT = true.
float x
int angularspecselector_
angular emitted gluon spectrum selection
Log< level::Warning, false > LogWarning
bool declareStableParticles(const std::vector< int > &)
bool rotate_
Switch to rotate event plane.
double phi0_
Event plane angle.
HepMC::FourVector * fVertex_
Event signal vertex.
def move(src, dest)
Definition: eostools.py:511
bool get_particles(HepMC::GenEvent *evt)
#define LogDebug(id)