CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
gen::HydjetHadronizer Class Reference

#include <HydjetHadronizer.h>

Inheritance diagram for gen::HydjetHadronizer:
gen::BaseHadronizer

Public Member Functions

const char * classname () const
 
bool decay ()
 
bool declareSpecialSettings (const std::vector< std::string > &)
 
bool declareStableParticles (const std::vector< int > &)
 
void finalizeEvent ()
 
bool generatePartonsAndHadronize ()
 
bool hadronize ()
 
 HydjetHadronizer (const edm::ParameterSet &)
 
bool initializeForExternalPartons ()
 
bool initializeForInternalPartons ()
 
bool readSettings (int)
 
bool residualDecay ()
 
void statistics ()
 
virtual ~HydjetHadronizer ()
 
- Public Member Functions inherited from gen::BaseHadronizer
 BaseHadronizer (edm::ParameterSet const &ps)
 
void cleanLHE ()
 
void generateLHE (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine, unsigned int ncpu)
 
edm::EventgetEDMEvent () const
 
HepMC::GenEvent * getGenEvent ()
 
GenEventInfoProductgetGenEventInfo ()
 
virtual GenLumiInfoHeadergetGenLumiInfoHeader () const
 
GenRunInfoProductgetGenRunInfo ()
 
const boost::shared_ptr< lhef::LHERunInfo > & getLHERunInfo () const
 
const std::string & gridpackPath () const
 
int randomIndex () const
 
const std::string & randomInitConfigDescription () const
 
void randomizeIndex (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine)
 
void resetEvent (HepMC::GenEvent *event)
 
void resetEventInfo (GenEventInfoProduct *eventInfo)
 
virtual bool select (HepMC::GenEvent *) const
 
void setEDMEvent (edm::Event &event)
 
void setLHEEvent (lhef::LHEEvent *event)
 
void setLHERunInfo (lhef::LHERunInfo *runInfo)
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
 ~BaseHadronizer ()
 

Private Member Functions

void add_heavy_ion_rec (HepMC::GenEvent *evt)
 
HepMC::GenParticle * build_hyjet (int index, int barcode)
 
HepMC::GenVertex * build_hyjet_vertex (int i, int id)
 
bool call_hyinit (double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh)
 
virtual void doSetRandomEngine (CLHEP::HepRandomEngine *v) override
 
virtual std::vector< std::string > const & doSharedResources () const override
 
bool get_particles (HepMC::GenEvent *evt)
 
bool hydjet_init (const edm::ParameterSet &pset)
 
double nuclear_radius () const
 
void rotateEvtPlane ()
 

Private Attributes

double abeamtarget_
 
int angularspecselector_
 
double bfixed_
 
double bmax_
 
double bmin_
 
int cflag_
 
double comenergy
 
double cosphi0_
 
bool docollisionalenloss_
 DEFAULT = true. More...
 
bool doradiativeenloss_
 
bool embedding_
 
HepMC::GenEvent * evt
 
double fracsoftmult_
 DEFAULT = true. More...
 
double hadfreeztemp_
 
std::string hymode_
 
unsigned int maxEventsToPrint_
 
double maxlongy_
 
double maxtrany_
 
int nhard_
 
int nmultiplicity_
 
unsigned int nquarkflavor_
 
int nsoft_
 
int nsub_
 
double phi0_
 
edm::ParameterSet pset_
 
Pythia6Servicepythia6Service_
 
unsigned int pythiaPylistVerbosity_
 
double qgpt0_
 
double qgptau0_
 
bool rotate_
 
unsigned int shadowingswitch_
 
double signn_
 
double sinphi0_
 
edm::InputTag src_
 

Static Private Attributes

static const std::vector< std::string > theSharedResources
 

Additional Inherited Members

- Protected Member Functions inherited from gen::BaseHadronizer
std::auto_ptr< HepMC::GenEvent > & event ()
 
std::auto_ptr< GenEventInfoProduct > & eventInfo ()
 
lhef::LHEEventlheEvent ()
 
lhef::LHERunInfolheRunInfo ()
 
GenRunInfoProductrunInfo ()
 
- Protected Attributes inherited from gen::BaseHadronizer
std::string lheFile_
 
int randomIndex_
 

Detailed Description

Definition at line 37 of file HydjetHadronizer.h.

Constructor & Destructor Documentation

HydjetHadronizer::HydjetHadronizer ( const edm::ParameterSet pset)

Definition at line 59 of file HydjetHadronizer.cc.

References embedding_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, maxEventsToPrint_, pythiaPylistVerbosity_, and src_.

59  :
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")),
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 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
HepMC::GenEvent * evt
bool docollisionalenloss_
DEFAULT = true.
BaseHadronizer(edm::ParameterSet const &ps)
unsigned int maxEventsToPrint_
Pythia6Service * pythia6Service_
double fracsoftmult_
DEFAULT = true.
unsigned int pythiaPylistVerbosity_
edm::ParameterSet pset_
unsigned int shadowingswitch_
HydjetHadronizer::~HydjetHadronizer ( )
virtual

Definition at line 112 of file HydjetHadronizer.cc.

References pythia6Service_.

113 {
114  // destructor
115  call_pystat(1);
116  delete pythia6Service_;
117 }
Pythia6Service * pythia6Service_

Member Function Documentation

void HydjetHadronizer::add_heavy_ion_rec ( HepMC::GenEvent *  evt)
private

Definition at line 128 of file HydjetHadronizer.cc.

References hyfpar, hyjpar, nsub_, nuclear_radius(), and phi0_.

Referenced by generatePartonsAndHadronize().

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 }
HepMC::GenEvent * evt
#define hyfpar
Definition: HydjetWrapper.h:54
#define hyjpar
Definition: HydjetWrapper.h:99
double nuclear_radius() const
HepMC::GenParticle * HydjetHadronizer::build_hyjet ( int  index,
int  barcode 
)
private

Definition at line 156 of file HydjetHadronizer.cc.

References crabWrap::convertStatus(), cosphi0_, GenParticle::GenParticle, hyjets, diffTreeTool::index, gen::p, sinphi0_, and x().

Referenced by get_particles().

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 }
double p[5][pyjets_maxn]
#define hyjets
Definition: HydjetWrapper.h:84
def convertStatus(status)
Definition: crabWrap.py:174
HepMC::GenVertex * HydjetHadronizer::build_hyjet_vertex ( int  i,
int  id 
)
private

Definition at line 181 of file HydjetHadronizer.cc.

References cosphi0_, hyjets, mps_fire::i, sinphi0_, lumiQTWidget::t, and x().

Referenced by get_particles().

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 }
#define hyjets
Definition: HydjetWrapper.h:84
bool HydjetHadronizer::call_hyinit ( double  energy,
double  a,
int  ifb,
double  bmin,
double  bmax,
double  bfix,
int  nh 
)
private

Definition at line 351 of file HydjetHadronizer.cc.

References HYINIT.

Referenced by initializeForInternalPartons().

353 {
354  // initialize hydjet
355 
356  pydatr.mrpy[2]=1;
357  HYINIT(energy,a,ifb,bmin,bmax,bfix,nh);
358  return true;
359 }
#define HYINIT
Definition: HydjetWrapper.h:18
double a
Definition: hdecay.h:121
const char * HydjetHadronizer::classname ( ) const

Definition at line 509 of file HydjetHadronizer.cc.

510 {
511  return "gen::HydjetHadronizer";
512 }
bool HydjetHadronizer::decay ( )

Definition at line 491 of file HydjetHadronizer.cc.

492 {
493  return true;
494 }
bool gen::HydjetHadronizer::declareSpecialSettings ( const std::vector< std::string > &  )
inline

Definition at line 50 of file HydjetHadronizer.h.

References data-class-funcs::classname, myMessageLogger_cff::statistics, and findQualityFiles::v.

50 { return true; }
bool HydjetHadronizer::declareStableParticles ( const std::vector< int > &  _pdg)

Definition at line 462 of file HydjetHadronizer.cc.

References gen::call_pygive(), gather_cfg::cout, mps_fire::i, and gen::pycomp_().

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 }
bool call_pygive(const std::string &line)
int pycomp_(int &)
void HydjetHadronizer::doSetRandomEngine ( CLHEP::HepRandomEngine *  v)
overrideprivatevirtual

Reimplemented from gen::BaseHadronizer.

Definition at line 121 of file HydjetHadronizer.cc.

References pythia6Service_, and gen::Pythia6Service::setRandomEngine().

122 {
124 }
double v[5][pyjets_maxn]
Pythia6Service * pythia6Service_
void setRandomEngine(CLHEP::HepRandomEngine *v)
virtual std::vector<std::string> const& gen::HydjetHadronizer::doSharedResources ( ) const
inlineoverrideprivatevirtual

Reimplemented from gen::BaseHadronizer.

Definition at line 59 of file HydjetHadronizer.h.

59 { return theSharedResources; }
static const std::vector< std::string > theSharedResources
void HydjetHadronizer::finalizeEvent ( )

Definition at line 501 of file HydjetHadronizer.cc.

502 {
503 }
bool HydjetHadronizer::generatePartonsAndHadronize ( )

Definition at line 198 of file HydjetHadronizer.cc.

References add_heavy_ion_rec(), bfixed_, cflag_, funct::cos(), cosphi0_, MillePedeFileConverter_cfg::e, embedding_, gen::BaseHadronizer::event(), edm::errors::EventCorruption, evt, get_particles(), edm::Event::getByLabel(), gen::BaseHadronizer::getEDMEvent(), edm::HepMCProduct::GetEvent(), HYEVNT, hyflow, hyfpar, hyjpar, input, nhard_, nsoft_, nsub_, phi0_, pypars, pyqpar, pythia6Service_, rotate_, rotateEvtPlane(), funct::sin(), sinphi0_, and src_.

199 {
200  Pythia6Service::InstanceWrapper guard(pythia6Service_);
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 }
#define HYEVNT
Definition: HydjetWrapper.h:26
HepMC::GenEvent * evt
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void add_heavy_ion_rec(HepMC::GenEvent *evt)
#define hyfpar
Definition: HydjetWrapper.h:54
std::auto_ptr< HepMC::GenEvent > & event()
static std::string const input
Definition: EdmProvDump.cc:44
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:416
Pythia6Service * pythia6Service_
#define hyjpar
Definition: HydjetWrapper.h:99
#define pypars
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
#define pyqpar
#define hyflow
Definition: HydjetWrapper.h:43
bool get_particles(HepMC::GenEvent *evt)
edm::Event & getEDMEvent() const
bool HydjetHadronizer::get_particles ( HepMC::GenEvent *  evt)
private

Definition at line 268 of file HydjetHadronizer.cc.

References build_hyjet(), build_hyjet_vertex(), GenParticle::GenParticle, hyjets, hyjpar, mps_fire::i, LogDebug, nhard_, nsoft_, nsub_, and HadronAndPartonSelector_cfi::particles.

Referenced by generatePartonsAndHadronize().

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 }
#define LogDebug(id)
HepMC::GenEvent * evt
HepMC::GenParticle * build_hyjet(int index, int barcode)
#define hyjpar
Definition: HydjetWrapper.h:99
HepMC::GenVertex * build_hyjet_vertex(int i, int id)
part
Definition: HCALResponse.h:20
#define hyjets
Definition: HydjetWrapper.h:84
bool HydjetHadronizer::hadronize ( )

Definition at line 486 of file HydjetHadronizer.cc.

487 {
488  return false;
489 }
bool HydjetHadronizer::hydjet_init ( const edm::ParameterSet pset)
private

Definition at line 363 of file HydjetHadronizer.cc.

References angularspecselector_, docollisionalenloss_, doradiativeenloss_, fracsoftmult_, hadfreeztemp_, hyflow, hyjpar, hymode_, maxlongy_, maxtrany_, nquarkflavor_, pyqpar, qgpt0_, qgptau0_, shadowingswitch_, and signn_.

Referenced by initializeForInternalPartons().

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 }
bool docollisionalenloss_
DEFAULT = true.
#define hyjpar
Definition: HydjetWrapper.h:99
double fracsoftmult_
DEFAULT = true.
#define pyqpar
#define hyflow
Definition: HydjetWrapper.h:43
unsigned int shadowingswitch_
bool gen::HydjetHadronizer::initializeForExternalPartons ( )
bool HydjetHadronizer::initializeForInternalPartons ( )

Definition at line 441 of file HydjetHadronizer.cc.

References abeamtarget_, bfixed_, bmax_, bmin_, call_hyinit(), cflag_, comenergy, hydjet_init(), nmultiplicity_, nuclear_radius(), pset_, and pythia6Service_.

441  {
442 
443  Pythia6Service::InstanceWrapper guard(pythia6Service_);
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 }
bool hydjet_init(const edm::ParameterSet &pset)
Pythia6Service * pythia6Service_
double nuclear_radius() const
bool call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh)
edm::ParameterSet pset_
double HydjetHadronizer::nuclear_radius ( ) const
inlineprivate

Definition at line 132 of file HydjetHadronizer.h.

References funct::pow().

Referenced by add_heavy_ion_rec(), and initializeForInternalPartons().

133  {
134  // Return the nuclear radius derived from the
135  // beam/target atomic mass number.
136 
137  return 1.15 * pow((double)abeamtarget_, 1./3.);
138  }
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
bool HydjetHadronizer::readSettings ( int  )

Definition at line 430 of file HydjetHadronizer.cc.

References pythia6Service_, and gen::Pythia6Service::setGeneralParams().

430  {
431 
432  Pythia6Service::InstanceWrapper guard(pythia6Service_);
434 
435  return true;
436 
437 }
Pythia6Service * pythia6Service_
bool HydjetHadronizer::residualDecay ( )

Definition at line 496 of file HydjetHadronizer.cc.

497 {
498  return true;
499 }
void HydjetHadronizer::rotateEvtPlane ( )
private

Definition at line 476 of file HydjetHadronizer.cc.

References funct::cos(), cosphi0_, phi0_, pi, gen::pyr_(), funct::sin(), and sinphi0_.

Referenced by generatePartonsAndHadronize().

477 {
478  const double pi = 3.14159265358979;
479  phi0_ = 2.*pi*gen::pyr_(0) - pi;
480  sinphi0_ = sin(phi0_);
481  cosphi0_ = cos(phi0_);
482 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const Double_t pi
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double pyr_(int *idummy)
void HydjetHadronizer::statistics ( )

Definition at line 505 of file HydjetHadronizer.cc.

506 {
507 }

Member Data Documentation

double gen::HydjetHadronizer::abeamtarget_
private

Definition at line 75 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

int gen::HydjetHadronizer::angularspecselector_
private

Definition at line 76 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::bfixed_
private

Definition at line 77 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and initializeForInternalPartons().

double gen::HydjetHadronizer::bmax_
private

Definition at line 78 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

double gen::HydjetHadronizer::bmin_
private

Definition at line 80 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

int gen::HydjetHadronizer::cflag_
private

Definition at line 82 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and initializeForInternalPartons().

double gen::HydjetHadronizer::comenergy
private

Definition at line 86 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

double gen::HydjetHadronizer::cosphi0_
private
bool gen::HydjetHadronizer::docollisionalenloss_
private

DEFAULT = true.

Definition at line 88 of file HydjetHadronizer.h.

Referenced by hydjet_init().

bool gen::HydjetHadronizer::doradiativeenloss_
private

Definition at line 87 of file HydjetHadronizer.h.

Referenced by hydjet_init().

bool gen::HydjetHadronizer::embedding_
private

Definition at line 85 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and HydjetHadronizer().

HepMC::GenEvent* gen::HydjetHadronizer::evt
private

Definition at line 73 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize().

double gen::HydjetHadronizer::fracsoftmult_
private

DEFAULT = true.

Definition at line 89 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::hadfreeztemp_
private

Definition at line 95 of file HydjetHadronizer.h.

Referenced by hydjet_init().

std::string gen::HydjetHadronizer::hymode_
private

Definition at line 97 of file HydjetHadronizer.h.

Referenced by hydjet_init().

unsigned int gen::HydjetHadronizer::maxEventsToPrint_
private

Definition at line 98 of file HydjetHadronizer.h.

Referenced by HydjetHadronizer().

double gen::HydjetHadronizer::maxlongy_
private

Definition at line 99 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::maxtrany_
private

Definition at line 102 of file HydjetHadronizer.h.

Referenced by hydjet_init().

int gen::HydjetHadronizer::nhard_
private

Definition at line 106 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and get_particles().

int gen::HydjetHadronizer::nmultiplicity_
private

Definition at line 107 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

unsigned int gen::HydjetHadronizer::nquarkflavor_
private

Definition at line 110 of file HydjetHadronizer.h.

Referenced by hydjet_init().

int gen::HydjetHadronizer::nsoft_
private

Definition at line 109 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and get_particles().

int gen::HydjetHadronizer::nsub_
private
double gen::HydjetHadronizer::phi0_
private
edm::ParameterSet gen::HydjetHadronizer::pset_
private

Definition at line 74 of file HydjetHadronizer.h.

Referenced by initializeForInternalPartons().

Pythia6Service* gen::HydjetHadronizer::pythia6Service_
private
unsigned int gen::HydjetHadronizer::pythiaPylistVerbosity_
private

number of active quark flavors in qgp DEFAULT=0; allowed values: 0,1,2,3.

Definition at line 112 of file HydjetHadronizer.h.

Referenced by HydjetHadronizer().

double gen::HydjetHadronizer::qgpt0_
private

Definition at line 113 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::qgptau0_
private

Definition at line 115 of file HydjetHadronizer.h.

Referenced by hydjet_init().

bool gen::HydjetHadronizer::rotate_
private

Definition at line 121 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize().

unsigned int gen::HydjetHadronizer::shadowingswitch_
private

Definition at line 123 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::signn_
private

Definition at line 125 of file HydjetHadronizer.h.

Referenced by hydjet_init().

double gen::HydjetHadronizer::sinphi0_
private
edm::InputTag gen::HydjetHadronizer::src_
private

Definition at line 129 of file HydjetHadronizer.h.

Referenced by generatePartonsAndHadronize(), and HydjetHadronizer().

const std::vector< std::string > HydjetHadronizer::theSharedResources
staticprivate