26 #include "HepMC/IO_HEPEVT.h"
27 #include "HepMC/PythiaWrapper.h"
41 abeamtarget_(pset.getParameter<double>(
"aBeamTarget")),
42 angularspecselector_(pset.getParameter<int>(
"angularSpectrumSelector")),
43 bmin_(pset.getParameter<double>(
"bMin")),
44 bmax_(pset.getParameter<double>(
"bMax")),
45 bfixed_(pset.getParameter<double>(
"bFixed")),
46 cflag_(pset.getParameter<int>(
"cFlag")),
47 comenergy(pset.getParameter<double>(
"comEnergy")),
48 doquench_(pset.getParameter<bool>(
"doQuench")),
49 doradiativeenloss_(pset.getParameter<bool>(
"doRadiativeEnLoss")),
50 docollisionalenloss_(pset.getParameter<bool>(
"doCollisionalEnLoss")),
51 doIsospin_(pset.getParameter<bool>(
"doIsospin")),
52 protonSide_(pset.getUntrackedParameter<int>(
"protonSide",0)),
53 embedding_(pset.getParameter<bool>(
"embeddingMode")),
55 nquarkflavor_(pset.getParameter<int>(
"qgpNumQuarkFlavor")),
56 qgpt0_(pset.getParameter<double>(
"qgpInitialTemperature")),
57 qgptau0_(pset.getParameter<double>(
"qgpProperTimeFormation")),
58 maxEventsToPrint_(pset.getUntrackedParameter<int>(
"maxEventsToPrint",1)),
59 pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>(
"pythiaHepMCVerbosity",
false)),
60 pythiaPylistVerbosity_(pset.getUntrackedParameter<int>(
"pythiaPylistVerbosity",0)),
62 filterType_(pset.getUntrackedParameter<
string>(
"filterType",
"None"))
105 HepMC::HeavyIon *hi =
new HepMC::HeavyIon(
121 evt->set_heavy_ion(*hi);
141 const HepMC::HeavyIon* hi = inev->heavy_ion();
143 bfixed_ = hi->impact_parameter();
146 LogWarning(
"EventEmbedding")<<
"Background event does not have heavy ion record!";
159 call_pyinit(
"CMS", projN.data(), targN.data(),
comenergy);
169 edm::LogInfo(
"PYQUENinAction") <<
"##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
181 evt->set_signal_process_id(
pypars.msti[0]);
182 evt->set_event_scale(
pypars.pari[16]);
230 string sHadOff(
"MSTP(111)=0");
247 edm::LogInfo(
"PYQUENinEnLoss") <<
"##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
250 edm::LogInfo(
"PYQUENinRad") <<
"##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
253 edm::LogInfo(
"PYQUENinColl") <<
"##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
256 edm::LogInfo(
"PYQUENinEnLoss") <<
"##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
274 if(random >
pfrac_) nuc =
"n";
282 double sinphi0 =
sin(angle);
283 double cosphi0 =
cos(angle);
286 vt!=evt->vertices_end(); ++vt )
289 double x0 = (*vt)->position().x();
290 double y0 = (*vt)->position().y();
291 double z = (*vt)->position().z();
292 double t = (*vt)->position().t();
294 double x = x0*cosphi0-y0*sinphi0;
295 double y = y0*cosphi0+x0*sinphi0;
297 (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;
300 for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
301 vt!=evt->particles_end(); ++vt )
304 double x0 = (*vt)->momentum().x();
305 double y0 = (*vt)->momentum().y();
306 double z = (*vt)->momentum().z();
307 double t = (*vt)->momentum().t();
309 double x = x0*cosphi0-y0*sinphi0;
310 double y = y0*cosphi0+x0*sinphi0;
312 (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
318 std::vector<int> pdg = _pdg;
319 for (
size_t i=0;
i < pdg.size();
i++ )
322 std::ostringstream pyCard ;
323 pyCard <<
"MDCY(" << pyCode <<
",1)=0";
360 return "gen::PyquenHadronizer";
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned int maxEventsToPrint_
bool pyquen_init(const edm::ParameterSet &pset)
PyquenHadronizer(const edm::ParameterSet &)
bool declareStableParticles(const std::vector< int > &)
Sin< T >::type sin(const T &t)
bool call_pygive(const std::string &line)
bool doradiativeenloss_
if true perform quenching (default = true)
bool docollisionalenloss_
DEFAULT = true.
std::auto_ptr< HepMC::GenEvent > & event()
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
bool initializeForInternalPartons()
static std::string const input
double bmax_
min impact param (fm); valid only if cflag_!=0
def gen
run2 Cosmic #### Run 256259 @ 0T 2015C### Run 272133 @ 3.8T 2016B###
static const std::string kPythia6
T x() const
Cartesian x coordinate.
int cflag_
fixed impact param (fm); valid only if cflag_=0
unsigned int pythiaPylistVerbosity_
HepMC verbosity flag.
bool doIsospin_
DEFAULT = true.
static BaseHiGenEvtSelector * get(std::string, const edm::ParameterSet &)
unsigned int angularspecselector_
beam/target atomic mass number
bool generatePartonsAndHadronize()
Cos< T >::type cos(const T &t)
edm::InputTag src_
Pythia PYLIST Verbosity flag.
void add_heavy_ion_rec(HepMC::GenEvent *evt)
Pythia6Service * pythia6Service_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
static const std::string kFortranInstance
bool pyqpythia_init(const edm::ParameterSet &pset)
unsigned int nquarkflavor_
Proton fraction in the nucleus.
bool doquench_
collision energy
const char * classname() const
double comenergy
centrality flag =0 fixed impact param, <>0 minbias
bool pythiaHepMCVerbosity_
Events to print if verbosity.
virtual ~PyquenHadronizer()
void setRandomEngine(CLHEP::HepRandomEngine *v)
static const std::vector< std::string > theSharedResources
BaseHiGenEvtSelector * selector_
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
double bfixed_
max impact param (fm); valid only if cflag_!=0
volatile std::atomic< bool > shutdown_flag false
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
int protonSide_
Run n&p with proper ratios; if false, only p+p collisions.
Power< A, B >::type pow(const A &a, const B &b)
edm::Event & getEDMEvent() const
HepMC::IO_HEPEVT hepevtio
T angle(T x1, T y1, T z1, T x2, T y2, T z2)