CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
14 
16 
22 
24 
25 #include "HepMC/IO_HEPEVT.h"
26 #include "HepMC/PythiaWrapper.h"
27 
28 using namespace gen;
29 using namespace edm;
30 using namespace std;
31 
32 HepMC::IO_HEPEVT hepevtio;
33 
35  BaseHadronizer(pset),
36  pset_(pset),
37 abeamtarget_(pset.getParameter<double>("aBeamTarget")),
38 angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
39 bmin_(pset.getParameter<double>("bMin")),
40 bmax_(pset.getParameter<double>("bMax")),
41 bfixed_(pset.getParameter<double>("bFixed")),
42 cflag_(pset.getParameter<int>("cFlag")),
43 comenergy(pset.getParameter<double>("comEnergy")),
44 doquench_(pset.getParameter<bool>("doQuench")),
45 doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
46 docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
47 doIsospin_(pset.getParameter<bool>("doIsospin")),
48  protonSide_(pset.getUntrackedParameter<int>("protonSide",0)),
49 embedding_(pset.getParameter<bool>("embeddingMode")),
50 evtPlane_(0),
51 nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
52 qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
53 qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
54 maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
55 pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
56 pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
57 pythia6Service_(new Pythia6Service(pset)),
58 filterType_(pset.getUntrackedParameter<string>("filterType","None"))
59 {
60  // Verbosity Level
61  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
62  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl;
63  // HepMC event verbosity Level
64  pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false);
65  LogDebug("HepMCverbosity") << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl;
66 
67  //Max number of events printed on verbosity level
68  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
69  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl;
70 
71  if(embedding_){
72  cflag_ = 0;
73  src_ = pset.getParameter<InputTag>("backgroundLabel");
74  }
76 
77 }
78 
79 
80 //_____________________________________________________________________
82 {
83  // distructor
84  call_pystat(1);
85 
86  delete pythia6Service_;
87 
88 }
89 
90 
91 //_____________________________________________________________________
92 void PyquenHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
93 {
94  HepMC::HeavyIon *hi = new HepMC::HeavyIon(
95  1, // Ncoll_hard
96  -1, // Npart_proj
97  -1, // Npart_targ
98  1, // Ncoll
99  -1, // spectator_neutrons
100  -1, // spectator_protons
101  -1, // N_Nwounded_collisions
102  -1, // Nwounded_N_collisions
103  -1, // Nwounded_Nwounded_collisions
104  plfpar.bgen, // impact_parameter in [fm]
105  evtPlane_, // event_plane_angle
106  0, // eccentricity
107  -1 // sigma_inel_NN
108  );
109 
110  evt->set_heavy_ion(*hi);
111 
112  delete hi;
113 }
114 
115 //_____________________________________________________________________
117 {
119 
120  // Not possible to retrieve impact paramter and event plane info
121  // at this part, need to overwrite filter() in
122  // PyquenGeneratorFilter
123 
124  const edm::Event& e = getEDMEvent();
125 
126  if(embedding_){
128  e.getByLabel(src_,input);
129  const HepMC::GenEvent * inev = input->GetEvent();
130  const HepMC::HeavyIon* hi = inev->heavy_ion();
131  if(hi){
132  bfixed_ = hi->impact_parameter();
133  evtPlane_ = hi->event_plane_angle();
134  }else{
135  LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
136  }
137  }
138 
139  // Generate PYQUEN event
140  // generate single partonic PYTHIA jet event
141 
142  // Take into account whether it's a nn or pp or pn interaction
143  if(doIsospin_){
144  string projN = "p";
145  string targN = "p";
146  if(protonSide_ != 1) projN = nucleon();
147  if(protonSide_ != 2) targN = nucleon();
148  call_pyinit("CMS", projN.data(), targN.data(), comenergy);
149  }
150  call_pyevnt();
151 
152  // call PYQUEN to apply parton rescattering and energy loss
153  // if doQuench=FALSE, it is pure PYTHIA
154  if( doquench_ ){
156  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
157  } else {
158  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
159  }
160 
161  // call PYTHIA to finish the hadronization
162  pyexec_();
163 
164  // fill the HEPEVT with the PYJETS event record
165  call_pyhepc(1);
166 
167  // event information
168  HepMC::GenEvent* evt = hepevtio.read_next_event();
169 
170  evt->set_signal_process_id(pypars.msti[0]); // type of the process
171  evt->set_event_scale(pypars.pari[16]); // Q^2
172 
174  add_heavy_ion_rec(evt);
175 
176  event().reset(evt);
177 
178  return true;
179 }
180 
182 {
183 
187 
188  //Proton to Nucleon fraction
189  pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
190 
191  //initialize pythia
193 
194  //initilize pyquen
196 
197  return true;
198 
199 }
200 
202 
203 
205 
206  // Call PYTHIA
207  call_pyinit("CMS", "p", "p", comenergy);
208 
209  return true;
210 }
211 
212 
213 //_____________________________________________________________________
215 {
216 
217  //Turn Hadronization Off whether or not there is quenching
218  // PYEXEC is called later anyway
219  string sHadOff("MSTP(111)=0");
220  gen::call_pygive(sHadOff);
221 
222  return true;
223 }
224 
225 
226 //_________________________________________________________________
228 {
229  // PYQUEN initialization
230 
231  // angular emitted gluon spectrum selection
232  pyqpar.ianglu = angularspecselector_;
233 
234  // type of medium induced partonic energy loss
236  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
237  pyqpar.ienglu = 0;
238  } else if ( doradiativeenloss_ ) {
239  edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
240  pyqpar.ienglu = 1;
241  } else if ( docollisionalenloss_ ) {
242  edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
243  pyqpar.ienglu = 2;
244  } else {
245  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
246  pyqpar.ienglu = 0;
247  }
248 
249  // number of active quark flavors in qgp
250  pyqpar.nfu = nquarkflavor_;
251  // initial temperature of QGP
252  pyqpar.T0u = qgpt0_;
253  // proper time of QGP formation
254  pyqpar.tau0u = qgptau0_;
255 
256  return true;
257 }
258 
260  int* dummy = 0;
261  double random = gen::pyr_(dummy);
262  const char* nuc = 0;
263  if(random > pfrac_) nuc = "n";
264  else nuc = "p";
265 
266  return nuc;
267 }
268 
269 void PyquenHadronizer::rotateEvtPlane(HepMC::GenEvent* evt, double angle){
270 
271  double sinphi0 = sin(angle);
272  double cosphi0 = cos(angle);
273 
274  for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
275  vt!=evt->vertices_end(); ++vt )
276  {
277 
278  double x0 = (*vt)->position().x();
279  double y0 = (*vt)->position().y();
280  double z = (*vt)->position().z();
281  double t = (*vt)->position().t();
282 
283  double x = x0*cosphi0-y0*sinphi0;
284  double y = y0*cosphi0+x0*sinphi0;
285 
286  (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;
287  }
288 
289  for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
290  vt!=evt->particles_end(); ++vt )
291  {
292 
293  double x0 = (*vt)->momentum().x();
294  double y0 = (*vt)->momentum().y();
295  double z = (*vt)->momentum().z();
296  double t = (*vt)->momentum().t();
297 
298  double x = x0*cosphi0-y0*sinphi0;
299  double y = y0*cosphi0+x0*sinphi0;
300 
301  (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
302  }
303 }
304 
305 bool PyquenHadronizer::declareStableParticles(const std::vector<int>& _pdg )
306 {
307  std::vector<int> pdg = _pdg;
308  for ( size_t i=0; i < pdg.size(); i++ )
309  {
310  int pyCode = pycomp_( pdg[i] );
311  std::ostringstream pyCard ;
312  pyCard << "MDCY(" << pyCode << ",1)=0";
313  std::cout << pyCard.str() << std::endl;
314  call_pygive( pyCard.str() );
315  }
316 
317  return true;
318 
319 }
320 
321 
322 
323 //____________________________________________________________________
324 
326 {
327  return false;
328 }
329 
331 {
332  return true;
333 }
334 
336 {
337  return true;
338 }
339 
341 
342 }
343 
345 }
346 
347 const char* PyquenHadronizer::classname() const
348 {
349  return "gen::PyquenHadronizer";
350 }
351 
352 
353 
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
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
float float float z
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
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:390
edm::ParameterSet pset_
int pycomp_(int &)
#define pypars
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, &lt;&gt;0 minbias
bool pythiaHepMCVerbosity_
Events to print if verbosity.
#define pyqpar
BaseHiGenEvtSelector * selector_
tuple cout
Definition: gather_cfg.py:121
#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
volatile std::atomic< bool > shutdown_flag false
Definition: DDAxes.h:10
double pyr_(int *idummy)
int protonSide_
Run n&amp;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
HepMC::IO_HEPEVT hepevtio
void pyexec_()
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11