CMS 3D CMS Logo

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