#include <PyquenHadronizer.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 () |
bool | initializeForExternalPartons () |
bool | initializeForInternalPartons () |
PyquenHadronizer (const edm::ParameterSet &) | |
bool | readSettings (int) |
bool | residualDecay () |
virtual bool | select (HepMC::GenEvent *evtTry) const |
void | statistics () |
virtual | ~PyquenHadronizer () |
Private Member Functions | |
void | add_heavy_ion_rec (HepMC::GenEvent *evt) |
const char * | nucleon () |
bool | pyqpythia_init (const edm::ParameterSet &pset) |
bool | pyquen_init (const edm::ParameterSet &pset) |
void | rotateEvtPlane (HepMC::GenEvent *evt, double angle) |
Private Attributes | |
double | abeamtarget_ |
unsigned int | angularspecselector_ |
beam/target atomic mass number | |
double | bfixed_ |
max impact param (fm); valid only if cflag_!=0 | |
double | bmax_ |
min impact param (fm); valid only if cflag_!=0 | |
double | bmin_ |
int | cflag_ |
fixed impact param (fm); valid only if cflag_=0 | |
double | comenergy |
centrality flag =0 fixed impact param, <>0 minbias | |
bool | docollisionalenloss_ |
DEFAULT = true. | |
bool | doIsospin_ |
DEFAULT = true. | |
bool | doquench_ |
collision energy | |
bool | doradiativeenloss_ |
if true perform quenching (default = true) | |
bool | embedding_ |
double | evtPlane_ |
std::string | filterType_ |
unsigned int | maxEventsToPrint_ |
unsigned int | nquarkflavor_ |
Proton fraction in the nucleus. | |
double | pfrac_ |
int | protonSide_ |
Run n&p with proper ratios; if false, only p+p collisions. | |
edm::ParameterSet | pset_ |
Pythia6Service * | pythia6Service_ |
bool | pythiaHepMCVerbosity_ |
Events to print if verbosity. | |
unsigned int | pythiaPylistVerbosity_ |
HepMC verbosity flag. | |
double | qgpt0_ |
double | qgptau0_ |
BaseHiGenEvtSelector * | selector_ |
edm::InputTag | src_ |
Pythia PYLIST Verbosity flag. |
Definition at line 25 of file PyquenHadronizer.h.
PyquenHadronizer::PyquenHadronizer | ( | const edm::ParameterSet & | pset | ) |
Definition at line 35 of file PyquenHadronizer.cc.
References cflag_, embedding_, filterType_, reco::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, maxEventsToPrint_, pythiaHepMCVerbosity_, pythiaPylistVerbosity_, selector_, and src_.
: BaseHadronizer(pset), pset_(pset), abeamtarget_(pset.getParameter<double>("aBeamTarget")), angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")), bmin_(pset.getParameter<double>("bMin")), bmax_(pset.getParameter<double>("bMax")), bfixed_(pset.getParameter<double>("bFixed")), cflag_(pset.getParameter<int>("cFlag")), comenergy(pset.getParameter<double>("comEnergy")), doquench_(pset.getParameter<bool>("doQuench")), doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")), docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")), doIsospin_(pset.getParameter<bool>("doIsospin")), protonSide_(pset.getUntrackedParameter<int>("protonSide",0)), embedding_(pset.getParameter<bool>("embeddingMode")), evtPlane_(0), nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")), qgpt0_(pset.getParameter<double>("qgpInitialTemperature")), qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")), maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint",1)), pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)), pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)), pythia6Service_(new Pythia6Service(pset)), filterType_(pset.getUntrackedParameter<string>("filterType","None")) { // Verbosity Level // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13 LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl; // HepMC event verbosity Level pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false); LogDebug("HepMCverbosity") << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl; //Max number of events printed on verbosity level maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0); LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl; if(embedding_){ cflag_ = 0; src_ = pset.getParameter<InputTag>("backgroundLabel"); } selector_ = HiGenEvtSelectorFactory::get(filterType_,pset); }
PyquenHadronizer::~PyquenHadronizer | ( | ) | [virtual] |
Definition at line 82 of file PyquenHadronizer.cc.
References pythia6Service_.
{ // distructor call_pystat(1); delete pythia6Service_; }
void PyquenHadronizer::add_heavy_ion_rec | ( | HepMC::GenEvent * | evt | ) | [private] |
Definition at line 93 of file PyquenHadronizer.cc.
References evtPlane_, and plfpar.
Referenced by generatePartonsAndHadronize().
{ HepMC::HeavyIon *hi = new HepMC::HeavyIon( 1, // Ncoll_hard -1, // Npart_proj -1, // Npart_targ 1, // Ncoll -1, // spectator_neutrons -1, // spectator_protons -1, // N_Nwounded_collisions -1, // Nwounded_N_collisions -1, // Nwounded_Nwounded_collisions plfpar.bgen, // impact_parameter in [fm] evtPlane_, // event_plane_angle 0, // eccentricity -1 // sigma_inel_NN ); evt->set_heavy_ion(*hi); delete hi; }
const char * PyquenHadronizer::classname | ( | ) | const |
Definition at line 347 of file PyquenHadronizer.cc.
{ return "gen::PyquenHadronizer"; }
bool PyquenHadronizer::decay | ( | ) |
Definition at line 330 of file PyquenHadronizer.cc.
{ return true; }
bool gen::PyquenHadronizer::declareSpecialSettings | ( | const std::vector< std::string > | ) | [inline] |
Definition at line 39 of file PyquenHadronizer.h.
{ return true; }
bool PyquenHadronizer::declareStableParticles | ( | const std::vector< int > | pdg | ) |
Definition at line 306 of file PyquenHadronizer.cc.
References gen::call_pygive(), gather_cfg::cout, i, and gen::pycomp_().
void PyquenHadronizer::finalizeEvent | ( | ) |
Definition at line 340 of file PyquenHadronizer.cc.
{ }
bool PyquenHadronizer::generatePartonsAndHadronize | ( | ) |
Definition at line 117 of file PyquenHadronizer.cc.
References abeamtarget_, add_heavy_ion_rec(), bfixed_, bmax_, bmin_, cflag_, comenergy, doIsospin_, doquench_, alignCSCRings::e, embedding_, gen::BaseHadronizer::event(), evtPlane_, edm::Event::getByLabel(), gen::BaseHadronizer::getEDMEvent(), hepevtio, LaserDQM_cfg::input, nucleon(), protonSide_, gen::pyexec_(), pypars, PYQUEN, pythia6Service_, rotateEvtPlane(), and src_.
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); // Not possible to retrieve impact paramter and event plane info // at this part, need to overwrite filter() in // PyquenGeneratorFilter const edm::Event& e = getEDMEvent(); if(embedding_){ 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(); evtPlane_ = hi->event_plane_angle(); }else{ LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!"; } } // Generate PYQUEN event // generate single partonic PYTHIA jet event // Take into account whether it's a nn or pp or pn interaction if(doIsospin_){ string projN = "p"; string targN = "p"; if(protonSide_ != 1) projN = nucleon(); if(protonSide_ != 2) targN = nucleon(); call_pyinit("CMS", projN.data(), targN.data(), comenergy); } call_pyevnt(); // call PYQUEN to apply parton rescattering and energy loss // if doQuench=FALSE, it is pure PYTHIA if( doquench_ ){ PYQUEN(abeamtarget_,cflag_,bfixed_,bmin_,bmax_); edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####"; } else { edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####"; } // call PYTHIA to finish the hadronization pyexec_(); // fill the HEPEVT with the PYJETS event record call_pyhepc(1); // event information HepMC::GenEvent* evt = hepevtio.read_next_event(); evt->set_signal_process_id(pypars.msti[0]); // type of the process evt->set_event_scale(pypars.pari[16]); // Q^2 if(embedding_) rotateEvtPlane(evt,evtPlane_); add_heavy_ion_rec(evt); event().reset(evt); return true; }
bool PyquenHadronizer::hadronize | ( | ) |
Definition at line 325 of file PyquenHadronizer.cc.
{ return false; }
bool gen::PyquenHadronizer::initializeForExternalPartons | ( | ) |
bool PyquenHadronizer::initializeForInternalPartons | ( | ) |
Definition at line 202 of file PyquenHadronizer.cc.
References comenergy, and pythia6Service_.
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); // Call PYTHIA call_pyinit("CMS", "p", "p", comenergy); return true; }
const char * PyquenHadronizer::nucleon | ( | ) | [private] |
Definition at line 260 of file PyquenHadronizer.cc.
References pfrac_, gen::pyr_(), and random.
Referenced by generatePartonsAndHadronize().
bool PyquenHadronizer::pyqpythia_init | ( | const edm::ParameterSet & | pset | ) | [private] |
Definition at line 215 of file PyquenHadronizer.cc.
References gen::call_pygive().
Referenced by readSettings().
{ //Turn Hadronization Off whether or not there is quenching // PYEXEC is called later anyway string sHadOff("MSTP(111)=0"); gen::call_pygive(sHadOff); return true; }
bool PyquenHadronizer::pyquen_init | ( | const edm::ParameterSet & | pset | ) | [private] |
Definition at line 228 of file PyquenHadronizer.cc.
References angularspecselector_, docollisionalenloss_, doradiativeenloss_, nquarkflavor_, pyqpar, qgpt0_, and qgptau0_.
Referenced by readSettings().
{ // PYQUEN initialization // angular emitted gluon spectrum selection pyqpar.ianglu = angularspecselector_; // type of medium induced partonic energy loss if( doradiativeenloss_ && docollisionalenloss_ ){ edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####"; pyqpar.ienglu = 0; } else if ( doradiativeenloss_ ) { edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####"; pyqpar.ienglu = 1; } else if ( docollisionalenloss_ ) { edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####"; pyqpar.ienglu = 2; } else { edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####"; pyqpar.ienglu = 0; } // 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_; return true; }
bool PyquenHadronizer::readSettings | ( | int | ) |
Definition at line 182 of file PyquenHadronizer.cc.
References abeamtarget_, pfrac_, funct::pow(), pset_, pyqpythia_init(), pyquen_init(), pythia6Service_, gen::Pythia6Service::setCSAParams(), and gen::Pythia6Service::setGeneralParams().
{ Pythia6Service::InstanceWrapper guard(pythia6Service_); pythia6Service_->setGeneralParams(); pythia6Service_->setCSAParams(); //Proton to Nucleon fraction pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3)); //initialize pythia pyqpythia_init(pset_); //initilize pyquen pyquen_init(pset_); return true; }
bool PyquenHadronizer::residualDecay | ( | ) |
Definition at line 335 of file PyquenHadronizer.cc.
{ return true; }
void PyquenHadronizer::rotateEvtPlane | ( | HepMC::GenEvent * | evt, |
double | angle | ||
) | [private] |
Definition at line 270 of file PyquenHadronizer.cc.
References funct::cos(), funct::sin(), lumiQTWidget::t, x, detailsBasic3DVector::y, and z.
Referenced by generatePartonsAndHadronize().
{ double sinphi0 = sin(angle); double cosphi0 = cos(angle); for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin(); vt!=evt->vertices_end(); ++vt ) { double x0 = (*vt)->position().x(); double y0 = (*vt)->position().y(); double z = (*vt)->position().z(); double t = (*vt)->position().t(); double x = x0*cosphi0-y0*sinphi0; double y = y0*cosphi0+x0*sinphi0; (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ; } for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin(); vt!=evt->particles_end(); ++vt ) { double x0 = (*vt)->momentum().x(); double y0 = (*vt)->momentum().y(); double z = (*vt)->momentum().z(); double t = (*vt)->momentum().t(); double x = x0*cosphi0-y0*sinphi0; double y = y0*cosphi0+x0*sinphi0; (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ; } }
virtual bool gen::PyquenHadronizer::select | ( | HepMC::GenEvent * | evtTry | ) | const [inline, virtual] |
Reimplemented from gen::BaseHadronizer.
Definition at line 40 of file PyquenHadronizer.h.
References BaseHiGenEvtSelector::filter(), and selector_.
void PyquenHadronizer::statistics | ( | ) |
Definition at line 344 of file PyquenHadronizer.cc.
{ }
double gen::PyquenHadronizer::abeamtarget_ [private] |
Definition at line 54 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), and readSettings().
unsigned int gen::PyquenHadronizer::angularspecselector_ [private] |
beam/target atomic mass number
Definition at line 55 of file PyquenHadronizer.h.
Referenced by pyquen_init().
double gen::PyquenHadronizer::bfixed_ [private] |
max impact param (fm); valid only if cflag_!=0
Definition at line 61 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
double gen::PyquenHadronizer::bmax_ [private] |
min impact param (fm); valid only if cflag_!=0
Definition at line 60 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
double gen::PyquenHadronizer::bmin_ [private] |
angular emitted gluon spectrum selection DEFAULT= 0 -- small angular emitted gluon spectrum = 1 -- broad angular emitted gluon spectrum = 2 -- collinear angular emitted gluon spectrum
Definition at line 59 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
int gen::PyquenHadronizer::cflag_ [private] |
fixed impact param (fm); valid only if cflag_=0
Definition at line 62 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), and PyquenHadronizer().
double gen::PyquenHadronizer::comenergy [private] |
centrality flag =0 fixed impact param, <>0 minbias
Definition at line 63 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), and initializeForInternalPartons().
bool gen::PyquenHadronizer::docollisionalenloss_ [private] |
bool gen::PyquenHadronizer::doIsospin_ [private] |
DEFAULT = true.
Definition at line 67 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
bool gen::PyquenHadronizer::doquench_ [private] |
collision energy
Definition at line 64 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
bool gen::PyquenHadronizer::doradiativeenloss_ [private] |
if true perform quenching (default = true)
Definition at line 65 of file PyquenHadronizer.h.
Referenced by pyquen_init().
bool gen::PyquenHadronizer::embedding_ [private] |
Definition at line 69 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), and PyquenHadronizer().
double gen::PyquenHadronizer::evtPlane_ [private] |
Definition at line 70 of file PyquenHadronizer.h.
Referenced by add_heavy_ion_rec(), and generatePartonsAndHadronize().
std::string gen::PyquenHadronizer::filterType_ [private] |
Definition at line 86 of file PyquenHadronizer.h.
Referenced by PyquenHadronizer().
unsigned int gen::PyquenHadronizer::maxEventsToPrint_ [private] |
proper time of QGP formation DEFAULT = 0.1 fm/c; allowed range [0.01,10.0]fm/c;
Definition at line 79 of file PyquenHadronizer.h.
Referenced by PyquenHadronizer().
unsigned int gen::PyquenHadronizer::nquarkflavor_ [private] |
Proton fraction in the nucleus.
Definition at line 73 of file PyquenHadronizer.h.
Referenced by pyquen_init().
double gen::PyquenHadronizer::pfrac_ [private] |
Definition at line 71 of file PyquenHadronizer.h.
Referenced by nucleon(), and readSettings().
int gen::PyquenHadronizer::protonSide_ [private] |
Run n&p with proper ratios; if false, only p+p collisions.
Definition at line 68 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize().
Definition at line 53 of file PyquenHadronizer.h.
Referenced by readSettings().
Definition at line 85 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), initializeForInternalPartons(), readSettings(), and ~PyquenHadronizer().
bool gen::PyquenHadronizer::pythiaHepMCVerbosity_ [private] |
Events to print if verbosity.
Definition at line 80 of file PyquenHadronizer.h.
Referenced by PyquenHadronizer().
unsigned int gen::PyquenHadronizer::pythiaPylistVerbosity_ [private] |
HepMC verbosity flag.
Definition at line 81 of file PyquenHadronizer.h.
Referenced by PyquenHadronizer().
double gen::PyquenHadronizer::qgpt0_ [private] |
number of active quark flavors in qgp DEFAULT=0; allowed values: 0,1,2,3.
Definition at line 75 of file PyquenHadronizer.h.
Referenced by pyquen_init().
double gen::PyquenHadronizer::qgptau0_ [private] |
initial temperature of QGP DEFAULT = 1GeV; allowed range [0.2,2.0]GeV;
Definition at line 77 of file PyquenHadronizer.h.
Referenced by pyquen_init().
Definition at line 87 of file PyquenHadronizer.h.
Referenced by PyquenHadronizer(), and select().
edm::InputTag gen::PyquenHadronizer::src_ [private] |
Pythia PYLIST Verbosity flag.
Definition at line 84 of file PyquenHadronizer.h.
Referenced by generatePartonsAndHadronize(), and PyquenHadronizer().