CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AMPTHadronizer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 
4 #include "boost/lexical_cast.hpp"
5 
14 
17 
21 
22 #include "HepMC/GenEvent.h"
23 #include "HepMC/HeavyIon.h"
24 #include "HepMC/SimpleVector.h"
25 #include "CLHEP/Random/RandomEngine.h"
26 
27 static const double pi = 3.14159265358979;
28 
29 using namespace edm;
30 using namespace std;
31 using namespace gen;
32 
33 CLHEP::HepRandomEngine* _amptRandomEngine;
34 
35 extern "C"
36 {
37  float gen::ranart_(int *idummy)
38  {
39  if(0) idummy = idummy;
40  float rannum = _amptRandomEngine->flat();
41  return rannum;
42  }
43 }
44 
45 extern "C"
46 {
47  float gen::ran1_(int *idummy)
48  {
49  if(0) idummy = idummy;
50  return _amptRandomEngine->flat();
51  }
52 }
53 
54 AMPTHadronizer::AMPTHadronizer(const ParameterSet &pset) :
55  BaseHadronizer(pset),
56  evt(0),
57  pset_(pset),
58  bmax_(pset.getParameter<double>("bMax")),
59  bmin_(pset.getParameter<double>("bMin")),
60  efrm_(pset.getParameter<double>("comEnergy")),
61  frame_(pset.getParameter<string>("frame")),
62  proj_(pset.getParameter<string>("proj")),
63  targ_(pset.getParameter<string>("targ")),
64  iap_(pset.getParameter<int>("iap")),
65  izp_(pset.getParameter<int>("izp")),
66  iat_(pset.getParameter<int>("iat")),
67  izt_(pset.getParameter<int>("izt")),
68  amptmode_(pset.getParameter<int>("amptmode")),
69  ntmax_(pset.getParameter<int>("ntmax")),
70  dt_(pset.getParameter<double>("dt")),
71  stringFragA_(pset.getParameter<double>("stringFragA")),
72  stringFragB_(pset.getParameter<double>("stringFragB")),
73  popcornmode_(pset.getParameter<bool>("popcornmode")),
74  popcornpar_(pset.getParameter<double>("popcornpar")),
75  shadowingmode_(pset.getParameter<bool>("shadowingmode")),
76  quenchingmode_(pset.getParameter<bool>("quenchingmode")),
77  quenchingpar_(pset.getParameter<double>("quenchingpar")),
78  pthard_(pset.getParameter<double>("pthard")),
79  mu_(pset.getParameter<double>("mu")),
80  izpc_(pset.getParameter<int>("izpc")),
81  alpha_(pset.getParameter<double>("alpha")),
82  dpcoal_(pset.getParameter<double>("dpcoal")),
83  drcoal_(pset.getParameter<double>("drcoal")),
84  ks0decay_(pset.getParameter<bool>("ks0decay")),
85  phidecay_(pset.getParameter<bool>("phidecay")),
86  deuteronmode_(pset.getParameter<int>("deuteronmode")),
87  deuteronfactor_(pset.getParameter<int>("deuteronfactor")),
88  deuteronxsec_(pset.getParameter<int>("deuteronxsec")),
89  minijetpt_(pset.getParameter<double>("minijetpt")),
90  maxmiss_(pset.getParameter<int>("maxmiss")),
91  doInitialAndFinalRadiation_(pset.getParameter<int>("doInitialAndFinalRadiation")),
92  ktkick_(pset.getParameter<int>("ktkick")),
93  diquarkembedding_(pset.getParameter<int>("diquarkembedding")),
94  diquarkpx_(pset.getParameter<double>("diquarkpx")),
95  diquarkpy_(pset.getParameter<double>("diquarkpy")),
96  diquarkx_(pset.getParameter<double>("diquarkx")),
97  diquarky_(pset.getParameter<double>("diquarky")),
98  phi0_(0.),
99  sinphi0_(0.),
100  cosphi0_(1.),
101  rotate_(pset.getParameter<bool>("rotateEventPlane"))
102 {
103  // Default constructor
105  _amptRandomEngine = &(rng->getEngine());
106 }
107 
108 
109 //_____________________________________________________________________
111 {
112 }
113 
114 //_____________________________________________________________________
115 void AMPTHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
116 {
117  // heavy ion record in the final CMSSW Event
118  HepMC::HeavyIon* hi = new HepMC::HeavyIon(
119  hmain1.jatt, // Ncoll_hard/N of SubEvents
120  hmain1.np, // Npart_proj
121  hmain1.nt, // Npart_targ
122  hmain1.n0+hmain1.n01+hmain1.n10+hmain1.n11, // Ncoll
123  0, // spectator_neutrons
124  0, // spectator_protons
125  hmain1.n01, // N_Nwounded_collisions
126  hmain1.n10, // Nwounded_N_collisions
127  hmain1.n11, // Nwounded_Nwounded_collisions
128  hparnt.hint1[18], // impact_parameter in [fm]
129  phi0_, // event_plane_angle
130  0, // eccentricity
131  hparnt.hint1[11] // sigma_inel_NN
132  );
133  evt->set_heavy_ion(*hi);
134  delete hi;
135 }
136 
137 //___________________________________________________________________
139 {
140  // Build particle object corresponding to index in ampt
141 
142  float px0 = hbt.plast[index][0];
143  float py0 = hbt.plast[index][1];
144  float pz0 = hbt.plast[index][2];
145  float m = hbt.plast[index][3];
146 
147  float px = px0*cosphi0_-py0*sinphi0_;
148  float py = py0*cosphi0_+px0*sinphi0_;
149  float pz = pz0;
150  float e = sqrt(px*px+py*py+pz*pz+m*m);
151  int status = 1;
152 
154  HepMC::FourVector(px,
155  py,
156  pz,
157  e),
158  INVFLV(hbt.lblast[index]),// id
159  status // status
160  );
161 
162  p->suggest_barcode(barcode);
163  return p;
164 }
165 
166 //___________________________________________________________________
167 HepMC::GenVertex* AMPTHadronizer::build_ampt_vertex(int i,int id)
168 {
169  // build verteces for the ampt stored events
170  HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(0,0,0,0),id);
171  return vertex;
172 }
173 //_____________________________________________________________________
175 {
176  // generate single event
177  if(rotate_) rotateEvtPlane();
178 
179  // generate a AMPT event
180  AMPT(frame_.data(), bmin_, bmax_, strlen(frame_.data()));
181 
182  // event information
183  HepMC::GenEvent *evt = new HepMC::GenEvent();
184  get_particles(evt);
185 
186  add_heavy_ion_rec(evt);
187 
188  event().reset(evt);
189 
190 
191  return true;
192 }
193 
194 //_____________________________________________________________________
195 bool AMPTHadronizer::get_particles(HepMC::GenEvent *evt )
196 {
197  HepMC::GenVertex* vertice;
198 
199  vector<HepMC::GenParticle*> particles;
200  vector<int> mother_ids;
201  vector<HepMC::GenVertex*> prods;
202 
203  vertice = new HepMC::GenVertex(HepMC::FourVector(0,0,0,0),0);
204  evt->add_vertex(vertice);
205  if(!evt->signal_process_vertex()) evt->set_signal_process_vertex(vertice);
206 
207  const unsigned int knumpart = hbt.nlast;
208  for (unsigned int ipart = 0; ipart<knumpart; ipart++) {
209  int mid = 0;
210  particles.push_back(build_ampt(ipart,ipart+1));
211  prods.push_back(build_ampt_vertex(ipart,0));
212  mother_ids.push_back(mid);
213  LogDebug("DecayChain")<<"Mother index : "<<mid;
214  }
215 
216  LogDebug("AMPT")<<"Number of particles in vector "<<particles.size();
217 
218  for (unsigned int ipart = 0; ipart<particles.size(); ipart++) {
219  HepMC::GenParticle* part = particles[ipart];
220 
221  int mid = mother_ids[ipart];
222  LogDebug("DecayChain")<<"Particle "<<ipart;
223  LogDebug("DecayChain")<<"Mother's ID "<<mid;
224  LogDebug("DecayChain")<<"Particle's PDG ID "<<part->pdg_id();
225 
226  if(mid <= 0){
227  vertice->add_particle_out(part);
228  continue;
229  }
230 
231  if(mid > 0){
232  HepMC::GenParticle* mother = particles[mid];
233  LogDebug("DecayChain")<<"Mother's PDG ID "<<mother->pdg_id();
234 
235  HepMC::GenVertex* prod_vertex = mother->end_vertex();
236  if(!prod_vertex){
237  prod_vertex = prods[ipart];
238  prod_vertex->add_particle_in(mother);
239  evt->add_vertex(prod_vertex);
240  prods[ipart]=0; // mark to protect deletion
241 
242  }
243  prod_vertex->add_particle_out(part);
244  }
245  }
246 
247  // cleanup vertices not assigned to evt
248  for (unsigned int i = 0; i<prods.size(); i++) {
249  if(prods[i]) delete prods[i];
250  }
251 
252  return true;
253 }
254 
255 //_____________________________________________________________________
256 bool AMPTHadronizer::call_amptset(double efrm, std::string frame, std::string proj, std::string targ, int iap, int izp, int iat, int izt)
257 {
258  // initialize hydjet
259  AMPTSET(efrm,frame.data(),proj.data(),targ.data(),iap,izp,iat,izt,strlen(frame.data()),strlen(proj.data()),strlen(targ.data()));
260  return true;
261 }
262 //______________________________________________________________________
264 {
265  anim.isoft=amptmode_;
266  input2.ntmax=ntmax_;
267  input1.dt=dt_;
268  ludat1.parj[40]=stringFragA_;
269  ludat1.parj[41]=stringFragB_;
270  popcorn.ipop=popcornmode_;
271  ludat1.parj[4]=popcornpar_;
272  hparnt.ihpr2[5]=shadowingmode_;
273  hparnt.ihpr2[3]=quenchingmode_;
274  hparnt.hipr1[13]=quenchingpar_;
275  hparnt.hipr1[7]=pthard_;
276  para2.xmu=mu_;
277  anim.izpc=izpc_;
278  para2.alpha=alpha_;
279  coal.dpcoal=dpcoal_;
280  coal.drcoal=drcoal_;
281  resdcy.iksdcy=ks0decay_;
282  phidcy.iphidcy=phidecay_;
283  para8.idpert=deuteronmode_;
284  para8.npertd=deuteronfactor_;
285  para8.idxsec=deuteronxsec_;
286  phidcy.pttrig=minijetpt_;
287  phidcy.maxmiss=maxmiss_;
289  hparnt.ihpr2[4]=ktkick_;
290  embed.iembed=diquarkembedding_;
291  embed.pxqembd=diquarkpx_;
292  embed.pyqembd=diquarkpy_;
293  embed.xembd=diquarkx_;
294  embed.yembd=diquarky_;
295 
296  return true;
297 }
298 
299 //_____________________________________________________________________
301 
302  // ampt running options
303  ampt_init(pset_);
304 
305  // initialize ampt
306  LogInfo("AMPTinAction") << "##### Calling AMPTSET(" << efrm_ << "," <<frame_<<","<<proj_<<","<<targ_<<","<<iap_<<","<<izp_<<","<<iat_<<","<<izt_<<") ####";
307 
309 
310  return true;
311 }
312 
313 bool AMPTHadronizer::declareStableParticles( const std::vector<int>& pdg )
314 {
315  return true;
316 }
317 
318 //________________________________________________________________
320  int zero = 0;
321  double test = (double)gen::ranart_(&zero);
322  phi0_ = 2.*pi*test - pi;
323  sinphi0_ = sin(phi0_);
324  cosphi0_ = cos(phi0_);
325 }
326 
327 //________________________________________________________________
329 {
330  return false;
331 }
332 
334 {
335  return true;
336 }
337 
339 {
340  return true;
341 }
342 
344  return;
345 }
346 
348  return;
349 }
350 
351 const char* AMPTHadronizer::classname() const
352 {
353  return "gen::AMPTHadronizer";
354 }
355 
#define LogDebug(id)
#define popcorn
Definition: AMPTWrapper.h:185
float ran1_(int *)
int i
Definition: DBlmapReader.cc:9
#define ludat1
Definition: AMPTWrapper.h:70
bool initializeForInternalPartons()
#define phidcy
Definition: AMPTWrapper.h:167
#define AMPT
Definition: AMPTWrapper.h:19
HepMC::GenEvent * evt
CLHEP::HepRandomEngine * _amptRandomEngine
void add_heavy_ion_rec(HepMC::GenEvent *evt)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool generatePartonsAndHadronize()
#define hparnt
Definition: AMPTWrapper.h:60
#define anim
Definition: AMPTWrapper.h:80
#define input2
Definition: AMPTWrapper.h:149
std::auto_ptr< HepMC::GenEvent > & event()
bool declareStableParticles(const std::vector< int > &)
#define embed
Definition: AMPTWrapper.h:178
const char * classname() const
#define AMPTSET
Definition: AMPTWrapper.h:14
#define resdcy
Definition: AMPTWrapper.h:157
bool call_amptset(double efrm, std::string frame, std::string proj, std::string targ, int iap, int izp, int iat, int izt)
#define hbt
Definition: AMPTWrapper.h:50
double p[5][pyjets_maxn]
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define hmain1
Definition: AMPTWrapper.h:40
#define para2
Definition: AMPTWrapper.h:100
#define input1
Definition: AMPTWrapper.h:129
virtual CLHEP::HepRandomEngine & getEngine() const =0
Use this to get the random number engine, this is the only function most users should call...
edm::ParameterSet pset_
float ranart_(int *)
part
Definition: HCALResponse.h:20
#define para8
Definition: AMPTWrapper.h:118
HepMC::GenParticle * build_ampt(int index, int barcode)
double pi
tuple status
Definition: ntuplemaker.py:245
bool ampt_init(const edm::ParameterSet &pset)
#define INVFLV
Definition: AMPTWrapper.h:24
HepMC::GenVertex * build_ampt_vertex(int i, int id)
#define coal
Definition: AMPTWrapper.h:89
bool get_particles(HepMC::GenEvent *evt)