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  * $Id: PyquenHadronizer.cc,v 1.15 2013/01/20 23:55:21 mnguyen Exp $
6 */
7 
8 #include <iostream>
9 #include "time.h"
10 
15 
17 
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 
36  BaseHadronizer(pset),
37  pset_(pset),
38 abeamtarget_(pset.getParameter<double>("aBeamTarget")),
39 angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
40 bmin_(pset.getParameter<double>("bMin")),
41 bmax_(pset.getParameter<double>("bMax")),
42 bfixed_(pset.getParameter<double>("bFixed")),
43 cflag_(pset.getParameter<int>("cFlag")),
44 comenergy(pset.getParameter<double>("comEnergy")),
45 doquench_(pset.getParameter<bool>("doQuench")),
46 doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
47 docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
48 doIsospin_(pset.getParameter<bool>("doIsospin")),
49  protonSide_(pset.getUntrackedParameter<int>("protonSide",0)),
50 embedding_(pset.getParameter<bool>("embeddingMode")),
51 evtPlane_(0),
52 nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
53 qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
54 qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
55 maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
56 pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
57 pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
58 pythia6Service_(new Pythia6Service(pset)),
59 filterType_(pset.getUntrackedParameter<string>("filterType","None"))
60 {
61  // Verbosity Level
62  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
63  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl;
64  // HepMC event verbosity Level
65  pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false);
66  LogDebug("HepMCverbosity") << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl;
67 
68  //Max number of events printed on verbosity level
69  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
70  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl;
71 
72  if(embedding_){
73  cflag_ = 0;
74  src_ = pset.getParameter<InputTag>("backgroundLabel");
75  }
77 
78 }
79 
80 
81 //_____________________________________________________________________
83 {
84  // distructor
85  call_pystat(1);
86 
87  delete pythia6Service_;
88 
89 }
90 
91 
92 //_____________________________________________________________________
93 void PyquenHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
94 {
95  HepMC::HeavyIon *hi = new HepMC::HeavyIon(
96  1, // Ncoll_hard
97  -1, // Npart_proj
98  -1, // Npart_targ
99  1, // Ncoll
100  -1, // spectator_neutrons
101  -1, // spectator_protons
102  -1, // N_Nwounded_collisions
103  -1, // Nwounded_N_collisions
104  -1, // Nwounded_Nwounded_collisions
105  plfpar.bgen, // impact_parameter in [fm]
106  evtPlane_, // event_plane_angle
107  0, // eccentricity
108  -1 // sigma_inel_NN
109  );
110 
111  evt->set_heavy_ion(*hi);
112 
113  delete hi;
114 }
115 
116 //_____________________________________________________________________
118 {
120 
121  // Not possible to retrieve impact paramter and event plane info
122  // at this part, need to overwrite filter() in
123  // PyquenGeneratorFilter
124 
125  const edm::Event& e = getEDMEvent();
126 
127  if(embedding_){
129  e.getByLabel(src_,input);
130  const HepMC::GenEvent * inev = input->GetEvent();
131  const HepMC::HeavyIon* hi = inev->heavy_ion();
132  if(hi){
133  bfixed_ = hi->impact_parameter();
134  evtPlane_ = hi->event_plane_angle();
135  }else{
136  LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
137  }
138  }
139 
140  // Generate PYQUEN event
141  // generate single partonic PYTHIA jet event
142 
143  // Take into account whether it's a nn or pp or pn interaction
144  if(doIsospin_){
145  string projN = "p";
146  string targN = "p";
147  if(protonSide_ != 1) projN = nucleon();
148  if(protonSide_ != 2) targN = nucleon();
149  call_pyinit("CMS", projN.data(), targN.data(), comenergy);
150  }
151  call_pyevnt();
152 
153  // call PYQUEN to apply parton rescattering and energy loss
154  // if doQuench=FALSE, it is pure PYTHIA
155  if( doquench_ ){
157  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
158  } else {
159  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
160  }
161 
162  // call PYTHIA to finish the hadronization
163  pyexec_();
164 
165  // fill the HEPEVT with the PYJETS event record
166  call_pyhepc(1);
167 
168  // event information
169  HepMC::GenEvent* evt = hepevtio.read_next_event();
170 
171  evt->set_signal_process_id(pypars.msti[0]); // type of the process
172  evt->set_event_scale(pypars.pari[16]); // Q^2
173 
175  add_heavy_ion_rec(evt);
176 
177  event().reset(evt);
178 
179  return true;
180 }
181 
183 {
184 
188 
189  //Proton to Nucleon fraction
190  pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
191 
192  //initialize pythia
194 
195  //initilize pyquen
197 
198  return true;
199 
200 }
201 
203 
204 
206 
207  // Call PYTHIA
208  call_pyinit("CMS", "p", "p", comenergy);
209 
210  return true;
211 }
212 
213 
214 //_____________________________________________________________________
216 {
217 
218  //Turn Hadronization Off whether or not there is quenching
219  // PYEXEC is called later anyway
220  string sHadOff("MSTP(111)=0");
221  gen::call_pygive(sHadOff);
222 
223  return true;
224 }
225 
226 
227 //_________________________________________________________________
229 {
230  // PYQUEN initialization
231 
232  // angular emitted gluon spectrum selection
233  pyqpar.ianglu = angularspecselector_;
234 
235  // type of medium induced partonic energy loss
237  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
238  pyqpar.ienglu = 0;
239  } else if ( doradiativeenloss_ ) {
240  edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
241  pyqpar.ienglu = 1;
242  } else if ( docollisionalenloss_ ) {
243  edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
244  pyqpar.ienglu = 2;
245  } else {
246  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
247  pyqpar.ienglu = 0;
248  }
249 
250  // number of active quark flavors in qgp
251  pyqpar.nfu = nquarkflavor_;
252  // initial temperature of QGP
253  pyqpar.T0u = qgpt0_;
254  // proper time of QGP formation
255  pyqpar.tau0u = qgptau0_;
256 
257  return true;
258 }
259 
261  int* dummy = 0;
262  double random = gen::pyr_(dummy);
263  const char* nuc = 0;
264  if(random > pfrac_) nuc = "n";
265  else nuc = "p";
266 
267  return nuc;
268 }
269 
270 void PyquenHadronizer::rotateEvtPlane(HepMC::GenEvent* evt, double angle){
271 
272  double sinphi0 = sin(angle);
273  double cosphi0 = cos(angle);
274 
275  for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
276  vt!=evt->vertices_end(); ++vt )
277  {
278 
279  double x0 = (*vt)->position().x();
280  double y0 = (*vt)->position().y();
281  double z = (*vt)->position().z();
282  double t = (*vt)->position().t();
283 
284  double x = x0*cosphi0-y0*sinphi0;
285  double y = y0*cosphi0+x0*sinphi0;
286 
287  (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;
288  }
289 
290  for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
291  vt!=evt->particles_end(); ++vt )
292  {
293 
294  double x0 = (*vt)->momentum().x();
295  double y0 = (*vt)->momentum().y();
296  double z = (*vt)->momentum().z();
297  double t = (*vt)->momentum().t();
298 
299  double x = x0*cosphi0-y0*sinphi0;
300  double y = y0*cosphi0+x0*sinphi0;
301 
302  (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
303  }
304 }
305 
306 bool PyquenHadronizer::declareStableParticles( std::vector<int> pdg )
307 {
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:25
bool pyquen_init(const edm::ParameterSet &pset)
PyquenHadronizer(const edm::ParameterSet &)
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 double double z
std::auto_ptr< HepMC::GenEvent > & event()
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
double bmax_
min impact param (fm); valid only if cflag_!=0
#define pypars
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.
bool declareStableParticles(const std::vector< int >)
void add_heavy_ion_rec(HepMC::GenEvent *evt)
Pythia6Service * pythia6Service_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
edm::ParameterSet pset_
int pycomp_(int &)
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:18
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
x
Definition: VDTMath.h:216
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