#include <HydjetHadronizer.h>
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 () |
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) |
bool | get_particles (HepMC::GenEvent *evt) |
bool | hydjet_init (const edm::ParameterSet &pset) |
double | nuclear_radius () const |
void | rotateEvtPlane () |
Private Attributes | |
double | abeamtarget_ |
double | bfixed_ |
double | bmax_ |
double | bmin_ |
int | cflag_ |
double | comenergy |
double | cosphi0_ |
bool | docollisionalenloss_ |
DEFAULT = true. | |
bool | doradiativeenloss_ |
bool | embedding_ |
HepMC::GenEvent * | evt |
double | fracsoftmult_ |
DEFAULT = true. | |
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_ |
Pythia6Service * | pythia6Service_ |
unsigned int | pythiaPylistVerbosity_ |
double | qgpt0_ |
double | qgptau0_ |
bool | rotate_ |
unsigned int | shadowingswitch_ |
double | signn_ |
double | sinphi0_ |
edm::InputTag | src_ |
Definition at line 34 of file HydjetHadronizer.h.
HydjetHadronizer::HydjetHadronizer | ( | const edm::ParameterSet & | pset | ) |
Definition at line 54 of file HydjetHadronizer.cc.
References embedding_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, maxEventsToPrint_, pythiaPylistVerbosity_, and src_.
: BaseHadronizer(pset), evt(0), pset_(pset), abeamtarget_(pset.getParameter<double>("aBeamTarget")), bfixed_(pset.getParameter<double>("bFixed")), bmax_(pset.getParameter<double>("bMax")), bmin_(pset.getParameter<double>("bMin")), cflag_(pset.getParameter<int>("cFlag")), embedding_(pset.getParameter<bool>("embeddingMode")), comenergy(pset.getParameter<double>("comEnergy")), doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")), docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")), fracsoftmult_(pset.getParameter<double>("fracSoftMultiplicity")), hadfreeztemp_(pset.getParameter<double>("hadronFreezoutTemperature")), hymode_(pset.getParameter<string>("hydjetMode")), maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)), maxlongy_(pset.getParameter<double>("maxLongitudinalRapidity")), maxtrany_(pset.getParameter<double>("maxTransverseRapidity")), nsub_(0), nhard_(0), nmultiplicity_(pset.getParameter<int>("nMultiplicity")), nsoft_(0), nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")), pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)), qgpt0_(pset.getParameter<double>("qgpInitialTemperature")), qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")), phi0_(0.), sinphi0_(0.), cosphi0_(1.), rotate_(pset.getParameter<bool>("rotateEventPlane")), shadowingswitch_(pset.getParameter<int>("shadowingSwitch")), signn_(pset.getParameter<double>("sigmaInelNN")), pythia6Service_(new Pythia6Service(pset)) { // Default constructor // PYLIST Verbosity Level // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13 pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0); LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_; //Max number of events printed on verbosity level maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0); LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_; if(embedding_) src_ = pset.getParameter<edm::InputTag>("backgroundLabel"); }
HydjetHadronizer::~HydjetHadronizer | ( | ) | [virtual] |
Definition at line 106 of file HydjetHadronizer.cc.
References pythia6Service_.
{ // destructor call_pystat(1); delete pythia6Service_; }
void HydjetHadronizer::add_heavy_ion_rec | ( | HepMC::GenEvent * | evt | ) | [private] |
Definition at line 115 of file HydjetHadronizer.cc.
References hyfpar, hyjpar, npart, nsub_, nuclear_radius(), and phi0_.
Referenced by generatePartonsAndHadronize().
{ // heavy ion record in the final CMSSW Event double npart = hyfpar.npart; int nproj = static_cast<int>(npart / 2); int ntarg = static_cast<int>(npart - nproj); HepMC::HeavyIon* hi = new HepMC::HeavyIon( nsub_, // Ncoll_hard/N of SubEvents nproj, // Npart_proj ntarg, // Npart_targ static_cast<int>(hyfpar.nbcol), // Ncoll 0, // spectator_neutrons 0, // spectator_protons 0, // N_Nwounded_collisions 0, // Nwounded_N_collisions 0, // Nwounded_Nwounded_collisions hyfpar.bgen * nuclear_radius(), // impact_parameter in [fm] phi0_, // event_plane_angle 0, // eccentricity hyjpar.sigin // sigma_inel_NN ); evt->set_heavy_ion(*hi); delete hi; }
HepMC::GenParticle * HydjetHadronizer::build_hyjet | ( | int | index, |
int | barcode | ||
) | [private] |
Definition at line 143 of file HydjetHadronizer.cc.
References crabWrap::convertStatus(), cosphi0_, configurableAnalysis::GenParticle, hyjets, getHLTprescales::index, gen::p, sinphi0_, x, and detailsBasic3DVector::y.
Referenced by get_particles().
{ // Build particle object corresponding to index in hyjets (soft+hard) double x0 = hyjets.phj[0][index]; double y0 = hyjets.phj[1][index]; double x = x0*cosphi0_-y0*sinphi0_; double y = y0*cosphi0_+x0*sinphi0_; HepMC::GenParticle* p = new HepMC::GenParticle( HepMC::FourVector(x, // px y, // py hyjets.phj[2][index], // pz hyjets.phj[3][index]), // E hyjets.khj[1][index], // id convertStatus(hyjets.khj[0][index] // status ) ); p->suggest_barcode(barcode); return p; }
HepMC::GenVertex * HydjetHadronizer::build_hyjet_vertex | ( | int | i, |
int | id | ||
) | [private] |
Definition at line 168 of file HydjetHadronizer.cc.
References cosphi0_, hyjets, i, sinphi0_, lumiQTWidget::t, x, detailsBasic3DVector::y, and z.
Referenced by get_particles().
{ // build verteces for the hyjets stored events double x0=hyjets.vhj[0][i]; double y0=hyjets.vhj[1][i]; double x = x0*cosphi0_-y0*sinphi0_; double y = y0*cosphi0_+x0*sinphi0_; double z=hyjets.vhj[2][i]; double t=hyjets.vhj[4][i]; HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t),id); return vertex; }
bool HydjetHadronizer::call_hyinit | ( | double | energy, |
double | a, | ||
int | ifb, | ||
double | bmin, | ||
double | bmax, | ||
double | bfix, | ||
int | nh | ||
) | [private] |
Definition at line 338 of file HydjetHadronizer.cc.
References HYINIT.
Referenced by initializeForInternalPartons().
const char * HydjetHadronizer::classname | ( | ) | const |
Definition at line 492 of file HydjetHadronizer.cc.
{ return "gen::HydjetHadronizer"; }
bool HydjetHadronizer::decay | ( | ) |
Definition at line 474 of file HydjetHadronizer.cc.
{ return true; }
bool gen::HydjetHadronizer::declareSpecialSettings | ( | const std::vector< std::string > | ) | [inline] |
Definition at line 47 of file HydjetHadronizer.h.
{ return true; }
bool HydjetHadronizer::declareStableParticles | ( | const std::vector< int > | pdg | ) |
Definition at line 446 of file HydjetHadronizer.cc.
References gen::call_pygive(), gather_cfg::cout, i, and gen::pycomp_().
void HydjetHadronizer::finalizeEvent | ( | ) |
Definition at line 484 of file HydjetHadronizer.cc.
{ }
bool HydjetHadronizer::generatePartonsAndHadronize | ( | ) |
Definition at line 185 of file HydjetHadronizer.cc.
References add_heavy_ion_rec(), bfixed_, cflag_, funct::cos(), cosphi0_, alignCSCRings::e, embedding_, gen::BaseHadronizer::event(), edm::errors::EventCorruption, evt, except, get_particles(), edm::Event::getByLabel(), gen::BaseHadronizer::getEDMEvent(), HYEVNT, hyflow, hyfpar, hyjpar, LaserDQM_cfg::input, nhard_, nsoft_, nsub_, phi0_, pypars, pyqpar, pythia6Service_, rotate_, rotateEvtPlane(), funct::sin(), sinphi0_, and src_.
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); // generate single event if(embedding_){ cflag_ = 0; const edm::Event& e = getEDMEvent(); Handle<HepMCProduct> input; e.getByLabel(src_,input); const HepMC::GenEvent * inev = input->GetEvent(); const HepMC::HeavyIon* hi = inev->heavy_ion(); if(hi){ bfixed_ = hi->impact_parameter(); phi0_ = hi->event_plane_angle(); sinphi0_ = sin(phi0_); cosphi0_ = cos(phi0_); }else{ LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!"; } }else if(rotate_) rotateEvtPlane(); nsoft_ = 0; nhard_ = 0; edm::LogInfo("HYDJETmode") << "##### HYDJET nhsel = " << hyjpar.nhsel; edm::LogInfo("HYDJETfpart") << "##### HYDJET fpart = " << hyflow.fpart; edm::LogInfo("HYDJETtf") << "##### HYDJET hadron freez-out temp, Tf = " << hyflow.Tf; edm::LogInfo("HYDJETinTemp") << "##### HYDJET: QGP init temperature, T0 ="<<pyqpar.T0u; edm::LogInfo("HYDJETinTau") << "##### HYDJET: QGP formation time,tau0 ="<<pyqpar.tau0u; // generate a HYDJET event int ntry = 0; while(nsoft_ == 0 && nhard_ == 0){ if(ntry > 100){ edm::LogError("HydjetEmptyEvent") << "##### HYDJET: No Particles generated, Number of tries ="<<ntry; // Throw an exception. Use the EventCorruption exception since it maps onto SkipEvent // which is what we want to do here. std::ostringstream sstr; sstr << "HydjetHadronizerProducer: No particles generated after " << ntry << " tries.\n"; edm::Exception except(edm::errors::EventCorruption, sstr.str()); throw except; } else { HYEVNT(); nsoft_ = hyfpar.nhyd; nsub_ = hyjpar.njet; nhard_ = hyfpar.npyt; ++ntry; } } if(hyjpar.nhsel < 3) nsub_++; // event information HepMC::GenEvent *evt = new HepMC::GenEvent(); if(nhard_>0 || nsoft_>0) get_particles(evt); evt->set_signal_process_id(pypars.msti[0]); // type of the process evt->set_event_scale(pypars.pari[16]); // Q^2 add_heavy_ion_rec(evt); event().reset(evt); return true; }
bool HydjetHadronizer::get_particles | ( | HepMC::GenEvent * | evt | ) | [private] |
Definition at line 255 of file HydjetHadronizer.cc.
References build_hyjet(), build_hyjet_vertex(), configurableAnalysis::GenParticle, hyjets, hyjpar, i, LogDebug, nhard_, nsoft_, and nsub_.
Referenced by generatePartonsAndHadronize().
{ // Hard particles. The first nhard_ lines from hyjets array. // Pythia/Pyquen sub-events (sub-collisions) for a given event // Return T/F if success/failure // Create particles from lujet entries, assign them into vertices and // put the vertices in the GenEvent, for each SubEvent // The SubEvent information is kept by storing indeces of main vertices // of subevents as a vector in GenHIEvent. LogDebug("SubEvent")<< "Number of sub events "<<nsub_; LogDebug("Hydjet")<<"Number of hard events "<<hyjpar.njet; LogDebug("Hydjet")<<"Number of hard particles "<<nhard_; LogDebug("Hydjet")<<"Number of soft particles "<<nsoft_; vector<HepMC::GenVertex*> sub_vertices(nsub_); int ihy = 0; for(int isub=0;isub<nsub_;isub++){ LogDebug("SubEvent") <<"Sub Event ID : "<<isub; int sub_up = (isub+1)*50000; // Upper limit in mother index, determining the range of Sub-Event vector<HepMC::GenParticle*> particles; vector<int> mother_ids; vector<HepMC::GenVertex*> prods; sub_vertices[isub] = new HepMC::GenVertex(HepMC::FourVector(0,0,0,0),isub); evt->add_vertex(sub_vertices[isub]); if(!evt->signal_process_vertex()) evt->set_signal_process_vertex(sub_vertices[isub]); while(ihy<nhard_+nsoft_ && (hyjets.khj[2][ihy] < sub_up || ihy > nhard_ )){ particles.push_back(build_hyjet(ihy,ihy+1)); prods.push_back(build_hyjet_vertex(ihy,isub)); mother_ids.push_back(hyjets.khj[2][ihy]); LogDebug("DecayChain")<<"Mother index : "<<hyjets.khj[2][ihy]; ihy++; } //Produce Vertices and add them to the GenEvent. Remember that GenParticles are adopted by //GenVertex and GenVertex is adopted by GenEvent. LogDebug("Hydjet")<<"Number of particles in vector "<<particles.size(); for (unsigned int i = 0; i<particles.size(); i++) { HepMC::GenParticle* part = particles[i]; //The Fortran code is modified to preserve mother id info, by seperating the beginning //mother indices of successive subevents by 5000 int mid = mother_ids[i]-isub*50000-1; LogDebug("DecayChain")<<"Particle "<<i; LogDebug("DecayChain")<<"Mother's ID "<<mid; LogDebug("DecayChain")<<"Particle's PDG ID "<<part->pdg_id(); if(mid <= 0){ sub_vertices[isub]->add_particle_out(part); continue; } if(mid > 0){ HepMC::GenParticle* mother = particles[mid]; LogDebug("DecayChain")<<"Mother's PDG ID "<<mother->pdg_id(); HepMC::GenVertex* prod_vertex = mother->end_vertex(); if(!prod_vertex){ prod_vertex = prods[i]; prod_vertex->add_particle_in(mother); evt->add_vertex(prod_vertex); prods[i]=0; // mark to protect deletion } prod_vertex->add_particle_out(part); } } // cleanup vertices not assigned to evt for (unsigned int i = 0; i<prods.size(); i++) { if(prods[i]) delete prods[i]; } } return true; }
bool HydjetHadronizer::hadronize | ( | ) |
Definition at line 469 of file HydjetHadronizer.cc.
{ return false; }
bool HydjetHadronizer::hydjet_init | ( | const edm::ParameterSet & | pset | ) | [private] |
Definition at line 350 of file HydjetHadronizer.cc.
References docollisionalenloss_, doradiativeenloss_, fracsoftmult_, hadfreeztemp_, hyflow, hyjpar, hymode_, maxlongy_, maxtrany_, nquarkflavor_, pyqpar, qgpt0_, qgptau0_, shadowingswitch_, and signn_.
Referenced by initializeForInternalPartons().
{ // set hydjet options // hydjet running mode mode // kHydroOnly --- nhsel=0 jet production off (pure HYDRO event), nhsel=0 // kHydroJets --- nhsle=1 jet production on, jet quenching off (HYDRO+njet*PYTHIA events) // kHydroQJet --- nhsel=2 jet production & jet quenching on (HYDRO+njet*PYQUEN events) // kJetsOnly --- nhsel=3 jet production on, jet quenching off, HYDRO off (njet*PYTHIA events) // kQJetsOnly --- nhsel=4 jet production & jet quenching on, HYDRO off (njet*PYQUEN events) if(hymode_ == "kHydroOnly") hyjpar.nhsel=0; else if ( hymode_ == "kHydroJets") hyjpar.nhsel=1; else if ( hymode_ == "kHydroQJets") hyjpar.nhsel=2; else if ( hymode_ == "kJetsOnly") hyjpar.nhsel=3; else if ( hymode_ == "kQJetsOnly") hyjpar.nhsel=4; else hyjpar.nhsel=2; // fraction of soft hydro induced multiplicity hyflow.fpart = fracsoftmult_; // hadron freez-out temperature hyflow.Tf = hadfreeztemp_; // maximum longitudinal collective rapidity hyflow.ylfl = maxlongy_; // maximum transverse collective rapidity hyflow.ytfl = maxtrany_; // shadowing on=1, off=0 hyjpar.ishad = shadowingswitch_; // set inelastic nucleon-nucleon cross section hyjpar.sigin = signn_; // number of active quark flavors in qgp pyqpar.nfu = nquarkflavor_; // initial temperature of QGP pyqpar.T0u = qgpt0_; // proper time of QGP formation pyqpar.tau0u = qgptau0_; // type of medium induced partonic energy loss if( doradiativeenloss_ && docollisionalenloss_ ){ edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####"; pyqpar.ienglu = 0; } else if ( doradiativeenloss_ ) { edm::LogInfo("HydjetenLoss") << "##### Only RADIATIVE partonic energy loss ON ####"; pyqpar.ienglu = 1; } else if ( docollisionalenloss_ ) { edm::LogInfo("HydjetEnLoss") << "##### Only COLLISIONAL partonic energy loss ON ####"; pyqpar.ienglu = 2; } else { edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####"; pyqpar.ienglu = 0; } return true; }
bool gen::HydjetHadronizer::initializeForExternalPartons | ( | ) |
bool HydjetHadronizer::initializeForInternalPartons | ( | ) |
Definition at line 425 of file HydjetHadronizer.cc.
References abeamtarget_, bfixed_, bmax_, bmin_, call_hyinit(), cflag_, comenergy, hydjet_init(), nmultiplicity_, nuclear_radius(), pset_, and pythia6Service_.
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); // pythia6Service_->setGeneralParams(); // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage const float ra = nuclear_radius(); LogInfo("RAScaling")<<"Nuclear radius(RA) = "<<ra; bmin_ /= ra; bmax_ /= ra; bfixed_ /= ra; // hydjet running options hydjet_init(pset_); // initialize hydjet LogInfo("HYDJETinAction") << "##### Calling HYINIT("<<comenergy<<","<<abeamtarget_<<"," <<cflag_<<","<<bmin_<<","<<bmax_<<","<<bfixed_<<","<<nmultiplicity_<<") ####"; call_hyinit(comenergy,abeamtarget_,cflag_,bmin_,bmax_,bfixed_,nmultiplicity_); return true; }
double HydjetHadronizer::nuclear_radius | ( | ) | const [inline, private] |
Definition at line 122 of file HydjetHadronizer.h.
References abeamtarget_, and funct::pow().
Referenced by add_heavy_ion_rec(), and initializeForInternalPartons().
{ // Return the nuclear radius derived from the // beam/target atomic mass number. return 1.15 * pow((double)abeamtarget_, 1./3.); }
bool HydjetHadronizer::readSettings | ( | int | ) |
Definition at line 414 of file HydjetHadronizer.cc.
References pythia6Service_, and gen::Pythia6Service::setGeneralParams().
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); pythia6Service_->setGeneralParams(); return true; }
bool HydjetHadronizer::residualDecay | ( | ) |
Definition at line 479 of file HydjetHadronizer.cc.
{ return true; }
void HydjetHadronizer::rotateEvtPlane | ( | ) | [private] |
Definition at line 459 of file HydjetHadronizer.cc.
References funct::cos(), cosphi0_, phi0_, pi, gen::pyr_(), funct::sin(), and sinphi0_.
Referenced by generatePartonsAndHadronize().
void HydjetHadronizer::statistics | ( | ) |
Definition at line 488 of file HydjetHadronizer.cc.
{ }
double gen::HydjetHadronizer::abeamtarget_ [private] |
Definition at line 66 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons(), and nuclear_radius().
double gen::HydjetHadronizer::bfixed_ [private] |
Definition at line 67 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and initializeForInternalPartons().
double gen::HydjetHadronizer::bmax_ [private] |
Definition at line 68 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons().
double gen::HydjetHadronizer::bmin_ [private] |
Definition at line 70 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons().
int gen::HydjetHadronizer::cflag_ [private] |
Definition at line 72 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and initializeForInternalPartons().
double gen::HydjetHadronizer::comenergy [private] |
Definition at line 76 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons().
double gen::HydjetHadronizer::cosphi0_ [private] |
Definition at line 110 of file HydjetHadronizer.h.
Referenced by build_hyjet(), build_hyjet_vertex(), generatePartonsAndHadronize(), and rotateEvtPlane().
bool gen::HydjetHadronizer::docollisionalenloss_ [private] |
bool gen::HydjetHadronizer::doradiativeenloss_ [private] |
Definition at line 77 of file HydjetHadronizer.h.
Referenced by hydjet_init().
bool gen::HydjetHadronizer::embedding_ [private] |
Definition at line 75 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and HydjetHadronizer().
HepMC::GenEvent* gen::HydjetHadronizer::evt [private] |
Definition at line 64 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize().
double gen::HydjetHadronizer::fracsoftmult_ [private] |
double gen::HydjetHadronizer::hadfreeztemp_ [private] |
Definition at line 85 of file HydjetHadronizer.h.
Referenced by hydjet_init().
std::string gen::HydjetHadronizer::hymode_ [private] |
Definition at line 87 of file HydjetHadronizer.h.
Referenced by hydjet_init().
unsigned int gen::HydjetHadronizer::maxEventsToPrint_ [private] |
Definition at line 88 of file HydjetHadronizer.h.
Referenced by HydjetHadronizer().
double gen::HydjetHadronizer::maxlongy_ [private] |
Definition at line 89 of file HydjetHadronizer.h.
Referenced by hydjet_init().
double gen::HydjetHadronizer::maxtrany_ [private] |
Definition at line 92 of file HydjetHadronizer.h.
Referenced by hydjet_init().
int gen::HydjetHadronizer::nhard_ [private] |
Definition at line 96 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and get_particles().
int gen::HydjetHadronizer::nmultiplicity_ [private] |
Definition at line 97 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons().
unsigned int gen::HydjetHadronizer::nquarkflavor_ [private] |
Definition at line 100 of file HydjetHadronizer.h.
Referenced by hydjet_init().
int gen::HydjetHadronizer::nsoft_ [private] |
Definition at line 99 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and get_particles().
int gen::HydjetHadronizer::nsub_ [private] |
Definition at line 95 of file HydjetHadronizer.h.
Referenced by add_heavy_ion_rec(), generatePartonsAndHadronize(), and get_particles().
double gen::HydjetHadronizer::phi0_ [private] |
Definition at line 108 of file HydjetHadronizer.h.
Referenced by add_heavy_ion_rec(), generatePartonsAndHadronize(), and rotateEvtPlane().
Definition at line 65 of file HydjetHadronizer.h.
Referenced by initializeForInternalPartons().
Definition at line 118 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), initializeForInternalPartons(), readSettings(), and ~HydjetHadronizer().
unsigned int gen::HydjetHadronizer::pythiaPylistVerbosity_ [private] |
number of active quark flavors in qgp DEFAULT=0; allowed values: 0,1,2,3.
Definition at line 102 of file HydjetHadronizer.h.
Referenced by HydjetHadronizer().
double gen::HydjetHadronizer::qgpt0_ [private] |
Definition at line 103 of file HydjetHadronizer.h.
Referenced by hydjet_init().
double gen::HydjetHadronizer::qgptau0_ [private] |
Definition at line 105 of file HydjetHadronizer.h.
Referenced by hydjet_init().
bool gen::HydjetHadronizer::rotate_ [private] |
Definition at line 111 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize().
unsigned int gen::HydjetHadronizer::shadowingswitch_ [private] |
Definition at line 113 of file HydjetHadronizer.h.
Referenced by hydjet_init().
double gen::HydjetHadronizer::signn_ [private] |
Definition at line 115 of file HydjetHadronizer.h.
Referenced by hydjet_init().
double gen::HydjetHadronizer::sinphi0_ [private] |
Definition at line 109 of file HydjetHadronizer.h.
Referenced by build_hyjet(), build_hyjet_vertex(), generatePartonsAndHadronize(), and rotateEvtPlane().
edm::InputTag gen::HydjetHadronizer::src_ [private] |
Definition at line 119 of file HydjetHadronizer.h.
Referenced by generatePartonsAndHadronize(), and HydjetHadronizer().