CMS 3D CMS Logo

PyquenHadronizer.cc
Go to the documentation of this file.
1 /*
2  *
3  * Generates PYQUEN HepMC events
4  *
5 */
6 
7 #include <iostream>
8 #include <ctime>
9 
11 
16 
18 
23 
25 
26 #include "HepMC/IO_HEPEVT.h"
27 #include "HepMC/PythiaWrapper.h"
28 
29 using namespace gen;
30 using namespace edm;
31 using namespace std;
32 
33 HepMC::IO_HEPEVT hepevtio;
34 
37 
39  BaseHadronizer(pset),
40  pset_(pset),
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")),
54 evtPlane_(0),
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)),
61 pythia6Service_(new Pythia6Service(pset)),
62 filterType_(pset.getUntrackedParameter<string>("filterType","None"))
63 {
64  // Verbosity Level
65  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
66  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl;
67  // HepMC event verbosity Level
68  pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false);
69  LogDebug("HepMCverbosity") << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl;
70 
71  //Max number of events printed on verbosity level
72  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
73  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl;
74 
75  if(embedding_){
76  cflag_ = 0;
77  src_ = pset.getParameter<InputTag>("backgroundLabel");
78  }
80 
81 }
82 
83 
84 //_____________________________________________________________________
86 {
87  // distructor
88  call_pystat(1);
89 
90  delete pythia6Service_;
91 
92 }
93 
94 
95 //_____________________________________________________________________
96 void PyquenHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v)
97 {
99 }
100 
101 
102 //_____________________________________________________________________
104 {
105  HepMC::HeavyIon *hi = new HepMC::HeavyIon(
106  1, // Ncoll_hard
107  -1, // Npart_proj
108  -1, // Npart_targ
109  1, // Ncoll
110  -1, // spectator_neutrons
111  -1, // spectator_protons
112  -1, // N_Nwounded_collisions
113  -1, // Nwounded_N_collisions
114  -1, // Nwounded_Nwounded_collisions
115  plfpar.bgen, // impact_parameter in [fm]
116  evtPlane_, // event_plane_angle
117  0, // eccentricity
118  -1 // sigma_inel_NN
119  );
120 
121  evt->set_heavy_ion(*hi);
122 
123  delete hi;
124 }
125 
126 //_____________________________________________________________________
128 {
130 
131  // Not possible to retrieve impact paramter and event plane info
132  // at this part, need to overwrite filter() in
133  // PyquenGeneratorFilter
134 
135  const edm::Event& e = getEDMEvent();
136 
137  if(embedding_){
139  e.getByLabel(src_,input);
140  const HepMC::GenEvent * inev = input->GetEvent();
141  const HepMC::HeavyIon* hi = inev->heavy_ion();
142  if(hi){
143  bfixed_ = hi->impact_parameter();
144  evtPlane_ = hi->event_plane_angle();
145  }else{
146  LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
147  }
148  }
149 
150  // Generate PYQUEN event
151  // generate single partonic PYTHIA jet event
152 
153  // Take into account whether it's a nn or pp or pn interaction
154  if(doIsospin_){
155  string projN = "p";
156  string targN = "p";
157  if(protonSide_ != 1) projN = nucleon();
158  if(protonSide_ != 2) targN = nucleon();
159  call_pyinit("CMS", projN.data(), targN.data(), comenergy);
160  }
161  call_pyevnt();
162 
163  // call PYQUEN to apply parton rescattering and energy loss
164  // if doQuench=FALSE, it is pure PYTHIA
165  if( doquench_ ){
167  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
168  } else {
169  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
170  }
171 
172  // call PYTHIA to finish the hadronization
173  pyexec_();
174 
175  // fill the HEPEVT with the PYJETS event record
176  call_pyhepc(1);
177 
178  // event information
179  HepMC::GenEvent* evt = hepevtio.read_next_event();
180 
181  evt->set_signal_process_id(pypars.msti[0]); // type of the process
182  evt->set_event_scale(pypars.pari[16]); // Q^2
183 
185  add_heavy_ion_rec(evt);
186 
187  event().reset(evt);
188 
189  return true;
190 }
191 
193 {
194 
198 
199  //Proton to Nucleon fraction
200  pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
201 
202  //initialize pythia
204 
205  //initilize pyquen
207 
208  return true;
209 
210 }
211 
213 
214 
216 
217  // Call PYTHIA
218  call_pyinit("CMS", "p", "p", comenergy);
219 
220  return true;
221 }
222 
223 
224 //_____________________________________________________________________
226 {
227 
228  //Turn Hadronization Off whether or not there is quenching
229  // PYEXEC is called later anyway
230  string sHadOff("MSTP(111)=0");
231  gen::call_pygive(sHadOff);
232 
233  return true;
234 }
235 
236 
237 //_________________________________________________________________
239 {
240  // PYQUEN initialization
241 
242  // angular emitted gluon spectrum selection
243  pyqpar.ianglu = angularspecselector_;
244 
245  // type of medium induced partonic energy loss
247  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
248  pyqpar.ienglu = 0;
249  } else if ( doradiativeenloss_ ) {
250  edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
251  pyqpar.ienglu = 1;
252  } else if ( docollisionalenloss_ ) {
253  edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
254  pyqpar.ienglu = 2;
255  } else {
256  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
257  pyqpar.ienglu = 0;
258  }
259 
260  // number of active quark flavors in qgp
261  pyqpar.nfu = nquarkflavor_;
262  // initial temperature of QGP
263  pyqpar.T0u = qgpt0_;
264  // proper time of QGP formation
265  pyqpar.tau0u = qgptau0_;
266 
267  return true;
268 }
269 
271  int* dummy = 0;
272  double random = gen::pyr_(dummy);
273  const char* nuc = 0;
274  if(random > pfrac_) nuc = "n";
275  else nuc = "p";
276 
277  return nuc;
278 }
279 
281 
282  double sinphi0 = sin(angle);
283  double cosphi0 = cos(angle);
284 
285  for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
286  vt!=evt->vertices_end(); ++vt )
287  {
288 
289  double x0 = (*vt)->position().x();
290  double y0 = (*vt)->position().y();
291  double z = (*vt)->position().z();
292  double t = (*vt)->position().t();
293 
294  double x = x0*cosphi0-y0*sinphi0;
295  double y = y0*cosphi0+x0*sinphi0;
296 
297  (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;
298  }
299 
300  for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
301  vt!=evt->particles_end(); ++vt )
302  {
303 
304  double x0 = (*vt)->momentum().x();
305  double y0 = (*vt)->momentum().y();
306  double z = (*vt)->momentum().z();
307  double t = (*vt)->momentum().t();
308 
309  double x = x0*cosphi0-y0*sinphi0;
310  double y = y0*cosphi0+x0*sinphi0;
311 
312  (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
313  }
314 }
315 
316 bool PyquenHadronizer::declareStableParticles(const std::vector<int>& _pdg )
317 {
318  std::vector<int> pdg = _pdg;
319  for ( size_t i=0; i < pdg.size(); i++ )
320  {
321  int pyCode = pycomp_( pdg[i] );
322  std::ostringstream pyCard ;
323  pyCard << "MDCY(" << pyCode << ",1)=0";
324  std::cout << pyCard.str() << std::endl;
325  call_pygive( pyCard.str() );
326  }
327 
328  return true;
329 
330 }
331 
332 
333 
334 //____________________________________________________________________
335 
337 {
338  return false;
339 }
340 
342 {
343  return true;
344 }
345 
347 {
348  return true;
349 }
350 
352 
353 }
354 
356 }
357 
358 const char* PyquenHadronizer::classname() const
359 {
360  return "gen::PyquenHadronizer";
361 }
362 
363 
364 
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned int maxEventsToPrint_
#define plfpar
Definition: PyquenWrapper.h:24
bool pyquen_init(const edm::ParameterSet &pset)
PyquenHadronizer(const edm::ParameterSet &)
bool declareStableParticles(const std::vector< int > &)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool call_pygive(const std::string &line)
bool doradiativeenloss_
if true perform quenching (default = true)
bool docollisionalenloss_
DEFAULT = true.
TRandom random
Definition: MVATrainer.cc:138
double v[5][pyjets_maxn]
std::auto_ptr< HepMC::GenEvent > & event()
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
static std::string const input
Definition: EdmProvDump.cc:44
double bmax_
min impact param (fm); valid only if cflag_!=0
static const std::string kPythia6
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
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
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
Definition: Event.h:464
HepMC::IO_HEPEVT hepevtio
edm::ParameterSet pset_
static const std::string kFortranInstance
int pycomp_(int &)
#define pypars
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
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.
void setRandomEngine(CLHEP::HepRandomEngine *v)
#define pyqpar
HLT enums.
static const std::vector< std::string > theSharedResources
BaseHiGenEvtSelector * selector_
#define PYQUEN
Definition: PyquenWrapper.h:17
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19
double bfixed_
max impact param (fm); valid only if cflag_!=0
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
double pyr_(int *idummy)
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)
Definition: Power.h:40
edm::Event & getEDMEvent() const
void pyexec_()
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11