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  int cm=1, va, vb, vc;
82  PYQVER(cm,va,vb,vc);
83 
84 }
85 
86 
87 //_____________________________________________________________________
89 {
90  // distructor
91  call_pystat(1);
92 
93  delete pythia6Service_;
94 
95 }
96 
97 
98 //_____________________________________________________________________
99 void PyquenHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v)
100 {
102 }
103 
104 
105 //_____________________________________________________________________
107 {
108  HepMC::HeavyIon *hi = new HepMC::HeavyIon(
109  1, // Ncoll_hard
110  -1, // Npart_proj
111  -1, // Npart_targ
112  1, // Ncoll
113  -1, // spectator_neutrons
114  -1, // spectator_protons
115  -1, // N_Nwounded_collisions
116  -1, // Nwounded_N_collisions
117  -1, // Nwounded_Nwounded_collisions
118  plfpar.bgen, // impact_parameter in [fm]
119  evtPlane_, // event_plane_angle
120  0, // eccentricity
121  -1 // sigma_inel_NN
122  );
123 
124  evt->set_heavy_ion(*hi);
125 
126  delete hi;
127 }
128 
129 //_____________________________________________________________________
131 {
133 
134  // Not possible to retrieve impact paramter and event plane info
135  // at this part, need to overwrite filter() in
136  // PyquenGeneratorFilter
137 
138  const edm::Event& e = getEDMEvent();
139 
140  if(embedding_){
142  e.getByLabel(src_,input);
143  const HepMC::GenEvent * inev = input->GetEvent();
144  const HepMC::HeavyIon* hi = inev->heavy_ion();
145  if(hi){
146  bfixed_ = hi->impact_parameter();
147  evtPlane_ = hi->event_plane_angle();
148  }else{
149  LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
150  }
151  }
152 
153  // Generate PYQUEN event
154  // generate single partonic PYTHIA jet event
155 
156  // Take into account whether it's a nn or pp or pn interaction
157  if(doIsospin_){
158  string projN = "p";
159  string targN = "p";
160  if(protonSide_ != 1) projN = nucleon();
161  if(protonSide_ != 2) targN = nucleon();
162  call_pyinit("CMS", projN.data(), targN.data(), comenergy);
163  }
164  call_pyevnt();
165 
166  // call PYQUEN to apply parton rescattering and energy loss
167  // if doQuench=FALSE, it is pure PYTHIA
168  if( doquench_ ){
170  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
171  } else {
172  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
173  }
174 
175  // call PYTHIA to finish the hadronization
176  pyexec_();
177 
178  // fill the HEPEVT with the PYJETS event record
179  call_pyhepc(1);
180 
181  // event information
182  HepMC::GenEvent* evt = hepevtio.read_next_event();
183 
184  evt->set_signal_process_id(pypars.msti[0]); // type of the process
185  evt->set_event_scale(pypars.pari[16]); // Q^2
186 
188  add_heavy_ion_rec(evt);
189 
190  HepMC::HEPEVT_Wrapper::check_hepevt_consistency();
191  //std::cout<<"Entries number: "<<HepMC::HEPEVT_Wrapper::number_entries() <<" Max. entries "<<HepMC::HEPEVT_Wrapper::max_number_entries()<<std::endl;
192 
193  event().reset(evt);
194 
195  return true;
196 }
197 
199 {
200 
204 
205  //Proton to Nucleon fraction
206  pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
207 
208  //initialize pythia
210 
211  //initilize pyquen
213 
214  return true;
215 
216 }
217 
219 
220 
222 
223  // Call PYTHIA
224  call_pyinit("CMS", "p", "p", comenergy);
225 
226  return true;
227 }
228 
229 
230 //_____________________________________________________________________
232 {
233 
234  //Turn Hadronization Off whether or not there is quenching
235  // PYEXEC is called later anyway
236  string sHadOff("MSTP(111)=0");
237  gen::call_pygive(sHadOff);
238 
239  return true;
240 }
241 
242 
243 //_________________________________________________________________
245 {
246  // PYQUEN initialization
247 
248  // angular emitted gluon spectrum selection
249  pyqpar.ianglu = angularspecselector_;
250 
251  // type of medium induced partonic energy loss
253  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
254  pyqpar.ienglu = 0;
255  } else if ( doradiativeenloss_ ) {
256  edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
257  pyqpar.ienglu = 1;
258  } else if ( docollisionalenloss_ ) {
259  edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
260  pyqpar.ienglu = 2;
261  } else {
262  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
263  pyqpar.ienglu = 0;
264  }
265 
266  // number of active quark flavors in qgp
267  pyqpar.nfu = nquarkflavor_;
268  // initial temperature of QGP
269  pyqpar.T0u = qgpt0_;
270  // proper time of QGP formation
271  pyqpar.tau0u = qgptau0_;
272 
273  return true;
274 }
275 
277  int* dummy = nullptr;
278  double random = gen::pyr_(dummy);
279  const char* nuc = nullptr;
280  if(random > pfrac_) nuc = "n";
281  else nuc = "p";
282 
283  return nuc;
284 }
285 
287 
288  double sinphi0 = sin(angle);
289  double cosphi0 = cos(angle);
290 
291  for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
292  vt!=evt->vertices_end(); ++vt )
293  {
294 
295  double x0 = (*vt)->position().x();
296  double y0 = (*vt)->position().y();
297  double z = (*vt)->position().z();
298  double t = (*vt)->position().t();
299 
300  double x = x0*cosphi0-y0*sinphi0;
301  double y = y0*cosphi0+x0*sinphi0;
302 
303  (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;
304  }
305 
306  for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
307  vt!=evt->particles_end(); ++vt )
308  {
309 
310  double x0 = (*vt)->momentum().x();
311  double y0 = (*vt)->momentum().y();
312  double z = (*vt)->momentum().z();
313  double t = (*vt)->momentum().t();
314 
315  double x = x0*cosphi0-y0*sinphi0;
316  double y = y0*cosphi0+x0*sinphi0;
317 
318  (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
319  }
320 }
321 
322 bool PyquenHadronizer::declareStableParticles(const std::vector<int>& _pdg )
323 {
324  std::vector<int> pdg = _pdg;
325  for ( size_t i=0; i < pdg.size(); i++ )
326  {
327  int pyCode = pycomp_( pdg[i] );
328  std::ostringstream pyCard ;
329  pyCard << "MDCY(" << pyCode << ",1)=0";
330  std::cout << pyCard.str() << std::endl;
331  call_pygive( pyCard.str() );
332  }
333 
334  return true;
335 
336 }
337 
338 
339 
340 //____________________________________________________________________
341 
343 {
344  return false;
345 }
346 
348 {
349  return true;
350 }
351 
353 {
354  return true;
355 }
356 
358 
359 }
360 
362 }
363 
364 const char* PyquenHadronizer::classname() const
365 {
366  return "gen::PyquenHadronizer";
367 }
368 
369 
370 
#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:29
bool pyquen_init(const edm::ParameterSet &pset)
PyquenHadronizer(const edm::ParameterSet &)
bool declareStableParticles(const std::vector< int > &)
double pyr_(int *idummy)
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]
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
#define PYQVER
Definition: PyquenWrapper.h:22
static std::string const input
Definition: EdmProvDump.cc:48
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_
std::unique_ptr< HepMC::GenEvent > & event()
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
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
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)
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