CMS 3D CMS Logo

Hydjet2Hadronizer.cc
Go to the documentation of this file.
1 
7 #include <TLorentzVector.h>
8 #include <TMath.h>
9 #include <TVector3.h>
10 
12 #include <cmath>
13 #include <fstream>
14 #include <iostream>
15 
22 
25 
26 #include "HepMC/GenEvent.h"
27 #include "HepMC/HeavyIon.h"
28 #include "HepMC/IO_HEPEVT.h"
29 #include "HepMC/PythiaWrapper6_4.h"
30 #include "HepMC/SimpleVector.h"
31 
36 
37 CLHEP::HepRandomEngine *hjRandomEngine;
38 
39 using namespace edm;
40 using namespace std;
41 using namespace gen;
42 
43 int Hydjet2Hadronizer::convertStatusForComponents(int sta, int typ, int pySta) {
44  int st = -1;
45  if (typ == 0) //soft
46  st = 2 - sta;
47  else if (typ == 1)
48  st = convertStatus(pySta);
49 
50  if (st == -1)
51  throw cms::Exception("ConvertStatus") << "Wrong status code!" << endl;
52 
53  if (separateHydjetComponents_) {
54  if (st == 1 && typ == 0)
55  return 6;
56  if (st == 1 && typ == 1)
57  return 7;
58  if (st == 2 && typ == 0)
59  return 16;
60  if (st == 2 && typ == 1)
61  return 17;
62  }
63  return st;
64 }
65 
66 int Hydjet2Hadronizer::convertStatus(int st) {
67  if (st <= 0)
68  return 0;
69  if (st <= 10)
70  return 1;
71  if (st <= 20)
72  return 2;
73  if (st <= 30)
74  return 3;
75  else
76  return -1;
77 }
78 
79 const std::vector<std::string> Hydjet2Hadronizer::theSharedResources = {edm::SharedResourceNames::kPythia6};
80 
81 //____________________________________________________________________________________________
82 Hydjet2Hadronizer::Hydjet2Hadronizer(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
84  rotate_(pset.getParameter<bool>("rotateEventPlane")),
85  evt(nullptr),
86  nsub_(0),
87  nhard_(0),
88  nsoft_(0),
89  phi0_(0.),
90  sinphi0_(0.),
91  cosphi0_(1.),
92  fVertex_(nullptr),
93  pythia6Service_(new Pythia6Service(pset))
94 
95 {
96  fParams.doPrintInfo = false;
97  fParams.allowEmptyEvent = false;
98  fParams.fNevnt = 0; //not used in CMSSW
99  fParams.femb = pset.getParameter<int>("embeddingMode"); //
100  fParams.fSqrtS = pset.getParameter<double>("fSqrtS"); // C.m.s. energy per nucleon pair
101  fParams.fAw = pset.getParameter<double>("fAw"); // Atomic weigth of nuclei, fAw
102  fParams.fIfb = pset.getParameter<int>(
103  "fIfb"); // Flag of type of centrality generation, fBfix (=0 is fixed by fBfix, >0 distributed [fBfmin, fBmax])
104  fParams.fBmin = pset.getParameter<double>("fBmin"); // Minimum impact parameter in units of nuclear radius, fBmin
105  fParams.fBmax = pset.getParameter<double>("fBmax"); // Maximum impact parameter in units of nuclear radius, fBmax
106  fParams.fBfix = pset.getParameter<double>("fBfix"); // Fixed impact parameter in units of nuclear radius, fBfix
107  fParams.fT = pset.getParameter<double>("fT"); // Temperature at chemical freeze-out, fT [GeV]
108  fParams.fMuB = pset.getParameter<double>("fMuB"); // Chemical baryon potential per unit charge, fMuB [GeV]
109  fParams.fMuS = pset.getParameter<double>("fMuS"); // Chemical strangeness potential per unit charge, fMuS [GeV]
110  fParams.fMuC = pset.getParameter<double>(
111  "fMuC"); // Chemical charm potential per unit charge, fMuC [GeV] (used if charm production is turned on)
112  fParams.fMuI3 = pset.getParameter<double>("fMuI3"); // Chemical isospin potential per unit charge, fMuI3 [GeV]
113  fParams.fThFO = pset.getParameter<double>("fThFO"); // Temperature at thermal freeze-out, fThFO [GeV]
114  fParams.fMu_th_pip =
115  pset.getParameter<double>("fMu_th_pip"); // Chemical potential of pi+ at thermal freeze-out, fMu_th_pip [GeV]
116  fParams.fTau = pset.getParameter<double>(
117  "fTau"); // Proper time proper at thermal freeze-out for central collisions, fTau [fm/c]
118  fParams.fSigmaTau = pset.getParameter<double>(
119  "fSigmaTau"); // Duration of emission at thermal freeze-out for central collisions, fSigmaTau [fm/c]
120  fParams.fR = pset.getParameter<double>(
121  "fR"); // Maximal transverse radius at thermal freeze-out for central collisions, fR [fm]
122  fParams.fYlmax =
123  pset.getParameter<double>("fYlmax"); // Maximal longitudinal flow rapidity at thermal freeze-out, fYlmax
124  fParams.fUmax = pset.getParameter<double>(
125  "fUmax"); // Maximal transverse flow rapidity at thermal freeze-out for central collisions, fUmax
126  fParams.frhou2 = pset.getParameter<double>("fRhou2"); //parameter to swich ON/OFF = 0) rhou2
127  fParams.frhou3 = pset.getParameter<double>("fRhou3"); //parameter to swich ON/OFF(0) rhou3
128  fParams.frhou4 = pset.getParameter<double>("fRhou4"); //parameter to swich ON/OFF(0) rhou4
129  fParams.fDelta =
130  pset.getParameter<double>("fDelta"); // Momentum azimuthal anizotropy parameter at thermal freeze-out, fDelta
131  fParams.fEpsilon =
132  pset.getParameter<double>("fEpsilon"); // Spatial azimuthal anisotropy parameter at thermal freeze-out, fEpsilon
133  fParams.fv2 = pset.getParameter<double>("fKeps2"); //parameter to swich ON/OFF(0) epsilon2 fluctuations
134  fParams.fv3 = pset.getParameter<double>("fKeps3"); //parameter to swich ON/OFF(0) epsilon3 fluctuations
135  fParams.fIfDeltaEpsilon = pset.getParameter<double>(
136  "fIfDeltaEpsilon"); // Flag to specify fDelta and fEpsilon values, fIfDeltaEpsilon (=0 user's ones, >=1 calculated)
137  fParams.fDecay =
138  pset.getParameter<int>("fDecay"); // Flag to switch on/off hadron decays, fDecay (=0 decays off, >=1 decays on)
139  fParams.fWeakDecay = pset.getParameter<double>(
140  "fWeakDecay"); // Low decay width threshold fWeakDecay[GeV]: width<fWeakDecay decay off, width>=fDecayWidth decay on; can be used to switch off weak decays
141  fParams.fEtaType = pset.getParameter<double>(
142  "fEtaType"); // Flag to choose longitudinal flow rapidity distribution, fEtaType (=0 uniform, >0 Gaussian with the dispersion Ylmax)
143  fParams.fTMuType = pset.getParameter<double>(
144  "fTMuType"); // Flag to use calculated T_ch, mu_B and mu_S as a function of fSqrtS, fTMuType (=0 user's ones, >0 calculated)
145  fParams.fCorrS = pset.getParameter<double>(
146  "fCorrS"); // Strangeness supression factor gamma_s with fCorrS value (0<fCorrS <=1, if fCorrS <= 0 then it is calculated)
147  fParams.fCharmProd = pset.getParameter<int>(
148  "fCharmProd"); // Flag to include thermal charm production, fCharmProd (=0 no charm production, >=1 charm production)
149  fParams.fCorrC = pset.getParameter<double>(
150  "fCorrC"); // Charmness enhancement factor gamma_c with fCorrC value (fCorrC >0, if fCorrC<0 then it is calculated)
151  fParams.fNhsel = pset.getParameter<int>(
152  "fNhsel"); //Flag to include jet (J)/jet quenching (JQ) and hydro (H) state production, fNhsel (0 H on & J off, 1 H/J on & JQ off, 2 H/J/HQ on, 3 J on & H/JQ off, 4 H off & J/JQ on)
153  fParams.fPyhist = pset.getParameter<int>(
154  "fPyhist"); // Flag to suppress the output of particle history from PYTHIA, fPyhist (=1 only final state particles; =0 full particle history from PYTHIA)
155  fParams.fIshad = pset.getParameter<int>(
156  "fIshad"); // Flag to switch on/off nuclear shadowing, fIshad (0 shadowing off, 1 shadowing on)
157  fParams.fPtmin =
158  pset.getParameter<double>("fPtmin"); // Minimal pt of parton-parton scattering in PYTHIA event, fPtmin [GeV/c]
159  fParams.fT0 = pset.getParameter<double>(
160  "fT0"); // Initial QGP temperature for central Pb+Pb collisions in mid-rapidity, fT0 [GeV]
161  fParams.fTau0 = pset.getParameter<double>("fTau0"); // Proper QGP formation time in fm/c, fTau0 (0.01<fTau0<10)
162  fParams.fNf = pset.getParameter<int>("fNf"); // Number of active quark flavours in QGP, fNf (0, 1, 2 or 3)
163  fParams.fIenglu = pset.getParameter<int>(
164  "fIenglu"); // Flag to fix type of partonic energy loss, fIenglu (0 radiative and collisional loss, 1 radiative loss only, 2 collisional loss only)
165  fParams.fIanglu = pset.getParameter<int>(
166  "fIanglu"); // Flag to fix type of angular distribution of in-medium emitted gluons, fIanglu (0 small-angular, 1 wide-angular, 2 collinear).
167 
168  edm::FileInPath f1("externals/hydjet2/particles.data");
169  strcpy(fParams.partDat, (f1.fullPath()).c_str());
170 
171  edm::FileInPath f2("externals/hydjet2/tabledecay.txt");
172  strcpy(fParams.tabDecay, (f2.fullPath()).c_str());
173 
174  fParams.fPythiaTune = false;
175 
176  if (pset.exists("signalVtx"))
177  signalVtx_ = pset.getUntrackedParameter<std::vector<double>>("signalVtx");
178 
179  if (signalVtx_.size() == 4) {
180  if (!fVertex_)
181  fVertex_ = new HepMC::FourVector();
182  LogDebug("EventSignalVertex") << "Setting event signal vertex "
183  << " x = " << signalVtx_.at(0) << " y = " << signalVtx_.at(1)
184  << " z= " << signalVtx_.at(2) << " t = " << signalVtx_.at(3) << endl;
185  fVertex_->set(signalVtx_.at(0), signalVtx_.at(1), signalVtx_.at(2), signalVtx_.at(3));
186  }
187 
188  // PYLIST Verbosity Level
189  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
190  pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
191  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_;
192  //Max number of events printed on verbosity level
193  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint", 0);
194  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_;
195  if (fParams.femb == 1) {
196  fParams.fIfb = 0;
197  src_ = iC.consumes<CrossingFrame<edm::HepMCProduct>>(
198  pset.getUntrackedParameter<edm::InputTag>("backgroundLabel", edm::InputTag("mix", "generatorSmeared")));
199  }
200 
201  separateHydjetComponents_ = pset.getUntrackedParameter<bool>("separateHydjetComponents", false);
202 }
203 //__________________________________________________________________________________________
205  call_pystat(1);
206  delete pythia6Service_;
207 }
208 
209 //_____________________________________________________________________
210 void Hydjet2Hadronizer::doSetRandomEngine(CLHEP::HepRandomEngine *v) {
212  hjRandomEngine = v;
213 }
214 
215 //______________________________________________________________________________________________________
219 
220  fParams.fSeed = hjRandomEngine->CLHEP::HepRandomEngine::getSeed();
221  LogInfo("Hydjet2Hadronizer|GenSeed") << "Seed for random number generation: "
222  << hjRandomEngine->CLHEP::HepRandomEngine::getSeed();
223 
224  return kTRUE;
225 }
226 
227 //______________________________________________________________________________________________________
230 
231  // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage
232  const double ra = nuclear_radius();
233  LogInfo("Hydjet2Hadronizer|RAScaling") << "Nuclear radius(RA) = " << ra;
234  fParams.fBmin /= ra;
235  fParams.fBmax /= ra;
236  fParams.fBfix /= ra;
237 
238  hj2 = new Hydjet2(fParams);
239 
240  return kTRUE;
241 }
242 
243 //__________________________________________________________________________________________
246 
247  // generate single event
248  if (fParams.femb == 1) {
249  const edm::Event &e = getEDMEvent();
250  HepMC::GenVertex *genvtx = nullptr;
251  const HepMC::GenEvent *inev = nullptr;
253  e.getByToken(src_, cf);
255  if (mix.size() < 1) {
256  throw cms::Exception("MatchVtx") << "Mixing has " << mix.size() << " sub-events, should have been at least 1"
257  << endl;
258  }
259  const HepMCProduct &bkg = mix.getObject(0);
260  if (!(bkg.isVtxGenApplied())) {
261  throw cms::Exception("MatchVtx") << "Input background does not have smeared vertex!" << endl;
262  } else {
263  inev = bkg.GetEvent();
264  }
265 
266  genvtx = inev->signal_process_vertex();
267 
268  if (!genvtx)
269  throw cms::Exception("MatchVtx") << "Input background does not have signal process vertex!" << endl;
270 
271  double aX, aY, aZ, aT;
272 
273  aX = genvtx->position().x();
274  aY = genvtx->position().y();
275  aZ = genvtx->position().z();
276  aT = genvtx->position().t();
277 
278  if (!fVertex_) {
279  fVertex_ = new HepMC::FourVector();
280  }
281  LogInfo("MatchVtx") << " setting vertex "
282  << " aX " << aX << " aY " << aY << " aZ " << aZ << " aT " << aT << endl;
283  fVertex_->set(aX, aY, aZ, aT);
284 
285  const HepMC::HeavyIon *hi = inev->heavy_ion();
286 
287  if (hi) {
288  fParams.fBfix = (hi->impact_parameter()) / nuclear_radius();
289  phi0_ = hi->event_plane_angle();
290  sinphi0_ = sin(phi0_);
291  cosphi0_ = cos(phi0_);
292  } else {
293  LogWarning("EventEmbedding") << "Background event does not have heavy ion record!";
294  }
295 
296  } else if (rotate_)
297  rotateEvtPlane();
298 
299  nsoft_ = 0;
300  nhard_ = 0;
301 
302  // generate one HYDJET event
303  int ntry = 0;
304 
305  while (nsoft_ == 0 && nhard_ == 0) {
306  if (ntry > 100) {
307  LogError("Hydjet2EmptyEvent") << "##### HYDJET2: No Particles generated, Number of tries =" << ntry;
308  // Throw an exception. Use the EventCorruption exception since it maps onto SkipEvent
309  // which is what we want to do here.
310  std::ostringstream sstr;
311  sstr << "Hydjet2HadronizerProducer: No particles generated after " << ntry << " tries.\n";
312  edm::Exception except(edm::errors::EventCorruption, sstr.str());
313  throw except;
314  } else {
315  hj2->GenerateEvent(fParams.fBfix);
316 
317  if (hj2->IsEmpty()) {
318  continue;
319  }
320 
321  nsoft_ = hj2->GetNhyd();
322  nsub_ = hj2->GetNjet();
323  nhard_ = hj2->GetNpyt();
324 
325  //100 trys
326  ++ntry;
327  }
328  }
329 
330  if (ev == 0) {
331  Sigin = hj2->GetSigin();
332  Sigjet = hj2->GetSigjet();
333  }
334  ev = true;
335 
336  if (fParams.fNhsel < 3)
337  nsub_++;
338 
339  // event information
340  std::unique_ptr<HepMC::GenEvent> evt = std::make_unique<HepMC::GenEvent>();
341  std::unique_ptr<edm::HepMCProduct> HepMCEvt = std::make_unique<edm::HepMCProduct>();
342 
343  if (nhard_ > 0 || nsoft_ > 0)
344  get_particles(evt.get());
345 
346  evt->set_signal_process_id(pypars.msti[0]); // type of the process
347  evt->set_event_scale(pypars.pari[16]); // Q^2
348  add_heavy_ion_rec(evt.get());
349 
350  if (fVertex_) {
351  // generate new vertex & apply the shift
352  // Copy the HepMC::GenEvent
353  HepMCEvt = std::make_unique<edm::HepMCProduct>(evt.get());
354  HepMCEvt->applyVtxGen(fVertex_);
355  evt = std::make_unique<HepMC::GenEvent>(*HepMCEvt->GetEvent());
356  }
357 
358  HepMC::HEPEVT_Wrapper::check_hepevt_consistency();
359  LogDebug("HEPEVT_info") << "Ev numb: " << HepMC::HEPEVT_Wrapper::event_number()
360  << " Entries number: " << HepMC::HEPEVT_Wrapper::number_entries() << " Max. entries "
361  << HepMC::HEPEVT_Wrapper::max_number_entries() << std::endl;
362 
363  event() = std::move(evt);
364  return kTRUE;
365 }
366 
367 //________________________________________________________________
368 bool Hydjet2Hadronizer::declareStableParticles(const std::vector<int> &_pdg) {
369  std::vector<int> pdg = _pdg;
370  for (size_t i = 0; i < pdg.size(); i++) {
371  int pyCode = pycomp_(pdg[i]);
372  std::ostringstream pyCard;
373  pyCard << "MDCY(" << pyCode << ",1)=0";
374  std::cout << pyCard.str() << std::endl;
375  call_pygive(pyCard.str());
376  }
377  return true;
378 }
379 //________________________________________________________________
380 bool Hydjet2Hadronizer::hadronize() { return false; }
381 bool Hydjet2Hadronizer::decay() { return true; }
382 bool Hydjet2Hadronizer::residualDecay() { return true; }
385 const char *Hydjet2Hadronizer::classname() const { return "gen::Hydjet2Hadronizer"; }
386 
387 //________________________________________________________________
389  const double pi = 3.14159265358979;
390  phi0_ = 2. * pi * gen::pyr_(nullptr) - pi;
391  sinphi0_ = sin(phi0_);
392  cosphi0_ = cos(phi0_);
393 }
394 
395 //_____________________________________________________________________
397  LogDebug("Hydjet2") << " Number of sub events " << nsub_;
398  LogDebug("Hydjet2") << " Number of hard events " << hj2->GetNjet();
399  LogDebug("Hydjet2") << " Number of hard particles " << nhard_;
400  LogDebug("Hydjet2") << " Number of soft particles " << nsoft_;
401  LogDebug("Hydjet2") << " nhard_ + nsoft_ = " << nhard_ + nsoft_ << " Ntot = " << hj2->GetNtot() << endl;
402 
403  int ihy = 0;
404  int isub_l = -1;
405  int stab = 0;
406 
407  vector<HepMC::GenParticle *> particle(hj2->GetNtot());
408  HepMC::GenVertex *sub_vertices = nullptr;
409 
410  while (ihy < hj2->GetNtot()) {
411  if ((hj2->GetiJet().at(ihy)) != isub_l) {
412  sub_vertices = new HepMC::GenVertex(HepMC::FourVector(0, 0, 0, 0), hj2->GetiJet().at(ihy));
413  evt->add_vertex(sub_vertices);
414  if (!evt->signal_process_vertex())
415  evt->set_signal_process_vertex(sub_vertices);
416  isub_l = hj2->GetiJet().at(ihy);
417  }
418 
420  (hj2->GetFinal()).at(ihy), (hj2->GetType()).at(ihy), (hj2->GetPythiaStatus().at(ihy)))) == 1)
421  stab++;
422  LogDebug("Hydjet2_array") << ihy << " MULTin ev.:" << hj2->GetNtot() << " SubEv.#" << hj2->GetiJet().at(ihy)
423  << " Part #" << ihy + 1 << ", PDG: " << hj2->GetPdg().at(ihy) << " (st. "
424  << convertStatus(hj2->GetPythiaStatus().at(ihy))
425  << ") mother=" << hj2->GetMotherIndex().at(ihy) + 1 << ", childs ("
426  << hj2->GetFirstDaughterIndex().at(ihy) + 1 << "-"
427  << hj2->GetLastDaughterIndex().at(ihy) + 1 << "), vtx (" << hj2->GetX().at(ihy) << ","
428  << hj2->GetY().at(ihy) << "," << hj2->GetZ().at(ihy) << ") " << std::endl;
429 
430  if ((hj2->GetMotherIndex().at(ihy)) <= 0) {
431  particle.at(ihy) = build_hyjet2(ihy, ihy + 1);
432  if (!sub_vertices)
433  LogError("Hydjet2_array") << "##### HYDJET2: Vertex not initialized!";
434  else
435  sub_vertices->add_particle_out(particle.at(ihy));
436  LogDebug("Hydjet2_array") << " ---> " << ihy + 1 << std::endl;
437  } else {
438  particle.at(ihy) = build_hyjet2(ihy, ihy + 1);
439  int mid = hj2->GetMotherIndex().at(ihy);
440 
441  while (((mid + 1) < ihy) && (std::abs(hj2->GetPdg().at(mid)) < 100) &&
442  ((hj2->GetFirstDaughterIndex().at(mid + 1)) <= ihy)) {
443  mid++;
444  LogDebug("Hydjet2_array") << "======== MID changed to " << mid
445  << " ======== PDG(mid) = " << hj2->GetPdg().at(mid) << std::endl;
446  }
447 
448  if (std::abs(hj2->GetPdg().at(mid)) < 100) {
449  mid = hj2->GetMotherIndex().at(ihy);
450  LogDebug("Hydjet2_array") << "======== MID changed BACK to " << mid
451  << " ======== PDG(mid) = " << hj2->GetPdg().at(mid) << std::endl;
452  }
453 
454  HepMC::GenParticle *mother = particle.at(mid);
455  HepMC::GenVertex *prod_vertex = mother->end_vertex();
456 
457  if (!prod_vertex) {
458  prod_vertex = build_hyjet2_vertex(ihy, (hj2->GetiJet().at(ihy)));
459  prod_vertex->add_particle_in(mother);
460  LogDebug("Hydjet2_array") << " <--- " << mid + 1 << std::endl;
461  evt->add_vertex(prod_vertex);
462  }
463  prod_vertex->add_particle_out(particle.at(ihy));
464  LogDebug("Hydjet2_array") << " ---" << mid + 1 << "---> " << ihy + 1 << std::endl;
465  }
466  ihy++;
467  }
468 
469  LogDebug("Hydjet2_array") << " MULTin ev.:" << hj2->GetNtot() << ", last index: " << ihy - 1
470  << ", stable particles: " << stab << std::endl;
471  return kTRUE;
472 }
473 
474 //___________________________________________________________________
476  // Build particle object corresponding to index in hyjets (soft+hard)
477  double px0 = (hj2->GetPx()).at(index);
478  double py0 = (hj2->GetPy()).at(index);
479 
480  double px = px0 * cosphi0_ - py0 * sinphi0_;
481  double py = py0 * cosphi0_ + px0 * sinphi0_;
482 
484  HepMC::FourVector(px, // px
485  py, // py
486  (hj2->GetPz()).at(index), // pz
487  (hj2->GetE()).at(index)), // E
488  (hj2->GetPdg()).at(index), // id
490  (hj2->GetFinal()).at(index), (hj2->GetType()).at(index), (hj2->GetPythiaStatus()).at(index)) // status
491  );
492 
493  p->suggest_barcode(barcode);
494  return p;
495 }
496 
497 //___________________________________________________________________
498 HepMC::GenVertex *Hydjet2Hadronizer::build_hyjet2_vertex(int i, int id) {
499  // build verteces for the hyjets stored events
500  double x0 = (hj2->GetX()).at(i);
501  double y0 = (hj2->GetY()).at(i);
502 
503  // convert to mm (as in PYTHIA6)
504  const double fm_to_mm = 1e-12;
505  double x = fm_to_mm * (x0 * cosphi0_ - y0 * sinphi0_);
506  double y = fm_to_mm * (y0 * cosphi0_ + x0 * sinphi0_);
507  double z = fm_to_mm * (hj2->GetZ()).at(i);
508  double t = fm_to_mm * (hj2->GetT()).at(i);
509 
510  HepMC::GenVertex *vertex = new HepMC::GenVertex(HepMC::FourVector(x, y, z, t), id);
511 
512  return vertex;
513 }
514 
515 //_____________________________________________________________________
517  // heavy ion record in the final CMSSW Event
518  int nproj = static_cast<int>((hj2->GetNpart()) / 2);
519  int ntarg = static_cast<int>((hj2->GetNpart()) - nproj);
520 
521  HepMC::HeavyIon *hi = new HepMC::HeavyIon(nsub_, // Ncoll_hard/N of SubEvents
522  nproj, // Npart_proj
523  ntarg, // Npart_targ
524  hj2->GetNbcol(), // Ncoll
525  0, // spectator_neutrons
526  0, // spectator_protons
527  0, // N_Nwounded_collisions
528  0, // Nwounded_N_collisions
529  0, // Nwounded_Nwounded_collisions
530  hj2->GetBgen() * nuclear_radius(), // impact_parameter in [fm]
531  phi0_, // event_plane_angle
532  hj2->GetPsiv3(), // eccentricity <<<---- psi for v3!!!
533  Sigin // sigma_inel_NN
534  );
535 
536  evt->set_heavy_ion(*hi);
537  delete hi;
538 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool declareStableParticles(const std::vector< int > &)
edm::ParameterSet pset
unsigned int maxEventsToPrint_
InitialParamsHydjet_t fParams
bool isVtxGenApplied() const
Definition: HepMCProduct.h:39
void add_heavy_ion_rec(HepMC::GenEvent *evt)
double pyr_(int *idummy)
unsigned int pythiaPylistVerbosity_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool exists(std::string const &parameterName) const
checks if a parameter exists
T const * product() const
Definition: Handle.h:70
bool call_pygive(const std::string &line)
HepMC::FourVector * fVertex_
std::vector< double > signalVtx_
Log< level::Error, false > LogError
double v[5][pyjets_maxn]
const Double_t pi
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > src_
const char * classname() const
T getUntrackedParameter(std::string const &, T const &) const
static const std::string kPythia6
HepMC::GenVertex * build_hyjet2_vertex(int i, int id)
double nuclear_radius() const
Pythia6Service * pythia6Service_
double p[5][pyjets_maxn]
edm::Event & getEDMEvent() const
CLHEP::HepRandomEngine * hjRandomEngine
Interface to the HYDJET++ (Hydjet2) generator (since core v. 2.4.3), produces HepMC events...
bool get_particles(HepMC::GenEvent *evt)
Definition: EPCuts.h:4
int convertStatusForComponents(int, int, int)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::unique_ptr< HepMC::GenEvent > & event()
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
Log< level::Info, false > LogInfo
int pycomp_(int &)
#define pypars
HepMC::GenParticle * build_hyjet2(int index, int barcode)
void setRandomEngine(CLHEP::HepRandomEngine *v)
HLT enums.
float x
Log< level::Warning, false > LogWarning
HepMC::GenEvent * evt
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)