CMS 3D CMS Logo

Hydjet2Hadronizer.cc
Go to the documentation of this file.
1 /*
2 
3 ##########################
4 # Hydjet2 #
5 # version: 2.2.0 patch1 #
6 ##########################
7 
8  Interface to the HYDJET++ generator, produces HepMC events
9 
10  Author: Andrey Belyaev (Andrey.Belyaev@cern.ch)
11 
12  Hydjet2Hadronizer is the modified InitialStateHydjet
13 
14  HYDJET++
15  version 2.2:
16  InitialStateHydjet is the modified InitialStateBjorken
17  The high-pt part related with PYTHIA-PYQUEN is included
18  InitialStateBjorken (FASTMC) was used.
19 
20 
21  InitialStateBjorken
22  version 2.0:
23  Ludmila Malinina malinina@lav01.sinp.msu.ru, SINP MSU/Moscow and JINR/Dubna
24  Ionut Arsene i.c.arsene@fys.uio.no, Oslo University
25  June 2007
26 
27  version 1.0:
28  Nikolai Amelin, Ludmila Malinina, Timur Pocheptsov (C) JINR/Dubna
29  amelin@sunhe.jinr.ru, malinina@sunhe.jinr.ru, pocheptsov@sunhe.jinr.ru
30  November. 2, 2005
31 */
32 
33 //expanding localy equilibated fireball with volume hadron radiation
34 //thermal part: Blast wave model, Bjorken-like parametrization
35 //high-pt: PYTHIA + jet quenching model PYQUEN
36 
37 #include <TLorentzVector.h>
38 #include <TVector3.h>
39 #include <TMath.h>
40 
51 
52 #include <iostream>
53 #include <fstream>
54 #include <cmath>
55 #include "boost/lexical_cast.hpp"
56 
63 
66 
67 #include "HepMC/PythiaWrapper6_4.h"
68 #include "HepMC/GenEvent.h"
69 #include "HepMC/HeavyIon.h"
70 #include "HepMC/SimpleVector.h"
71 
76 
78 
79 extern "C" void hyevnt_();
80 extern "C" void myini_();
81 
82 extern HYIPARCommon HYIPAR;
83 extern HYFPARCommon HYFPAR;
84 extern HYJPARCommon HYJPAR;
85 extern HYPARTCommon HYPART;
86 
87 using namespace edm;
88 using namespace std;
89 using namespace gen;
90 
92 
93 // definition of the static member fLastIndex
95 bool ev = false;
96 namespace {
97  int convertStatusForComponents(int sta, int typ) {
98  if (sta == 1 && typ == 0)
99  return 6;
100  if (sta == 1 && typ == 1)
101  return 7;
102  if (sta == 2 && typ == 0)
103  return 16;
104  if (sta == 2 && typ == 1)
105  return 17;
106 
107  else
108  return sta;
109  }
110 
111  int convertStatus(int st) {
112  if (st <= 0)
113  return 0;
114  if (st <= 10)
115  return 1;
116  if (st <= 20)
117  return 2;
118  if (st <= 30)
119  return 3;
120 
121  else
122  return st;
123  }
124 } // namespace
125 
126 const std::vector<std::string> Hydjet2Hadronizer::theSharedResources = {edm::SharedResourceNames::kPythia6,
128 
129 //____________________________________________________________________________________________
130 Hydjet2Hadronizer::Hydjet2Hadronizer(const edm::ParameterSet& pset)
131  : BaseHadronizer(pset),
132  fSqrtS(pset.getParameter<double>("fSqrtS")), // C.m.s. energy per nucleon pair
133  fAw(pset.getParameter<double>("fAw")), // Atomic weigth of nuclei, fAw
134  fIfb(pset.getParameter<int>(
135  "fIfb")), // Flag of type of centrality generation, fBfix (=0 is fixed by fBfix, >0 distributed [fBfmin, fBmax])
136  fBmin(pset.getParameter<double>("fBmin")), // Minimum impact parameter in units of nuclear radius, fBmin
137  fBmax(pset.getParameter<double>("fBmax")), // Maximum impact parameter in units of nuclear radius, fBmax
138  fBfix(pset.getParameter<double>("fBfix")), // Fixed impact parameter in units of nuclear radius, fBfix
139  fT(pset.getParameter<double>("fT")), // Temperature at chemical freeze-out, fT [GeV]
140  fMuB(pset.getParameter<double>("fMuB")), // Chemical baryon potential per unit charge, fMuB [GeV]
141  fMuS(pset.getParameter<double>("fMuS")), // Chemical strangeness potential per unit charge, fMuS [GeV]
142  fMuC(pset.getParameter<double>(
143  "fMuC")), // Chemical charm potential per unit charge, fMuC [GeV] (used if charm production is turned on)
144  fMuI3(pset.getParameter<double>("fMuI3")), // Chemical isospin potential per unit charge, fMuI3 [GeV]
145  fThFO(pset.getParameter<double>("fThFO")), // Temperature at thermal freeze-out, fThFO [GeV]
146  fMu_th_pip(pset.getParameter<double>(
147  "fMu_th_pip")), // Chemical potential of pi+ at thermal freeze-out, fMu_th_pip [GeV]
148  fTau(pset.getParameter<double>(
149  "fTau")), // Proper time proper at thermal freeze-out for central collisions, fTau [fm/c]
150  fSigmaTau(pset.getParameter<double>(
151  "fSigmaTau")), // Duration of emission at thermal freeze-out for central collisions, fSigmaTau [fm/c]
152  fR(pset.getParameter<double>(
153  "fR")), // Maximal transverse radius at thermal freeze-out for central collisions, fR [fm]
154  fYlmax(pset.getParameter<double>("fYlmax")), // Maximal longitudinal flow rapidity at thermal freeze-out, fYlmax
155  fUmax(pset.getParameter<double>(
156  "fUmax")), // Maximal transverse flow rapidity at thermal freeze-out for central collisions, fUmax
157  fDelta(pset.getParameter<double>(
158  "fDelta")), // Momentum azimuthal anizotropy parameter at thermal freeze-out, fDelta
159  fEpsilon(pset.getParameter<double>(
160  "fEpsilon")), // Spatial azimuthal anisotropy parameter at thermal freeze-out, fEpsilon
161  fIfDeltaEpsilon(pset.getParameter<double>(
162  "fIfDeltaEpsilon")), // Flag to specify fDelta and fEpsilon values, fIfDeltaEpsilon (=0 user's ones, >=1 calculated)
163  fDecay(pset.getParameter<int>(
164  "fDecay")), // Flag to switch on/off hadron decays, fDecay (=0 decays off, >=1 decays on)
165  fWeakDecay(pset.getParameter<double>(
166  "fWeakDecay")), // Low decay width threshold fWeakDecay[GeV]: width<fWeakDecay decay off, width>=fDecayWidth decay on; can be used to switch off weak decays
167  fEtaType(pset.getParameter<double>(
168  "fEtaType")), // Flag to choose longitudinal flow rapidity distribution, fEtaType (=0 uniform, >0 Gaussian with the dispersion Ylmax)
169  fTMuType(pset.getParameter<double>(
170  "fTMuType")), // Flag to use calculated T_ch, mu_B and mu_S as a function of fSqrtS, fTMuType (=0 user's ones, >0 calculated)
171  fCorrS(pset.getParameter<double>(
172  "fCorrS")), // Strangeness supression factor gamma_s with fCorrS value (0<fCorrS <=1, if fCorrS <= 0 then it is calculated)
173  fCharmProd(pset.getParameter<int>(
174  "fCharmProd")), // Flag to include thermal charm production, fCharmProd (=0 no charm production, >=1 charm production)
175  fCorrC(pset.getParameter<double>(
176  "fCorrC")), // Charmness enhancement factor gamma_c with fCorrC value (fCorrC >0, if fCorrC<0 then it is calculated)
177  fNhsel(pset.getParameter<int>(
178  "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)
179  fPyhist(pset.getParameter<int>(
180  "fPyhist")), // Flag to suppress the output of particle history from PYTHIA, fPyhist (=1 only final state particles; =0 full particle history from PYTHIA)
181  fIshad(pset.getParameter<int>(
182  "fIshad")), // Flag to switch on/off nuclear shadowing, fIshad (0 shadowing off, 1 shadowing on)
183  fPtmin(pset.getParameter<double>(
184  "fPtmin")), // Minimal pt of parton-parton scattering in PYTHIA event, fPtmin [GeV/c]
185  fT0(pset.getParameter<double>(
186  "fT0")), // Initial QGP temperature for central Pb+Pb collisions in mid-rapidity, fT0 [GeV]
187  fTau0(pset.getParameter<double>("fTau0")), // Proper QGP formation time in fm/c, fTau0 (0.01<fTau0<10)
188  fNf(pset.getParameter<int>("fNf")), // Number of active quark flavours in QGP, fNf (0, 1, 2 or 3)
189  fIenglu(pset.getParameter<int>(
190  "fIenglu")), // Flag to fix type of partonic energy loss, fIenglu (0 radiative and collisional loss, 1 radiative loss only, 2 collisional loss only)
191  fIanglu(pset.getParameter<int>(
192  "fIanglu")), // Flag to fix type of angular distribution of in-medium emitted gluons, fIanglu (0 small-angular, 1 wide-angular, 2 collinear).
193 
194  embedding_(pset.getParameter<bool>("embeddingMode")),
195  rotate_(pset.getParameter<bool>("rotateEventPlane")),
196  evt(nullptr),
197  nsub_(0),
198  nhard_(0),
199  nsoft_(0),
200  phi0_(0.),
201  sinphi0_(0.),
202  cosphi0_(1.),
203  pythia6Service_(new Pythia6Service(pset))
204 
205 {
206  // constructor
207  // PYLIST Verbosity Level
208  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
209  pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
210  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_;
211  //Max number of events printed on verbosity level
212  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint", 0);
213  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_;
214  if (embedding_)
215  src_ = pset.getParameter<edm::InputTag>("backgroundLabel");
216 }
217 
218 //__________________________________________________________________________________________
220  // destructor
221  call_pystat(1);
222  delete pythia6Service_;
223 }
224 
225 //_____________________________________________________________________
226 void Hydjet2Hadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v) {
228  hjRandomEngine = v;
229 }
230 
231 //______________________________________________________________________________________________________
235 
236  SERVICE.iseed_fromC = hjRandomEngine->CLHEP::HepRandomEngine::getSeed();
237  LogInfo("Hydjet2Hadronizer|GenSeed") << "Seed for random number generation: "
238  << hjRandomEngine->CLHEP::HepRandomEngine::getSeed();
239 
240  fNPartTypes = 0; //counter of hadron species
241 
242  return kTRUE;
243 }
244 
245 //______________________________________________________________________________________________________
248 
249  // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage
250  const float ra = nuclear_radius();
251  LogInfo("Hydjet2Hadronizer|RAScaling") << "Nuclear radius(RA) = " << ra;
252  fBmin /= ra;
253  fBmax /= ra;
254  fBfix /= ra;
255 
256  //check and redefine input parameters
257  if (fTMuType > 0 && fSqrtS > 2.24) {
258  if (fSqrtS < 2.24) {
259  LogError("Hydjet2Hadronizer|sqrtS") << "SqrtS<2.24 not allowed with fTMuType>0";
260  return false;
261  }
262 
263  //sqrt(s) = 2.24 ==> T_kin = 0.8 GeV
264  //see J. Cleymans, H. Oeschler, K. Redlich,S. Wheaton, Phys Rev. C73 034905 (2006)
265  fMuB = 1.308 / (1. + fSqrtS * 0.273);
266  fT = 0.166 - 0.139 * fMuB * fMuB - 0.053 * fMuB * fMuB * fMuB * fMuB;
267  fMuI3 = 0.;
268  fMuS = 0.;
269 
270  //create strange potential object and set strangeness density 0
272  psp->SetBaryonPotential(fMuB);
273  psp->SetTemperature(fT);
274 
275  //compute strangeness potential
276  if (fMuB > 0.01)
278  LogInfo("Hydjet2Hadronizer|Strange") << "fMuS = " << fMuS;
279 
280  //if user choose fYlmax larger then allowed by kinematics at the specified beam energy sqrt(s)
281  if (fYlmax > TMath::Log(fSqrtS / 0.94)) {
282  LogError("Hydjet2Hadronizer|Ylmax") << "fYlmax more then TMath::Log(fSqrtS vs 0.94)!!! ";
283  return false;
284  }
285 
286  if (fCorrS <= 0.) {
287  //see F. Becattini, J. Mannien, M. Gazdzicki, Phys Rev. C73 044905 (2006)
288  fCorrS = 1. - 0.386 * TMath::Exp(-1.23 * fT / fMuB);
289  LogInfo("Hydjet2Hadronizer|Strange")
290  << "The phenomenological f-la F. Becattini et al. PRC73 044905 (2006) for CorrS was used." << endl
291  << "Strangeness suppression parameter = " << fCorrS;
292  }
293  LogInfo("Hydjet2Hadronizer|Strange")
294  << "The phenomenological f-la J. Cleymans et al. PRC73 034905 (2006) for Tch mu_B was used." << endl
295  << "The simulation will be done with the calculated parameters:" << endl
296  << "Baryon chemical potential = " << fMuB << " [GeV]" << endl
297  << "Strangeness chemical potential = " << fMuS << " [GeV]" << endl
298  << "Isospin chemical potential = " << fMuI3 << " [GeV]" << endl
299  << "Strangeness suppression parameter = " << fCorrS << endl
300  << "Eta_max = " << fYlmax;
301  }
302 
303  LogInfo("Hydjet2Hadronizer|Param") << "Used eta_max = " << fYlmax << endl
304  << "maximal allowed eta_max TMath::Log(fSqrtS/0.94)= "
305  << TMath::Log(fSqrtS / 0.94);
306 
307  //initialisation of high-pt part
308  HYJPAR.nhsel = fNhsel;
309  HYJPAR.ptmin = fPtmin;
310  HYJPAR.ishad = fIshad;
311  HYJPAR.iPyhist = fPyhist;
312  HYIPAR.bminh = fBmin;
313  HYIPAR.bmaxh = fBmax;
314  HYIPAR.AW = fAw;
315 
316  HYPYIN.ifb = fIfb;
317  HYPYIN.bfix = fBfix;
318  HYPYIN.ene = fSqrtS;
319 
320  PYQPAR.T0 = fT0;
321  PYQPAR.tau0 = fTau0;
322  PYQPAR.nf = fNf;
323  PYQPAR.ienglu = fIenglu;
324  PYQPAR.ianglu = fIanglu;
325  myini_();
326 
327  // calculation of multiplicities of different particle species
328  // according to the grand canonical approach
329  GrandCanonical gc(15, fT, fMuB, fMuS, fMuI3, fMuC);
330  GrandCanonical gc_ch(15, fT, fMuB, fMuS, fMuI3, fMuC);
331  GrandCanonical gc_pi_th(15, fThFO, 0., 0., fMu_th_pip, fMuC);
332  GrandCanonical gc_th_0(15, fThFO, 0., 0., 0., 0.);
333 
334  // std::ofstream outMult("densities.txt");
335  // outMult<<"encoding particle density chemical potential "<<std::endl;
336 
337  double Nocth = 0; //open charm
338  double NJPsith = 0; //JPsi
339 
340  //effective volume for central
341  double dYl = 2 * fYlmax; //uniform distr. [-Ylmax; Ylmax]
342  if (fEtaType > 0)
343  dYl = TMath::Sqrt(2 * TMath::Pi()) * fYlmax; //Gaussian distr.
344  fVolEff = 2 * TMath::Pi() * fTau * dYl * (fR * fR) / TMath::Power((fUmax), 2) *
345  ((fUmax)*TMath::SinH((fUmax)) - TMath::CosH((fUmax)) + 1);
346  LogInfo("Hydjet2Hadronizer|Param") << "central Effective volume = " << fVolEff << " [fm^3]";
347 
348  double particleDensity_pi_ch = 0;
349  double particleDensity_pi_th = 0;
350  // double particleDensity_th_0=0;
351 
352  if (fThFO != fT && fThFO > 0) {
353  GrandCanonical gc_ch(15, fT, fMuB, fMuS, fMuI3, fMuC);
354  GrandCanonical gc_pi_th(15, fThFO, 0., 0., fMu_th_pip, fMuC);
355  GrandCanonical gc_th_0(15, fThFO, 0., 0., 0., 0.);
356  particleDensity_pi_ch = gc_ch.ParticleNumberDensity(fDatabase->GetPDGParticle(211));
357  particleDensity_pi_th = gc_pi_th.ParticleNumberDensity(fDatabase->GetPDGParticle(211));
358  }
359 
360  for (int particleIndex = 0; particleIndex < fDatabase->GetNParticles(); particleIndex++) {
361  ParticlePDG* currParticle = fDatabase->GetPDGParticleByIndex(particleIndex);
362  int encoding = currParticle->GetPDG();
363 
364  //strangeness supression
365  double gammaS = 1;
366  int S = int(currParticle->GetStrangeness());
367  if (encoding == 333)
368  S = 2;
369  if (fCorrS < 1. && S != 0)
370  gammaS = TMath::Power(fCorrS, -TMath::Abs(S));
371 
372  //average densities
373  double particleDensity = gc.ParticleNumberDensity(currParticle) / gammaS;
374 
375  //compute chemical potential for single f.o. mu==mu_ch
376  double mu = fMuB * int(currParticle->GetBaryonNumber()) + fMuS * int(currParticle->GetStrangeness()) +
377  fMuI3 * int(currParticle->GetElectricCharge()) + fMuC * int(currParticle->GetCharmness());
378 
379  //thermal f.o.
380  if (fThFO != fT && fThFO > 0) {
381  double particleDensity_ch = gc_ch.ParticleNumberDensity(currParticle);
382  double particleDensity_th_0 = gc_th_0.ParticleNumberDensity(currParticle);
383  double numb_dens_bolt = particleDensity_pi_th * particleDensity_ch / particleDensity_pi_ch;
384  mu = fThFO * TMath::Log(numb_dens_bolt / particleDensity_th_0);
385  if (abs(encoding) == 211 || encoding == 111)
386  mu = fMu_th_pip;
387  particleDensity = numb_dens_bolt;
388  }
389 
390  // set particle densities to zero for some particle codes
391  // pythia quark codes
392  if (abs(encoding) <= 9) {
393  particleDensity = 0;
394  }
395  // leptons
396  if (abs(encoding) > 10 && abs(encoding) < 19) {
397  particleDensity = 0;
398  }
399  // exchange bosons
400  if (abs(encoding) > 20 && abs(encoding) < 30) {
401  particleDensity = 0;
402  }
403  // pythia special codes (e.g. strings, clusters ...)
404  if (abs(encoding) > 80 && abs(encoding) < 100) {
405  particleDensity = 0;
406  }
407  // pythia di-quark codes
408  // Note: in PYTHIA all diquark codes have the tens digits equal to zero
409  if (abs(encoding) > 1000 && abs(encoding) < 6000) {
410  int tens = ((abs(encoding) - (abs(encoding) % 10)) / 10) % 10;
411  if (tens == 0) { // its a diquark;
412  particleDensity = 0;
413  }
414  }
415  // K0S and K0L
416  if (abs(encoding) == 130 || abs(encoding) == 310) {
417  particleDensity = 0;
418  }
419  // charmed particles
420 
421  if (encoding == 443)
422  NJPsith = particleDensity * fVolEff / dYl;
423 
424  // We generate thermo-statistically only J/psi(443), D_+(411), D_-(-411), D_0(421),
425  //Dbar_0(-421), D1_+(413), D1_-(-413), D1_0(423), D1bar_0(-423)
426  //Dcs(431) Lambdac(4122)
427  if (currParticle->GetCharmQNumber() != 0 || currParticle->GetCharmAQNumber() != 0) {
428  //ml if(abs(encoding)!=443 &&
429  //ml abs(encoding)!=411 && abs(encoding)!=421 &&
430  //ml abs(encoding)!=413 && abs(encoding)!=423 && abs(encoding)!=4122 && abs(encoding)!=431) {
431  //ml particleDensity=0; }
432 
433  if (abs(encoding) == 441 || abs(encoding) == 10441 || abs(encoding) == 10443 || abs(encoding) == 20443 ||
434  abs(encoding) == 445 || abs(encoding) == 4232 || abs(encoding) == 4322 || abs(encoding) == 4132 ||
435  abs(encoding) == 4312 || abs(encoding) == 4324 || abs(encoding) == 4314 || abs(encoding) == 4332 ||
436  abs(encoding) == 4334) {
437  particleDensity = 0;
438  } else {
439  if (abs(encoding) != 443) { //only open charm
440  Nocth = Nocth + particleDensity * fVolEff / dYl;
441  LogInfo("Hydjet2Hadronizer|Charam") << encoding << " Nochth " << Nocth;
442  // particleDensity=particleDensity*fCorrC;
443  // if(abs(encoding)==443)particleDensity=particleDensity*fCorrC;
444  }
445  }
446  }
447 
448  // bottom mesons
449  if ((abs(encoding) > 500 && abs(encoding) < 600) || (abs(encoding) > 10500 && abs(encoding) < 10600) ||
450  (abs(encoding) > 20500 && abs(encoding) < 20600) || (abs(encoding) > 100500 && abs(encoding) < 100600)) {
451  particleDensity = 0;
452  }
453  // bottom baryons
454  if (abs(encoding) > 5000 && abs(encoding) < 6000) {
455  particleDensity = 0;
456  }
458 
459  if (particleDensity > 0.) {
460  fPartEnc[fNPartTypes] = encoding;
461  fPartMult[2 * fNPartTypes] = particleDensity;
462  fPartMu[2 * fNPartTypes] = mu;
463  ++fNPartTypes;
464  if (fNPartTypes > 1000)
465  LogError("Hydjet2Hadronizer") << "fNPartTypes is too large" << fNPartTypes;
466 
467  //outMult<<encoding<<" "<<particleDensity*fVolEff/dYl <<" "<<mu<<std::endl;
468  }
469  }
470 
471  //put open charm number and cc number in Params
472  fNocth = Nocth;
473  fNccth = NJPsith;
474 
475  return kTRUE;
476 }
477 
478 //__________________________________________________________________________________________
481 
482  // Initialize the static "last index variable"
484 
485  //----- high-pt part------------------------------
486  TLorentzVector partJMom, partJPos, zeroVec;
487 
488  // generate single event
489  if (embedding_) {
490  fIfb = 0;
491  const edm::Event& e = getEDMEvent();
493  e.getByLabel(src_, input);
494  const HepMC::GenEvent* inev = input->GetEvent();
495  const HepMC::HeavyIon* hi = inev->heavy_ion();
496  if (hi) {
497  fBfix = hi->impact_parameter();
498  phi0_ = hi->event_plane_angle();
499  sinphi0_ = sin(phi0_);
500  cosphi0_ = cos(phi0_);
501  } else {
502  LogWarning("EventEmbedding") << "Background event does not have heavy ion record!";
503  }
504  } else if (rotate_)
505  rotateEvtPlane();
506 
507  nsoft_ = 0;
508  nhard_ = 0;
509  /*
510  edm::LogInfo("HYDJET2mode") << "##### HYDJET2 fNhsel = " << fNhsel;
511  edm::LogInfo("HYDJET2fpart") << "##### HYDJET2 fpart = " << hyflow.fpart;??
512  edm::LogInfo("HYDJET2tf") << "##### HYDJET2 hadron freez-out temp, Tf = " << hyflow.Tf;??
513  edm::LogInfo("HYDJET2tf") << "##### HYDJET2 hadron freez-out temp, Tf = " << hyflow.Tf;??
514  edm::LogInfo("HYDJET2inTemp") << "##### HYDJET2: QGP init temperature, fT0 ="<<fT0;
515  edm::LogInfo("HYDJET2inTau") << "##### HYDJET2: QGP formation time, fTau0 ="<<fTau0;
516  */
517  // generate a HYDJET event
518  int ntry = 0;
519  while (nsoft_ == 0 && nhard_ == 0) {
520  if (ntry > 100) {
521  LogError("Hydjet2EmptyEvent") << "##### HYDJET2: No Particles generated, Number of tries =" << ntry;
522  // Throw an exception. Use the EventCorruption exception since it maps onto SkipEvent
523  // which is what we want to do here.
524  std::ostringstream sstr;
525  sstr << "Hydjet2HadronizerProducer: No particles generated after " << ntry << " tries.\n";
526  edm::Exception except(edm::errors::EventCorruption, sstr.str());
527  throw except;
528  } else {
529  //generate non-equilibrated part event
530  hyevnt_();
531 
533  if (fNhsel != 0) {
534  //get number of particles in jets
535  int numbJetPart = HYPART.njp;
536 
537  for (int i = 0; i < numbJetPart; ++i) {
538  int pythiaStatus = int(HYPART.ppart[i][0]); // PYTHIA status code
539  int pdg = int(HYPART.ppart[i][1]); // PYTHIA species code
540  double px = HYPART.ppart[i][2]; // px
541  double py = HYPART.ppart[i][3]; // py
542  double pz = HYPART.ppart[i][4]; // pz
543  double e = HYPART.ppart[i][5]; // E
544  double vx = HYPART.ppart[i][6]; // x
545  double vy = HYPART.ppart[i][7]; // y
546  double vz = HYPART.ppart[i][8]; // z
547  double vt = HYPART.ppart[i][9]; // t
548  // particle line number in pythia are 1 based while we use a 0 based numbering
549  int mother_index = int(HYPART.ppart[i][10]) - 1; //line number of parent particle
550  int daughter_index1 = int(HYPART.ppart[i][11]) - 1; //line number of first daughter
551  int daughter_index2 = int(HYPART.ppart[i][12]) - 1; //line number of last daughter
552 
553  // For status codes 3, 13, 14 the first and last daughter indexes have a different meaning
554  // used for color flow in PYTHIA. So these indexes will be reset to zero.
555  if (TMath::Abs(daughter_index1) > numbJetPart || TMath::Abs(daughter_index2) > numbJetPart ||
556  TMath::Abs(daughter_index1) > TMath::Abs(daughter_index2)) {
557  daughter_index1 = -1;
558  daughter_index2 = -1;
559  }
560 
561  ParticlePDG* partDef = fDatabase->GetPDGParticle(pdg);
562 
563  int type = 1; //from jet
564  if (partDef) {
565  int motherPdg = int(HYPART.ppart[mother_index][1]);
566  if (motherPdg == 0)
567  motherPdg = -1;
568  partJMom.SetXYZT(px, py, pz, e);
569  partJPos.SetXYZT(vx, vy, vz, vt);
570  Particle particle(partDef, partJPos, partJMom, 0, 0, type, motherPdg, zeroVec, zeroVec);
571  int index = particle.SetIndex();
572  if (index != i) {
573  // LogWarning("Hydjet2Hadronizer") << " Allocated Hydjet2 index is not synchronized with the PYTHIA index !" << endl
574  // << " Collision history information is destroyed! It happens when a PYTHIA code is not" << endl
575  // << " implemented in Hydjet2 particle list particles.data! Check it out!";
576  }
577  particle.SetPythiaStatusCode(pythiaStatus);
578  particle.SetMother(mother_index);
579  particle.SetFirstDaughterIndex(daughter_index1);
580  particle.SetLastDaughterIndex(daughter_index2);
581  if (pythiaStatus != 1)
582  particle.SetDecayed();
583  allocator.AddParticle(particle, source);
584  } else {
585  LogWarning("Hydjet2Hadronizer")
586  << " PYTHIA particle of specie " << pdg << " is not in Hydjet2 particle list" << endl
587  << " Please define it in particles.data, otherwise the history information will be de-synchronized and "
588  "lost!";
589  }
590  }
591  } //nhsel !=0 not only hydro!
592 
593  //----------HYDRO part--------------------------------------
594 
595  // get impact parameter
596  double impactParameter = HYFPAR.bgen;
597 
598  // Sergey psiforv3
599  psiforv3 = 0.; //AS-ML Nov2012 epsilon3 //
600  double e3 = (0.2 / 5.5) * TMath::Power(impactParameter, 1. / 3.);
601  psiforv3 = TMath::TwoPi() * (-0.5 + CLHEP::RandFlat::shoot(hjRandomEngine)) / 3.;
602  SERVICEEV.psiv3 = -psiforv3;
603 
604  if (fNhsel < 3) {
605  const double weightMax = 2 * TMath::CosH(fUmax);
606  const int nBins = 100;
607  double probList[nBins];
608  RandArrayFunction arrayFunctDistE(nBins);
609  RandArrayFunction arrayFunctDistR(nBins);
610  TLorentzVector partPos, partMom, n1, p0;
611  TVector3 vec3;
612  const TLorentzVector zeroVec;
613  //set maximal hadron energy
614  const double eMax = 5.;
615  //-------------------------------------
616  // get impact parameter
617 
618  double Delta = fDelta;
619  double Epsilon = fEpsilon;
620 
621  if (fIfDeltaEpsilon > 0) {
622  double Epsilon0 = 0.5 * impactParameter; //e0=b/2Ra
623  double coeff = (HYIPAR.RA / fR) / 12.; //phenomenological coefficient
624  Epsilon = Epsilon0 * coeff;
625  double C = 5.6;
626  double A = C * Epsilon * (1 - Epsilon * Epsilon);
627  if (TMath::Abs(Epsilon) < 0.0001 || TMath::Abs(A) < 0.0001)
628  Delta = 0.0;
629  if (TMath::Abs(Epsilon) > 0.0001 && TMath::Abs(A) > 0.0001)
630  Delta = 0.5 * (TMath::Sqrt(1 + 4 * A * (Epsilon + A)) - 1) / A;
631  }
632 
633  //effective volume for central
634  double dYl = 2 * fYlmax; //uniform distr. [-Ylmax; Ylmax]
635  if (fEtaType > 0)
636  dYl = TMath::Sqrt(2 * TMath::Pi()) * fYlmax; //Gaussian distr.
637 
638  double VolEffcent = 2 * TMath::Pi() * fTau * dYl * (fR * fR) / TMath::Power((fUmax), 2) *
639  ((fUmax)*TMath::SinH((fUmax)) - TMath::CosH((fUmax)) + 1);
640 
641  //effective volume for non-central Simpson2
642  double VolEffnoncent = fTau * dYl * SimpsonIntegrator2(0., 2. * TMath::Pi(), Epsilon, Delta);
643 
644  fVolEff = VolEffcent * HYFPAR.npart / HYFPAR.npart0;
645 
646  double coeff_RB = TMath::Sqrt(VolEffcent * HYFPAR.npart / HYFPAR.npart0 / VolEffnoncent);
647  double coeff_R1 = HYFPAR.npart / HYFPAR.npart0;
648  coeff_R1 = TMath::Power(coeff_R1, 0.333333);
649 
650  double Veff = fVolEff;
651  //------------------------------------
652  //cycle on particles types
653 
654  double Nbcol = HYFPAR.nbcol;
655  double NccNN = SERVICE.charm;
656  double Ncc = Nbcol * NccNN / dYl;
657  double Nocth = fNocth;
658  double NJPsith = fNccth;
659 
660  double gammaC = 1.0;
661  if (fCorrC <= 0) {
662  gammaC = CharmEnhancementFactor(Ncc, Nocth, NJPsith, 0.001);
663  } else {
664  gammaC = fCorrC;
665  }
666 
667  LogInfo("Hydjet2Hadronizer|Param") << " gammaC = " << gammaC;
668 
669  for (int i = 0; i < fNPartTypes; ++i) {
670  double Mparam = fPartMult[2 * i] * Veff;
671  const int encoding = fPartEnc[i];
672 
673  //ml if(abs(encoding)==443)Mparam = Mparam * gammaC * gammaC;
674  //ml if(abs(encoding)==411 || abs(encoding)==421 ||abs(encoding)==413 || abs(encoding)==423
675  //ml || abs(encoding)==4122 || abs(encoding)==431)
676 
677  ParticlePDG* partDef0 = fDatabase->GetPDGParticle(encoding);
678 
679  if (partDef0->GetCharmQNumber() != 0 || partDef0->GetCharmAQNumber() != 0)
680  Mparam = Mparam * gammaC;
681  if (abs(encoding) == 443)
682  Mparam = Mparam * gammaC;
683 
684  LogInfo("Hydjet2Hadronizer|Param") << encoding << " " << Mparam / dYl;
685 
686  int multiplicity = CLHEP::RandPoisson::shoot(hjRandomEngine, Mparam);
687 
688  LogInfo("Hydjet2Hadronizer|Param")
689  << "specie: " << encoding << "; average mult: = " << Mparam << "; multiplicity = " << multiplicity;
690 
691  if (multiplicity > 0) {
692  ParticlePDG* partDef = fDatabase->GetPDGParticle(encoding);
693  if (!partDef) {
694  LogError("Hydjet2Hadronizer") << "No particle with encoding " << encoding;
695  continue;
696  }
697 
698  if (fCharmProd <= 0 && (partDef->GetCharmQNumber() != 0 || partDef->GetCharmAQNumber() != 0)) {
699  LogInfo("Hydjet2Hadronizer|Param") << "statistical charmed particle not allowed ! " << encoding;
700  continue;
701  }
702  if (partDef->GetCharmQNumber() != 0 || partDef->GetCharmAQNumber() != 0)
703  LogInfo("Hydjet2Hadronizer|Param") << " charm pdg generated " << encoding;
704 
705  //compute chemical potential for single f.o. mu==mu_ch
706  //compute chemical potential for thermal f.o.
707  double mu = fPartMu[2 * i];
708 
709  //choose Bose-Einstein or Fermi-Dirac statistics
710  const double d = !(int(2 * partDef->GetSpin()) & 1) ? -1 : 1;
711  const double mass = partDef->GetMass();
712 
713  //prepare histogram to sample hadron energy:
714  double h = (eMax - mass) / nBins;
715  double x = mass + 0.5 * h;
716  int i;
717  for (i = 0; i < nBins; ++i) {
718  if (x >= mu && fThFO > 0)
719  probList[i] = x * TMath::Sqrt(x * x - mass * mass) / (TMath::Exp((x - mu) / (fThFO)) + d);
720  if (x >= mu && fThFO <= 0)
721  probList[i] = x * TMath::Sqrt(x * x - mass * mass) / (TMath::Exp((x - mu) / (fT)) + d);
722  if (x < mu)
723  probList[i] = 0.;
724  x += h;
725  }
726  arrayFunctDistE.PrepareTable(probList);
727 
728  //prepare histogram to sample hadron transverse radius:
729  h = (fR) / nBins;
730  x = 0.5 * h;
731  double param = (fUmax) / (fR);
732  for (i = 0; i < nBins; ++i) {
733  probList[i] = x * TMath::CosH(param * x);
734  x += h;
735  }
736  arrayFunctDistR.PrepareTable(probList);
737 
738  //loop over hadrons, assign hadron coordinates and momenta
739  double weight = 0., yy = 0., px0 = 0., py0 = 0., pz0 = 0.;
740  double e = 0., x0 = 0., y0 = 0., z0 = 0., t0 = 0., etaF = 0.;
741  double r, RB, phiF;
742 
743  RB = fR * coeff_RB * coeff_R1 * TMath::Sqrt((1 + e3) / (1 - e3));
744 
745  for (int j = 0; j < multiplicity; ++j) {
746  do {
747  fEtaType <= 0 ? etaF = fYlmax * (2. * CLHEP::RandFlat::shoot(hjRandomEngine) - 1.)
748  : etaF = (fYlmax) * (CLHEP::RandGauss::shoot(hjRandomEngine));
749  n1.SetXYZT(0., 0., TMath::SinH(etaF), TMath::CosH(etaF));
750 
751  if (TMath::Abs(etaF) > 5.)
752  continue;
753 
754  //old
755  //double RBold = fR * TMath::Sqrt(1-fEpsilon);
756 
757  //RB = fR * coeff_RB * coeff_R1;
758 
759  //double impactParameter =HYFPAR.bgen;
760  //double e0 = 0.5*impactParameter;
761  //double RBold1 = fR * TMath::Sqrt(1-e0);
762 
763  double rho = TMath::Sqrt(CLHEP::RandFlat::shoot(hjRandomEngine));
764  double phi = TMath::TwoPi() * CLHEP::RandFlat::shoot(hjRandomEngine);
765  double Rx = TMath::Sqrt(1 - Epsilon) * RB;
766  double Ry = TMath::Sqrt(1 + Epsilon) * RB;
767 
768  x0 = Rx * rho * TMath::Cos(phi);
769  y0 = Ry * rho * TMath::Sin(phi);
770  r = TMath::Sqrt(x0 * x0 + y0 * y0);
771  phiF = TMath::Abs(TMath::ATan(y0 / x0));
772 
773  if (x0 < 0 && y0 > 0)
774  phiF = TMath::Pi() - phiF;
775  if (x0 < 0 && y0 < 0)
776  phiF = TMath::Pi() + phiF;
777  if (x0 > 0 && y0 < 0)
778  phiF = 2. * TMath::Pi() - phiF;
779 
780  //new Nov2012 AS-ML
781  if (r > RB * (1 + e3 * TMath::Cos(3 * (phiF + psiforv3))) / (1 + e3))
782  continue;
783 
784  //proper time with emission duration
785  double tau =
786  coeff_R1 * fTau + sqrt(2.) * fSigmaTau * coeff_R1 * (CLHEP::RandGauss::shoot(hjRandomEngine));
787  z0 = tau * TMath::SinH(etaF);
788  t0 = tau * TMath::CosH(etaF);
789  double rhou = fUmax * r / RB;
790  double rhou3 = 0.063 * TMath::Sqrt((0.5 * impactParameter) / 0.67);
791  double rhou4 = 0.023 * ((0.5 * impactParameter) / 0.67);
792  double rrcoeff = 1. / TMath::Sqrt(1. + Delta * TMath::Cos(2 * phiF));
793  //AS-ML Nov.2012
794  rhou3 = 0.; //rhou4=0.;
795  rhou = rhou * (1 + rrcoeff * rhou3 * TMath::Cos(3 * (phiF + psiforv3)) +
796  rrcoeff * rhou4 * TMath::Cos(4 * phiF)); //ML new suggestion of AS mar2012
797  double delta1 = 0.;
798  Delta = Delta * (1.0 + delta1 * TMath::Cos(phiF) - delta1 * TMath::Cos(3 * phiF));
799 
800  double uxf = TMath::SinH(rhou) * TMath::Sqrt(1 + Delta) * TMath::Cos(phiF);
801  double uyf = TMath::SinH(rhou) * TMath::Sqrt(1 - Delta) * TMath::Sin(phiF);
802  double utf = TMath::CosH(etaF) * TMath::CosH(rhou) *
803  TMath::Sqrt(1 + Delta * TMath::Cos(2 * phiF) * TMath::TanH(rhou) * TMath::TanH(rhou));
804  double uzf = TMath::SinH(etaF) * TMath::CosH(rhou) *
805  TMath::Sqrt(1 + Delta * TMath::Cos(2 * phiF) * TMath::TanH(rhou) * TMath::TanH(rhou));
806 
807  vec3.SetXYZ(uxf / utf, uyf / utf, uzf / utf);
808  n1.Boost(-vec3);
809 
810  yy = weightMax * CLHEP::RandFlat::shoot(hjRandomEngine);
811 
812  double php0 = TMath::TwoPi() * CLHEP::RandFlat::shoot(hjRandomEngine);
813  double ctp0 = 2. * CLHEP::RandFlat::shoot(hjRandomEngine) - 1.;
814  double stp0 = TMath::Sqrt(1. - ctp0 * ctp0);
815  e = mass + (eMax - mass) * arrayFunctDistE();
816  double pp0 = TMath::Sqrt(e * e - mass * mass);
817  px0 = pp0 * stp0 * TMath::Sin(php0);
818  py0 = pp0 * stp0 * TMath::Cos(php0);
819  pz0 = pp0 * ctp0;
820  p0.SetXYZT(px0, py0, pz0, e);
821 
822  //weight for rdr
823  weight = (n1 * p0) / e; // weight for rdr gammar: weight = (n1 * p0) / n1[3] / e;
824 
825  } while (yy >= weight);
826 
827  if (abs(z0) > 1000 || abs(x0) > 1000)
828  LogInfo("Hydjet2Hadronizer|Param") << " etaF = " << etaF << std::endl;
829 
830  partMom.SetXYZT(px0, py0, pz0, e);
831  partPos.SetXYZT(x0, y0, z0, t0);
832  partMom.Boost(vec3);
833 
834  int type = 0; //hydro
835  Particle particle(partDef, partPos, partMom, 0., 0, type, -1, zeroVec, zeroVec);
836  particle.SetIndex();
837  allocator.AddParticle(particle, source);
838  } //nhsel==4 , no hydro part
839  }
840  }
841  }
842 
844 
845  Npart = (int)HYFPAR.npart;
846  Bgen = HYFPAR.bgen;
847  Njet = (int)HYJPAR.njet;
848  Nbcol = (int)HYFPAR.nbcol;
849 
850  //if(source.empty()) {
851  // LogError("Hydjet2Hadronizer") << "Source is not initialized!! Trying again...";
852  //return ;
853  //}
854  //Run the decays
855  if (RunDecays())
857 
858  LPIT_t it;
859  LPIT_t e;
860 
861  //Fill the decayed arrays
862  Ntot = 0;
863  Nhyd = 0;
864  Npyt = 0;
865  for (it = source.begin(), e = source.end(); it != e; ++it) {
866  TVector3 pos(it->Pos().Vect());
867  TVector3 mom(it->Mom().Vect());
868  float m1 = it->TableMass();
869  pdg[Ntot] = it->Encoding();
870  Mpdg[Ntot] = it->GetLastMotherPdg();
871  Px[Ntot] = mom[0];
872  Py[Ntot] = mom[1];
873  Pz[Ntot] = mom[2];
874  E[Ntot] = TMath::Sqrt(mom.Mag2() + m1 * m1);
875  X[Ntot] = pos[0];
876  Y[Ntot] = pos[1];
877  Z[Ntot] = pos[2];
878  T[Ntot] = it->T();
879  type[Ntot] = it->GetType();
880  pythiaStatus[Ntot] = convertStatus(it->GetPythiaStatusCode());
881  Index[Ntot] = it->GetIndex();
882  MotherIndex[Ntot] = it->GetMother();
883  NDaughters[Ntot] = it->GetNDaughters();
884  FirstDaughterIndex[Ntot] = -1;
885  LastDaughterIndex[Ntot] = -1;
886  //index of first daughter
887  FirstDaughterIndex[Ntot] = it->GetFirstDaughterIndex();
888  //index of last daughter
889  LastDaughterIndex[Ntot] = it->GetLastDaughterIndex();
890  if (type[Ntot] == 1) { // jets
891  if (pythiaStatus[Ntot] == 1 && NDaughters[Ntot] == 0) { // code for final state particle in pythia
892  final[Ntot] = 1;
893  } else {
894  final[Ntot] = pythiaStatus[Ntot];
895  }
896  }
897  if (type[Ntot] == 0) { // hydro
898  if (NDaughters[Ntot] == 0)
899  final[Ntot] = 1;
900  else
901  final[Ntot] = 2;
902  }
903 
904  if (type[Ntot] == 0)
905  Nhyd++;
906  if (type[Ntot] == 1)
907  Npyt++;
908 
909  Ntot++;
910  if (Ntot > kMax)
911  LogError("Hydjet2Hadronizer") << "Ntot is too large" << Ntot;
912  }
913 
914  nsoft_ = Nhyd;
915  nsub_ = Njet;
916  nhard_ = Npyt;
917 
918  //100 trys
919 
920  ++ntry;
921  }
922  }
923 
924  if (ev == 0) {
925  Sigin = HYJPAR.sigin;
926  Sigjet = HYJPAR.sigjet;
927  }
928  ev = 1;
929 
930  if (fNhsel < 3)
931  nsub_++;
932 
933  // event information
935  if (nhard_ > 0 || nsoft_ > 0)
936  get_particles(evt);
937 
938  evt->set_signal_process_id(pypars.msti[0]); // type of the process
939  evt->set_event_scale(pypars.pari[16]); // Q^2
940  add_heavy_ion_rec(evt);
941 
942  event().reset(evt);
943 
945 
946  return kTRUE;
947 }
948 
949 //________________________________________________________________
950 bool Hydjet2Hadronizer::declareStableParticles(const std::vector<int>& _pdg) {
951  std::vector<int> pdg = _pdg;
952  for (size_t i = 0; i < pdg.size(); i++) {
953  int pyCode = pycomp_(pdg[i]);
954  std::ostringstream pyCard;
955  pyCard << "MDCY(" << pyCode << ",1)=0";
956  std::cout << pyCard.str() << std::endl;
957  call_pygive(pyCard.str());
958  }
959  return true;
960 }
961 //________________________________________________________________
962 bool Hydjet2Hadronizer::hadronize() { return false; }
963 bool Hydjet2Hadronizer::decay() { return true; }
964 bool Hydjet2Hadronizer::residualDecay() { return true; }
967 const char* Hydjet2Hadronizer::classname() const { return "gen::Hydjet2Hadronizer"; }
968 
969 //----------------------------------------------------------------------------------------------
970 
971 //______________________________________________________________________________________________
972 //f2=f(phi,r)
973 double Hydjet2Hadronizer::f2(double x, double y, double Delta) {
974  LogDebug("f2") << "in f2: "
975  << "delta" << Delta;
976  double RsB = fR; //test: podstavit' *coefff_RB
977  double rhou = fUmax * y / RsB;
978  double ff =
979  y * TMath::CosH(rhou) * TMath::Sqrt(1 + Delta * TMath::Cos(2 * x) * TMath::TanH(rhou) * TMath::TanH(rhou));
980  //n_mu u^mu f-la 20
981  return ff;
982 }
983 
984 //____________________________________________________________________________________________
985 double Hydjet2Hadronizer::SimpsonIntegrator(double a, double b, double phi, double Delta) {
986  LogDebug("SimpsonIntegrator") << "in SimpsonIntegrator"
987  << "delta - " << Delta;
988  int nsubIntervals = 100;
989  double h = (b - a) / nsubIntervals;
990  double s = f2(phi, a + 0.5 * h, Delta);
991  double t = 0.5 * (f2(phi, a, Delta) + f2(phi, b, Delta));
992  double x = a;
993  double y = a + 0.5 * h;
994  for (int i = 1; i < nsubIntervals; i++) {
995  x += h;
996  y += h;
997  s += f2(phi, y, Delta);
998  t += f2(phi, x, Delta);
999  }
1000  t += 2.0 * s;
1001  return t * h / 3.0;
1002 }
1003 
1004 //______________________________________________________________________________________________
1005 double Hydjet2Hadronizer::SimpsonIntegrator2(double a, double b, double Epsilon, double Delta) {
1006  LogInfo("SimpsonIntegrator2") << "in SimpsonIntegrator2: epsilon - " << Epsilon << " delta - " << Delta;
1007  int nsubIntervals = 10000;
1008  double h = (b - a) / nsubIntervals; //-1-pi, phi
1009  double s = 0;
1010  // double h2 = (fR)/nsubIntervals; //0-R maximal RB ?
1011 
1012  double x = 0; //phi
1013  for (int j = 1; j < nsubIntervals; j++) {
1014  x += h; // phi
1015  double e = Epsilon;
1016  double RsB = fR; //test: podstavit' *coefff_RB
1017  double RB = RsB * (TMath::Sqrt(1 - e * e) / TMath::Sqrt(1 + e * TMath::Cos(2 * x))); //f-la7 RB
1018  double sr = SimpsonIntegrator(0, RB, x, Delta);
1019  s += sr;
1020  }
1021  return s * h;
1022 }
1023 
1024 //___________________________________________________________________________________________________
1025 double Hydjet2Hadronizer::MidpointIntegrator2(double a, double b, double Delta, double Epsilon) {
1026  int nsubIntervals = 2000;
1027  int nsubIntervals2 = 1;
1028  double h = (b - a) / nsubIntervals; //0-pi , phi
1029  double h2 = (fR) / nsubIntervals; //0-R maximal RB ?
1030 
1031  double x = a + 0.5 * h;
1032  double y = 0;
1033 
1034  double t = f2(x, y, Delta);
1035 
1036  double e = Epsilon;
1037 
1038  for (int j = 1; j < nsubIntervals; j++) {
1039  x += h; // integr phi
1040 
1041  double RsB = fR; //test: podstavit' *coefff_RB
1042  double RB = RsB * (TMath::Sqrt(1 - e * e) / TMath::Sqrt(1 + e * TMath::Cos(2 * x))); //f-la7 RB
1043 
1044  nsubIntervals2 = int(RB / h2) + 1;
1045  // integr R
1046  y = 0;
1047  for (int i = 1; i < nsubIntervals2; i++)
1048  t += f2(x, (y += h2), Delta);
1049  }
1050  return t * h * h2;
1051 }
1052 
1053 //__________________________________________________________________________________________________________
1054 double Hydjet2Hadronizer::CharmEnhancementFactor(double Ncc, double Ndth, double NJPsith, double Epsilon) {
1055  double gammaC = 100.;
1056  double x1 = gammaC * Ndth;
1057  double var1 = Ncc - 0.5 * gammaC * Ndth * TMath::BesselI1(x1) / TMath::BesselI0(x1) - gammaC * gammaC * NJPsith;
1058  LogInfo("Charam") << "gammaC 20"
1059  << " var " << var1 << endl;
1060  gammaC = 1.;
1061  double x0 = gammaC * Ndth;
1062  double var0 = Ncc - 0.5 * gammaC * Ndth * TMath::BesselI1(x0) / TMath::BesselI0(x0) - gammaC * gammaC * NJPsith;
1063  LogInfo("Charam") << "gammaC 1"
1064  << " var " << var0;
1065 
1066  for (int i = 1; i < 1000; i++) {
1067  if (var1 * var0 < 0) {
1068  gammaC = gammaC + 0.01 * i;
1069  double x = gammaC * Ndth;
1070  var0 = Ncc - 0.5 * gammaC * Ndth * TMath::BesselI1(x) / TMath::BesselI0(x) - gammaC * gammaC * NJPsith;
1071  } else {
1072  LogInfo("Charam") << "gammaC " << gammaC << " var0 " << var0;
1073  return gammaC;
1074  }
1075  }
1076  LogInfo("Charam") << "gammaC not found ? " << gammaC << " var0 " << var0;
1077  return -100;
1078 }
1079 //----------------------------------------------------------------------------------------------
1080 
1081 //_____________________________________________________________________
1082 
1083 //________________________________________________________________
1085  const double pi = 3.14159265358979;
1086  phi0_ = 2. * pi * gen::pyr_(nullptr) - pi;
1087  sinphi0_ = sin(phi0_);
1088  cosphi0_ = cos(phi0_);
1089 }
1090 
1091 //_____________________________________________________________________
1093  // Hard particles. The first nhard_ lines from hyjets array.
1094  // Pythia/Pyquen sub-events (sub-collisions) for a given event
1095  // Return T/F if success/failure
1096  // Create particles from lujet entries, assign them into vertices and
1097  // put the vertices in the GenEvent, for each SubEvent
1098  // The SubEvent information is kept by storing indeces of main vertices
1099  // of subevents as a vector in GenHIEvent.
1100  LogDebug("SubEvent") << "Number of sub events " << nsub_;
1101  LogDebug("Hydjet2") << "Number of hard events " << Njet;
1102  LogDebug("Hydjet2") << "Number of hard particles " << nhard_;
1103  LogDebug("Hydjet2") << "Number of soft particles " << nsoft_;
1104 
1105  vector<HepMC::GenVertex*> sub_vertices(nsub_);
1106 
1107  int ihy = 0;
1108  for (int isub = 0; isub < nsub_; isub++) {
1109  LogDebug("SubEvent") << "Sub Event ID : " << isub;
1110 
1111  int sub_up = (isub + 1) * 150000; // Upper limit in mother index, determining the range of Sub-Event
1112  vector<HepMC::GenParticle*> particles;
1113  vector<int> mother_ids;
1114  vector<HepMC::GenVertex*> prods;
1115 
1116  sub_vertices[isub] = new HepMC::GenVertex(HepMC::FourVector(0, 0, 0, 0), isub);
1117  evt->add_vertex(sub_vertices[isub]);
1118 
1119  if (!evt->signal_process_vertex())
1120  evt->set_signal_process_vertex(sub_vertices[isub]);
1121 
1122  while (ihy < nhard_ + nsoft_ && (MotherIndex[ihy] < sub_up || ihy > nhard_)) {
1123  particles.push_back(build_hyjet2(ihy, ihy + 1));
1124  prods.push_back(build_hyjet2_vertex(ihy, isub));
1125  mother_ids.push_back(MotherIndex[ihy]);
1126  LogDebug("DecayChain") << "Mother index : " << MotherIndex[ihy];
1127  ihy++;
1128  }
1129  //Produce Vertices and add them to the GenEvent. Remember that GenParticles are adopted by
1130  //GenVertex and GenVertex is adopted by GenEvent.
1131  LogDebug("Hydjet2") << "Number of particles in vector " << particles.size();
1132 
1133  for (unsigned int i = 0; i < particles.size(); i++) {
1134  HepMC::GenParticle* part = particles[i];
1135  //The Fortran code is modified to preserve mother id info, by seperating the beginning
1136  //mother indices of successive subevents by 5000
1137  int mid = mother_ids[i] - isub * 150000 - 1;
1138  LogDebug("DecayChain") << "Particle " << i;
1139  LogDebug("DecayChain") << "Mother's ID " << mid;
1140  LogDebug("DecayChain") << "Particle's PDG ID " << part->pdg_id();
1141 
1142  if (mid <= 0) {
1143  sub_vertices[isub]->add_particle_out(part);
1144  continue;
1145  }
1146 
1147  if (mid > 0) {
1148  HepMC::GenParticle* mother = particles[mid];
1149  LogDebug("DecayChain") << "Mother's PDG ID " << mother->pdg_id();
1150  HepMC::GenVertex* prod_vertex = mother->end_vertex();
1151  if (!prod_vertex) {
1152  prod_vertex = prods[i];
1153  prod_vertex->add_particle_in(mother);
1154  evt->add_vertex(prod_vertex);
1155  prods[i] = nullptr; // mark to protect deletion
1156  }
1157 
1158  prod_vertex->add_particle_out(part);
1159  }
1160  }
1161 
1162  // cleanup vertices not assigned to evt
1163  for (unsigned int i = 0; i < prods.size(); i++) {
1164  if (prods[i])
1165  delete prods[i];
1166  }
1167  }
1168 
1169  return kTRUE;
1170 }
1171 
1172 //___________________________________________________________________
1174  // Build particle object corresponding to index in hyjets (soft+hard)
1175 
1176  double px0 = Px[index];
1177  double py0 = Py[index];
1178 
1179  double px = px0 * cosphi0_ - py0 * sinphi0_;
1180  double py = py0 * cosphi0_ + px0 * sinphi0_;
1181  // cout<< "status: "<<convertStatus(final[index], type[index])<<endl;
1182  HepMC::GenParticle* p = new HepMC::GenParticle(HepMC::FourVector(px, // px
1183  py, // py
1184  Pz[index], // pz
1185  E[index]), // E
1186  pdg[index], // id
1187  convertStatusForComponents(final[index], type[index]) // status
1188 
1189  );
1190 
1191  p->suggest_barcode(barcode);
1192  return p;
1193 }
1194 
1195 //___________________________________________________________________
1196 HepMC::GenVertex* Hydjet2Hadronizer::build_hyjet2_vertex(int i, int id) {
1197  // build verteces for the hyjets stored events
1198 
1199  double x0 = X[i];
1200  double y0 = Y[i];
1201  double x = x0 * cosphi0_ - y0 * sinphi0_;
1202  double y = y0 * cosphi0_ + x0 * sinphi0_;
1203  double z = Z[i];
1204  double t = T[i];
1205 
1206  HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x, y, z, t), id);
1207  return vertex;
1208 }
1209 
1210 //_____________________________________________________________________
1212  // heavy ion record in the final CMSSW Event
1213  double npart = Npart;
1214  int nproj = static_cast<int>(npart / 2);
1215  int ntarg = static_cast<int>(npart - nproj);
1216 
1217  HepMC::HeavyIon* hi = new HepMC::HeavyIon(nsub_, // Ncoll_hard/N of SubEvents
1218  nproj, // Npart_proj
1219  ntarg, // Npart_targ
1220  Nbcol, // Ncoll
1221  0, // spectator_neutrons
1222  0, // spectator_protons
1223  0, // N_Nwounded_collisions
1224  0, // Nwounded_N_collisions
1225  0, // Nwounded_Nwounded_collisions
1226  Bgen * nuclear_radius(), // impact_parameter in [fm]
1227  phi0_, // event_plane_angle
1228  0, // eccentricity
1229  Sigin // sigma_inel_NN
1230  );
1231 
1232  evt->set_heavy_ion(*hi);
1233  delete hi;
1234 }
#define LogDebug(id)
const double TwoPi
const double Pi
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
void SetBaryonPotential(double value)
bool declareStableParticles(const std::vector< int > &)
ParticlePDG * GetPDGParticle(int pdg)
Definition: DatabasePDG.cc:218
#define HYPYIN
Definition: HYJET_COMMONS.h:84
int GetNParticles(bool all=kFALSE)
Definition: DatabasePDG.cc:582
double ParticleNumberDensity(ParticlePDG *particle)
ParticleAllocator allocator
#define HYJPAR
Definition: HYJET_COMMONS.h:69
#define PYQPAR
TString RunInputHYDJETstr
int GetPDG()
Definition: ParticlePDG.h:69
void add_heavy_ion_rec(HepMC::GenEvent *evt)
#define HYPART
unsigned int pythiaPylistVerbosity_
double SimpsonIntegrator2(double, double, double, double)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double BesselI0(double x)
#define nullptr
bool call_pygive(const std::string &line)
double CalculateStrangePotential()
fCorrS
Strangeness suppression factor ###.
fIfDeltaEpsilon
Anizotropy parameter at thermal freeze-out ###.
Definition: weight.py:1
double npart
Definition: HydjetWrapper.h:46
bool ev
double MidpointIntegrator2(double, double, double, double)
double SimpsonIntegrator(double, double, double, double)
#define SERVICEEV
Definition: HYJET_COMMONS.h:52
double f2(double, double, double)
double v[5][pyjets_maxn]
double GetCharmAQNumber()
Definition: ParticlePDG.h:81
static std::string const input
Definition: EdmProvDump.cc:48
double nuclear_radius() const
#define HYFPAR
Definition: HYJET_COMMONS.h:99
const Double_t pi
static const std::string kPythia6
HepMC::GenVertex * build_hyjet2_vertex(int i, int id)
double GetSpin()
Definition: ParticlePDG.h:73
void PrepareTable(const double *aProbFunc)
fTMuType
Thermodinamic parameters at chemical freez-out ###.
Pythia6Service * pythia6Service_
double p[5][pyjets_maxn]
ParticlePDG * GetPDGParticleByIndex(int index)
Definition: DatabasePDG.cc:198
std::vector< vec2 > vec3
Definition: HCALResponse.h:17
double GetCharmness()
Definition: ParticlePDG.h:88
void FreeList(List_t &list)
Definition: Particle.cc:118
bool get_particles(HepMC::GenEvent *evt)
void hyevnt_()
double GetElectricCharge()
Definition: ParticlePDG.h:89
T sqrt(T t)
Definition: SSEVec.h:19
fTau
Volume parameters at thermal freeze-out ###.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T Abs(T a)
Definition: MathUtil.h:49
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double CharmEnhancementFactor(double, double, double, double)
double GetBaryonNumber()
Definition: ParticlePDG.h:82
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
std::unique_ptr< HepMC::GenEvent > & event()
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
d
Definition: ztail.py:151
CLHEP::HepRandomEngine * hjRandomEngine
double GetStrangeness()
Definition: ParticlePDG.h:87
#define HYIPAR
Definition: HYJET_COMMONS.h:26
virtual void Evolve(List_t &secondaries, ParticleAllocator &allocator, double weakDecayLimit)
Definition: InitialState.cc:13
static int fLastIndex
Definition: Particle.h:38
static void InitIndexing()
Definition: Particle.h:140
#define SERVICE
Definition: HYJET_COMMONS.h:38
const char * classname() const
static const std::string kFortranInstance
int pycomp_(int &)
#define pypars
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
fSigmaTau
Volume parameters at thermal freeze-out ###.
part
Definition: HCALResponse.h:20
double b
Definition: hdecay.h:118
void myini_()
HepMC::GenParticle * build_hyjet2(int index, int barcode)
void setRandomEngine(CLHEP::HepRandomEngine *v)
bool RunDecays() override
HLT enums.
double GetCharmQNumber()
Definition: ParticlePDG.h:80
DatabasePDG * fDatabase
Definition: InitialState.h:16
double GetWeakDecayLimit() override
void SetTemperature(double value)
double GetMass()
Definition: ParticlePDG.h:70
double BesselI1(double x)
double a
Definition: hdecay.h:119
double ppart[150000][20]
#define kMax
void AddParticle(const Particle &particle, List_t &list)
Definition: Particle.cc:107
std::list< Particle >::iterator LPIT_t
Definition: Particle.h:175
double pyr_(int *idummy)
long double T
fThFO
Thermodinamic parameters at thermal freez-out ###.
edm::Event & getEDMEvent() const
fYlmax
Maximal longitudinal flow rapidity at thermal freeze-out ###.