CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HydjetHadronizer.cc
Go to the documentation of this file.
1 /*
2 
3  ########################
4  # Hydjet1 #
5  # version: 1.9 patch1 #
6  ########################
7 
8  * Interface to the HYDJET generator, produces HepMC events
9  *
10  * Original Author: Camelia Mironov
11  */
12 
13 #include <iostream>
14 #include <cmath>
15 
16 #include "boost/lexical_cast.hpp"
17 
24 
30 
31 #include "HepMC/PythiaWrapper6_4.h"
32 #include "HepMC/GenEvent.h"
33 #include "HepMC/HeavyIon.h"
34 #include "HepMC/SimpleVector.h"
35 
40 
41 using namespace edm;
42 using namespace std;
43 using namespace gen;
44 
45 namespace {
46  int convertStatus(int st){
47  if(st<= 0) return 0;
48  if(st<=10) return 1;
49  if(st<=20) return 2;
50  if(st<=30) return 3;
51  else return st;
52  }
53 }
54 
55 const std::vector<std::string> HydjetHadronizer::theSharedResources = { edm::SharedResourceNames::kPythia6,
57 
58 //_____________________________________________________________________
59 HydjetHadronizer::HydjetHadronizer(const ParameterSet &pset) :
60  BaseHadronizer(pset),
61  evt(0),
62  pset_(pset),
63  abeamtarget_(pset.getParameter<double>("aBeamTarget")),
64  angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
65  bfixed_(pset.getParameter<double>("bFixed")),
66  bmax_(pset.getParameter<double>("bMax")),
67  bmin_(pset.getParameter<double>("bMin")),
68  cflag_(pset.getParameter<int>("cFlag")),
69  embedding_(pset.getParameter<bool>("embeddingMode")),
70  comenergy(pset.getParameter<double>("comEnergy")),
71  doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
72  docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
73  fracsoftmult_(pset.getParameter<double>("fracSoftMultiplicity")),
74  hadfreeztemp_(pset.getParameter<double>("hadronFreezoutTemperature")),
75  hymode_(pset.getParameter<string>("hydjetMode")),
76  maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)),
77  maxlongy_(pset.getParameter<double>("maxLongitudinalRapidity")),
78  maxtrany_(pset.getParameter<double>("maxTransverseRapidity")),
79  nsub_(0),
80  nhard_(0),
81  nmultiplicity_(pset.getParameter<int>("nMultiplicity")),
82  nsoft_(0),
83  nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
84  pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
85  qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
86  qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
87  phi0_(0.),
88  sinphi0_(0.),
89  cosphi0_(1.),
90  rotate_(pset.getParameter<bool>("rotateEventPlane")),
91  shadowingswitch_(pset.getParameter<int>("shadowingSwitch")),
92  signn_(pset.getParameter<double>("sigmaInelNN")),
93  pythia6Service_(new Pythia6Service(pset))
94 {
95  // Default constructor
96 
97  // PYLIST Verbosity Level
98  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
99  pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0);
100  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_;
101 
102  //Max number of events printed on verbosity level
103  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
104  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_;
105 
106  if(embedding_) src_ = pset.getParameter<edm::InputTag>("backgroundLabel");
107 
108 }
109 
110 
111 //_____________________________________________________________________
113 {
114  // destructor
115  call_pystat(1);
116  delete pythia6Service_;
117 }
118 
119 
120 //_____________________________________________________________________
121 void HydjetHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v)
122 {
124 }
125 
126 
127 //_____________________________________________________________________
128 void HydjetHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
129 {
130  // heavy ion record in the final CMSSW Event
131  double npart = hyfpar.npart;
132  int nproj = static_cast<int>(npart / 2);
133  int ntarg = static_cast<int>(npart - nproj);
134 
135  HepMC::HeavyIon* hi = new HepMC::HeavyIon(
136  nsub_, // Ncoll_hard/N of SubEvents
137  nproj, // Npart_proj
138  ntarg, // Npart_targ
139  static_cast<int>(hyfpar.nbcol), // Ncoll
140  0, // spectator_neutrons
141  0, // spectator_protons
142  0, // N_Nwounded_collisions
143  0, // Nwounded_N_collisions
144  0, // Nwounded_Nwounded_collisions
145  hyfpar.bgen * nuclear_radius(), // impact_parameter in [fm]
146  phi0_, // event_plane_angle
147  0,//hypsi3.psi3, // eccentricity
148  hyjpar.sigin // sigma_inel_NN
149  );
150 
151  evt->set_heavy_ion(*hi);
152  delete hi;
153 }
154 
155 //___________________________________________________________________
157 {
158  // Build particle object corresponding to index in hyjets (soft+hard)
159 
160  double x0 = hyjets.phj[0][index];
161  double y0 = hyjets.phj[1][index];
162 
163  double x = x0*cosphi0_-y0*sinphi0_;
164  double y = y0*cosphi0_+x0*sinphi0_;
165 
167  HepMC::FourVector(x, // px
168  y, // py
169  hyjets.phj[2][index], // pz
170  hyjets.phj[3][index]), // E
171  hyjets.khj[1][index], // id
172  convertStatus(hyjets.khj[0][index] // status
173  )
174  );
175 
176  p->suggest_barcode(barcode);
177  return p;
178 }
179 
180 //___________________________________________________________________
181 HepMC::GenVertex* HydjetHadronizer::build_hyjet_vertex(int i,int id)
182 {
183  // build verteces for the hyjets stored events
184 
185  double x0=hyjets.vhj[0][i];
186  double y0=hyjets.vhj[1][i];
187  double x = x0*cosphi0_-y0*sinphi0_;
188  double y = y0*cosphi0_+x0*sinphi0_;
189  double z=hyjets.vhj[2][i];
190  double t=hyjets.vhj[4][i];
191 
192  HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t),id);
193  return vertex;
194 }
195 
196 //___________________________________________________________________
197 
199 {
201 
202  // generate single event
203  if(embedding_){
204  cflag_ = 0;
205  const edm::Event& e = getEDMEvent();
207  e.getByLabel(src_,input);
208  const HepMC::GenEvent * inev = input->GetEvent();
209  const HepMC::HeavyIon* hi = inev->heavy_ion();
210  if(hi){
211  bfixed_ = hi->impact_parameter();
212  phi0_ = hi->event_plane_angle();
213  sinphi0_ = sin(phi0_);
214  cosphi0_ = cos(phi0_);
215  }else{
216  LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
217  }
218  }else if(rotate_) rotateEvtPlane();
219 
220  nsoft_ = 0;
221  nhard_ = 0;
222 
223  edm::LogInfo("HYDJETmode") << "##### HYDJET nhsel = " << hyjpar.nhsel;
224  edm::LogInfo("HYDJETfpart") << "##### HYDJET fpart = " << hyflow.fpart;
225  edm::LogInfo("HYDJETtf") << "##### HYDJET hadron freez-out temp, Tf = " << hyflow.Tf;
226  edm::LogInfo("HYDJETinTemp") << "##### HYDJET: QGP init temperature, T0 ="<<pyqpar.T0u;
227  edm::LogInfo("HYDJETinTau") << "##### HYDJET: QGP formation time,tau0 ="<<pyqpar.tau0u;
228 
229  // generate a HYDJET event
230  int ntry = 0;
231  while(nsoft_ == 0 && nhard_ == 0){
232  if(ntry > 100){
233  edm::LogError("HydjetEmptyEvent") << "##### HYDJET: No Particles generated, Number of tries ="<<ntry;
234 
235  // Throw an exception. Use the EventCorruption exception since it maps onto SkipEvent
236  // which is what we want to do here.
237 
238  std::ostringstream sstr;
239  sstr << "HydjetHadronizerProducer: No particles generated after " << ntry << " tries.\n";
240  edm::Exception except(edm::errors::EventCorruption, sstr.str());
241  throw except;
242  } else {
243  HYEVNT();
244  nsoft_ = hyfpar.nhyd;
245  nsub_ = hyjpar.njet;
246  nhard_ = hyfpar.npyt;
247  ++ntry;
248  }
249  }
250 
251  if(hyjpar.nhsel < 3) nsub_++;
252 
253  // event information
254  HepMC::GenEvent *evt = new HepMC::GenEvent();
255 
256  if(nhard_>0 || nsoft_>0) get_particles(evt);
257 
258  evt->set_signal_process_id(pypars.msti[0]); // type of the process
259  evt->set_event_scale(pypars.pari[16]); // Q^2
260  add_heavy_ion_rec(evt);
261 
262  event().reset(evt);
263  return true;
264 }
265 
266 
267 //_____________________________________________________________________
268 bool HydjetHadronizer::get_particles(HepMC::GenEvent *evt )
269 {
270  // Hard particles. The first nhard_ lines from hyjets array.
271  // Pythia/Pyquen sub-events (sub-collisions) for a given event
272  // Return T/F if success/failure
273  // Create particles from lujet entries, assign them into vertices and
274  // put the vertices in the GenEvent, for each SubEvent
275  // The SubEvent information is kept by storing indeces of main vertices
276  // of subevents as a vector in GenHIEvent.
277 
278  LogDebug("SubEvent")<< "Number of sub events "<<nsub_;
279  LogDebug("Hydjet")<<"Number of hard events "<<hyjpar.njet;
280  LogDebug("Hydjet")<<"Number of hard particles "<<nhard_;
281  LogDebug("Hydjet")<<"Number of soft particles "<<nsoft_;
282 
283  vector<HepMC::GenVertex*> sub_vertices(nsub_);
284 
285  int ihy = 0;
286  for(int isub=0;isub<nsub_;isub++){
287  LogDebug("SubEvent") <<"Sub Event ID : "<<isub;
288 
289  int sub_up = (isub+1)*50000; // Upper limit in mother index, determining the range of Sub-Event
290  vector<HepMC::GenParticle*> particles;
291  vector<int> mother_ids;
292  vector<HepMC::GenVertex*> prods;
293 
294  sub_vertices[isub] = new HepMC::GenVertex(HepMC::FourVector(0,0,0,0),isub);
295  evt->add_vertex(sub_vertices[isub]);
296  if(!evt->signal_process_vertex()) evt->set_signal_process_vertex(sub_vertices[isub]);
297 
298  while(ihy<nhard_+nsoft_ && (hyjets.khj[2][ihy] < sub_up || ihy > nhard_ )){
299  particles.push_back(build_hyjet(ihy,ihy+1));
300  prods.push_back(build_hyjet_vertex(ihy,isub));
301  mother_ids.push_back(hyjets.khj[2][ihy]);
302  LogDebug("DecayChain")<<"Mother index : "<<hyjets.khj[2][ihy];
303 
304  ihy++;
305  }
306 
307  //Produce Vertices and add them to the GenEvent. Remember that GenParticles are adopted by
308  //GenVertex and GenVertex is adopted by GenEvent.
309 
310  LogDebug("Hydjet")<<"Number of particles in vector "<<particles.size();
311 
312  for (unsigned int i = 0; i<particles.size(); i++) {
313  HepMC::GenParticle* part = particles[i];
314 
315  //The Fortran code is modified to preserve mother id info, by seperating the beginning
316  //mother indices of successive subevents by 5000
317  int mid = mother_ids[i]-isub*50000-1;
318  LogDebug("DecayChain")<<"Particle "<<i;
319  LogDebug("DecayChain")<<"Mother's ID "<<mid;
320  LogDebug("DecayChain")<<"Particle's PDG ID "<<part->pdg_id();
321 
322  if(mid <= 0){
323  sub_vertices[isub]->add_particle_out(part);
324  continue;
325  }
326 
327  if(mid > 0){
328  HepMC::GenParticle* mother = particles[mid];
329  LogDebug("DecayChain")<<"Mother's PDG ID "<<mother->pdg_id();
330 
331  HepMC::GenVertex* prod_vertex = mother->end_vertex();
332  if(!prod_vertex){
333  prod_vertex = prods[i];
334  prod_vertex->add_particle_in(mother);
335  evt->add_vertex(prod_vertex);
336  prods[i]=0; // mark to protect deletion
337  }
338  prod_vertex->add_particle_out(part);
339  }
340  }
341  // cleanup vertices not assigned to evt
342  for (unsigned int i = 0; i<prods.size(); i++) {
343  if(prods[i]) delete prods[i];
344  }
345  }
346  return true;
347 }
348 
349 
350 //______________________________________________________________
351 bool HydjetHadronizer::call_hyinit(double energy,double a, int ifb, double bmin,
352  double bmax,double bfix,int nh)
353 {
354  // initialize hydjet
355 
356  pydatr.mrpy[2]=1;
357  HYINIT(energy,a,ifb,bmin,bmax,bfix,nh);
358  return true;
359 }
360 
361 
362 //______________________________________________________________
364 {
365  // set hydjet options
366 
367  // hydjet running mode mode
368  // kHydroOnly --- nhsel=0 jet production off (pure HYDRO event), nhsel=0
369  // kHydroJets --- nhsle=1 jet production on, jet quenching off (HYDRO+njet*PYTHIA events)
370  // kHydroQJet --- nhsel=2 jet production & jet quenching on (HYDRO+njet*PYQUEN events)
371  // kJetsOnly --- nhsel=3 jet production on, jet quenching off, HYDRO off (njet*PYTHIA events)
372  // kQJetsOnly --- nhsel=4 jet production & jet quenching on, HYDRO off (njet*PYQUEN events)
373 
374  if(hymode_ == "kHydroOnly") hyjpar.nhsel=0;
375  else if ( hymode_ == "kHydroJets") hyjpar.nhsel=1;
376  else if ( hymode_ == "kHydroQJets") hyjpar.nhsel=2;
377  else if ( hymode_ == "kJetsOnly") hyjpar.nhsel=3;
378  else if ( hymode_ == "kQJetsOnly") hyjpar.nhsel=4;
379  else hyjpar.nhsel=2;
380 
381  // fraction of soft hydro induced multiplicity
382  hyflow.fpart = fracsoftmult_;
383 
384  // hadron freez-out temperature
385  hyflow.Tf = hadfreeztemp_;
386 
387  // maximum longitudinal collective rapidity
388  hyflow.ylfl = maxlongy_;
389 
390  // maximum transverse collective rapidity
391  hyflow.ytfl = maxtrany_;
392 
393  // shadowing on=1, off=0
394  hyjpar.ishad = shadowingswitch_;
395 
396  // set inelastic nucleon-nucleon cross section
397  hyjpar.sigin = signn_;
398 
399  // angular emitted gluon spectrum selection
400  pyqpar.ianglu = angularspecselector_;
401 
402  // number of active quark flavors in qgp
403  pyqpar.nfu = nquarkflavor_;
404 
405  // initial temperature of QGP
406  pyqpar.T0u = qgpt0_;
407 
408  // proper time of QGP formation
409  pyqpar.tau0u = qgptau0_;
410 
411  // type of medium induced partonic energy loss
413  edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
414  pyqpar.ienglu = 0;
415  } else if ( doradiativeenloss_ ) {
416  edm::LogInfo("HydjetenLoss") << "##### Only RADIATIVE partonic energy loss ON ####";
417  pyqpar.ienglu = 1;
418  } else if ( docollisionalenloss_ ) {
419  edm::LogInfo("HydjetEnLoss") << "##### Only COLLISIONAL partonic energy loss ON ####";
420  pyqpar.ienglu = 2;
421  } else {
422  edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
423  pyqpar.ienglu = 0;
424  }
425  return true;
426 }
427 
428 //_____________________________________________________________________
429 
431 
434 
435  return true;
436 
437 }
438 
439 //_____________________________________________________________________
440 
442 
444  // pythia6Service_->setGeneralParams();
445 
446  // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage
447  const float ra = nuclear_radius();
448  LogInfo("RAScaling")<<"Nuclear radius(RA) = "<<ra;
449  bmin_ /= ra;
450  bmax_ /= ra;
451  bfixed_ /= ra;
452 
453  // hydjet running options
455  // initialize hydjet
456  LogInfo("HYDJETinAction") << "##### Calling HYINIT("<<comenergy<<","<<abeamtarget_<<","
457  <<cflag_<<","<<bmin_<<","<<bmax_<<","<<bfixed_<<","<<nmultiplicity_<<") ####";
459  return true;
460 }
461 
462 bool HydjetHadronizer::declareStableParticles(const std::vector<int>& _pdg )
463 {
464  std::vector<int> pdg = _pdg;
465  for ( size_t i=0; i < pdg.size(); i++ ) {
466  int pyCode = pycomp_( pdg[i] );
467  std::ostringstream pyCard ;
468  pyCard << "MDCY(" << pyCode << ",1)=0";
469  std::cout << pyCard.str() << std::endl;
470  call_pygive( pyCard.str() );
471  }
472  return true;
473 }
474 
475 //________________________________________________________________
477 {
478  const double pi = 3.14159265358979;
479  phi0_ = 2.*pi*gen::pyr_(0) - pi;
480  sinphi0_ = sin(phi0_);
481  cosphi0_ = cos(phi0_);
482 }
483 
484 
485 //________________________________________________________________
487 {
488  return false;
489 }
490 
492 {
493  return true;
494 }
495 
497 {
498  return true;
499 }
500 
502 {
503 }
504 
506 {
507 }
508 
509 const char* HydjetHadronizer::classname() const
510 {
511  return "gen::HydjetHadronizer";
512 }
#define LogDebug(id)
#define HYEVNT
Definition: HydjetWrapper.h:26
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
HepMC::GenEvent * evt
bool docollisionalenloss_
DEFAULT = true.
bool hydjet_init(const edm::ParameterSet &pset)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool call_pygive(const std::string &line)
double npart
Definition: HydjetWrapper.h:49
void add_heavy_ion_rec(HepMC::GenEvent *evt)
#define hyfpar
Definition: HydjetWrapper.h:54
double v[5][pyjets_maxn]
std::auto_ptr< HepMC::GenEvent > & event()
static std::string const input
Definition: EdmProvDump.cc:43
const Double_t pi
static const std::string kPythia6
#define HYINIT
Definition: HydjetWrapper.h:18
T x() const
Cartesian x coordinate.
double p[5][pyjets_maxn]
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const char * classname() const
unsigned int maxEventsToPrint_
def convertStatus
Definition: crabWrap.py:174
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
HepMC::GenParticle * build_hyjet(int index, int barcode)
Pythia6Service * pythia6Service_
#define hyjpar
Definition: HydjetWrapper.h:99
static const std::string kFortranInstance
int pycomp_(int &)
#define pypars
double fracsoftmult_
DEFAULT = true.
HepMC::GenVertex * build_hyjet_vertex(int i, int id)
part
Definition: HCALResponse.h:20
#define hyjets
Definition: HydjetWrapper.h:84
double nuclear_radius() const
bool call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh)
unsigned int pythiaPylistVerbosity_
void setRandomEngine(CLHEP::HepRandomEngine *v)
#define pyqpar
#define hyflow
Definition: HydjetWrapper.h:43
double a
Definition: hdecay.h:121
edm::ParameterSet pset_
unsigned int shadowingswitch_
tuple cout
Definition: gather_cfg.py:121
double pyr_(int *idummy)
bool declareStableParticles(const std::vector< int > &)
bool get_particles(HepMC::GenEvent *evt)
edm::Event & getEDMEvent() const