CMS 3D CMS Logo

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