CMS 3D CMS Logo

HijingHadronizer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 
4 #include "boost/lexical_cast.hpp"
5 
12 
19 
20 #include "HepMC/GenEvent.h"
21 #include "HepMC/HeavyIon.h"
22 #include "HepMC/SimpleVector.h"
23 
24 #include "CLHEP/Random/RandomEngine.h"
25 
26 static const double pi = 3.14159265358979;
27 
28 using namespace edm;
29 using namespace std;
30 using namespace gen;
31 
32 static CLHEP::HepRandomEngine* hijRandomEngine;
33 
34 extern "C"
35 {
36  float gen::hijran_(int *idummy)
37  {
38  return hijRandomEngine->flat();
39  }
40 }
41 
42 extern "C" {
43  float ran_(unsigned int* iseed){
44  return hijRandomEngine->flat();
45  // return ranff_(iseed);
46  // return gen::pyr_(0);
47  }
48 }
49 
50 extern "C" {
51  float rlu_(unsigned int* iseed){
52  return hijRandomEngine->flat();
53  // return ranff_(iseed);
54  // return gen::pyr_(0);
55  }
56 }
57 
58 const std::vector<std::string> HijingHadronizer::theSharedResources = { edm::SharedResourceNames::kPythia6 };
59 
60 HijingHadronizer::HijingHadronizer(const ParameterSet &pset) :
61  BaseHadronizer(pset),
62  evt(0),
63  pset_(pset),
64  bmax_(pset.getParameter<double>("bMax")),
65  bmin_(pset.getParameter<double>("bMin")),
66  efrm_(pset.getParameter<double>("comEnergy")),
67  frame_(pset.getParameter<string>("frame")),
68  proj_(pset.getParameter<string>("proj")),
69  targ_(pset.getParameter<string>("targ")),
70  iap_(pset.getParameter<int>("iap")),
71  izp_(pset.getParameter<int>("izp")),
72  iat_(pset.getParameter<int>("iat")),
73  izt_(pset.getParameter<int>("izt")),
74  phi0_(0.),
75  sinphi0_(0.),
76  cosphi0_(1.),
77  rotate_(pset.getParameter<bool>("rotateEventPlane"))
78 {
79  // Default constructor
80 }
81 
82 
83 //_____________________________________________________________________
85 {
86  // destructor
87 }
88 
89 
90 //_____________________________________________________________________
91 void HijingHadronizer::doSetRandomEngine(CLHEP::HepRandomEngine* v)
92 {
94 }
95 
96 
97 //_____________________________________________________________________
99 {
100  // heavy ion record in the final CMSSW Event
101  HepMC::HeavyIon* hi = new HepMC::HeavyIon(
102  himain1.jatt, // Ncoll_hard/N of SubEvents
103  himain1.np, // Npart_proj
104  himain1.nt, // Npart_targ
105  himain1.n0+himain1.n01+himain1.n10+himain1.n11, // Ncoll
106  0, // spectator_neutrons
107  0, // spectator_protons
108  himain1.n01, // N_Nwounded_collisions
109  himain1.n10, // Nwounded_N_collisions
110  himain1.n11, // Nwounded_Nwounded_collisions
111  //gsfs Changed from 19 to 18 (Fortran counts from 1 , not 0)
112  hiparnt.hint1[18], // impact_parameter in [fm]
113  phi0_, // event_plane_angle
114  0, // eccentricity
115  //gsfs Changed from 12 to 11 (Fortran counts from 1 , not 0)
116  hiparnt.hint1[11] // sigma_inel_NN
117  );
118  evt->set_heavy_ion(*hi);
119  delete hi;
120 }
121 
122 //___________________________________________________________________
124 {
125  // Build particle object corresponding to index in hijing
126 
127  double x0 = himain2.patt[0][index];
128  double y0 = himain2.patt[1][index];
129 
130  double x = x0*cosphi0_-y0*sinphi0_;
131  double y = y0*cosphi0_+x0*sinphi0_;
132 
133  // Hijing gives V0's status=4, they need to have status=1 to be decayed in geant
134  // also change status=11 to status=2
135  if(himain2.katt[3][index]<=10 && himain2.katt[3][index]>0) himain2.katt[3][index]=1;
136  if(himain2.katt[3][index]<=20 && himain2.katt[3][index]>10) himain2.katt[3][index]=2;
137 
139  HepMC::FourVector(x, // px
140  y, // py
141  himain2.patt[2][index], // pz
142  himain2.patt[3][index]), // E
143  himain2.katt[0][index],// id
144  himain2.katt[3][index] // status
145  );
146  p->suggest_barcode(barcode);
147 
148  return p;
149 }
150 
151 //___________________________________________________________________
152 HepMC::GenVertex* HijingHadronizer::build_hijing_vertex(int i,int id)
153 {
154  // build verteces for the hijing stored events
155  double x0=himain2.vatt[0][i];
156  double y0=himain2.vatt[1][i];
157  double x = x0*cosphi0_-y0*sinphi0_;
158  double y = y0*cosphi0_+x0*sinphi0_;
159  double z=himain2.vatt[2][i];
160  double t=himain2.vatt[3][i];
161 
162  HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t),id);
163  return vertex;
164 }
165 
167 {
168  // generate single event
169  if(rotate_) rotateEvtPlane();
170 
171  // generate a HIJING event
172 
173  float f_bmin = bmin_;
174  float f_bmax = bmax_;
175  HIJING(frame_.data(), f_bmin, f_bmax, strlen(frame_.data()));
176 
177  // event information
178  HepMC::GenEvent *evt = new HepMC::GenEvent();
179  get_particles(evt);
180 
181  // evt->set_signal_process_id(pypars.msti[0]); // type of the process
182  // evt->set_event_scale(pypars.pari[16]); // Q^2
183  add_heavy_ion_rec(evt);
184 
185  event().reset(evt);
186 
187  return true;
188 }
189 
190 //_____________________________________________________________________
191 bool HijingHadronizer::get_particles(HepMC::GenEvent *evt )
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 = himain1.natt;
204 
205  for (unsigned int ipart = 0; ipart<knumpart; ipart++) {
206 
207  int mid = himain2.katt[2][ipart] - 1; // careful of fortan to c++ array index
208 
209  particles.push_back(build_hijing(ipart,ipart+1));
210  prods.push_back(build_hijing_vertex(ipart,0));
211  mother_ids.push_back(mid);
212  LogDebug("DecayChain")<<"Mother index : "<<mid;
213  }
214 
215  LogDebug("Hijing")<<"Number of particles in vector "<<particles.size();
216 
217  for (unsigned int ipart = 0; ipart<particles.size(); ipart++) {
218  HepMC::GenParticle* part = particles[ipart];
219 
220  int mid = mother_ids[ipart];
221  LogDebug("DecayChain")<<"Particle "<<ipart;
222  LogDebug("DecayChain")<<"Mother's ID "<<mid;
223  LogDebug("DecayChain")<<"Particle's PDG ID "<<part->pdg_id();
224 
225  // remove zero pT particles from list, protection for fastJet against pt=0 jets
226  if(part->status()==1&&sqrt(part->momentum().px()*part->momentum().px()+part->momentum().py()*part->momentum().py())==0)
227  continue;
228 
229  if(mid <= 0){
230  vertice->add_particle_out(part);
231  continue;
232  }
233 
234  if(mid > 0){
235  HepMC::GenParticle* mother = particles[mid];
236  LogDebug("DecayChain")<<"Mother's PDG ID "<<mother->pdg_id();
237  HepMC::GenVertex* prod_vertex = mother->end_vertex();
238  if(!prod_vertex){
239  prod_vertex = prods[ipart];
240  prod_vertex->add_particle_in(mother);
241 
242  evt->add_vertex(prod_vertex);
243  prods[ipart]=0; // mark to protect deletion
244 
245  }
246  prod_vertex->add_particle_out(part);
247  }
248  }
249 
250  // cleanup vertices not assigned to evt
251  for (unsigned int i = 0; i<prods.size(); i++) {
252  if(prods[i]) delete prods[i];
253  }
254 
255  return true;
256 }
257 
258 //_____________________________________________________________________
259 bool HijingHadronizer::call_hijset(double efrm, std::string frame, std::string proj, std::string targ, int iap, int izp, int iat, int izt)
260 {
261 
262  float ef = efrm;
263  // initialize hydjet
264  HIJSET(ef,frame.data(),proj.data(),targ.data(),iap,izp,iat,izt,strlen(frame.data()),strlen(proj.data()),strlen(targ.data()));
265  return true;
266 }
267 
268 //______________________________________________________________
270 
271  //initialize pythia5
272 
273  // std::string dumstr = "";
274  // call_pygive(dumstr);
275 
276  // initialize hijing
277  LogInfo("HIJINGinAction") << "##### Calling HIJSET(" << efrm_ << "," <<frame_<<","<<proj_<<","<<targ_<<","<<iap_<<","<<izp_<<","<<iat_<<","<<izt_<<") ####";
279 
280  return true;
281 
282 }
283 
284 bool HijingHadronizer::declareStableParticles( const std::vector<int>& pdg )
285 {
286  return true;
287 }
288 
289 //________________________________________________________________
291 
292  phi0_ = 2.*pi*gen::hijran_(0) - pi;
293  sinphi0_ = sin(phi0_);
294  cosphi0_ = cos(phi0_);
295 }
296 
297 //________________________________________________________________
299 {
300  return false;
301 }
302 
304 {
305  return true;
306 }
307 
309 {
310  return true;
311 }
312 
314  return;
315 }
316 
318  return;
319 }
320 
321 const char* HijingHadronizer::classname() const
322 {
323  return "gen::HijingHadronizer";
324 }
325 
#define LogDebug(id)
HepMC::GenEvent * evt
HepMC::GenParticle * build_hijing(int index, int barcode)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double v[5][pyjets_maxn]
std::auto_ptr< HepMC::GenEvent > & event()
static const std::string kPythia6
void add_heavy_ion_rec(HepMC::GenEvent *evt)
float rlu_(unsigned int *iseed)
double p[5][pyjets_maxn]
#define HIJSET
Definition: HijingWrapper.h:13
T sqrt(T t)
Definition: SSEVec.h:18
#define himain2
Definition: HijingWrapper.h:42
HepMC::GenVertex * build_hijing_vertex(int i, int id)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static CLHEP::HepRandomEngine * hijRandomEngine
static const double pi
const char * classname() const
bool get_particles(HepMC::GenEvent *evt)
float ran_(unsigned int *iseed)
int iseed
Definition: AMPTWrapper.h:124
#define HIJING
Definition: HijingWrapper.h:18
part
Definition: HCALResponse.h:20
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
HLT enums.
#define hiparnt
Definition: HijingWrapper.h:52
bool call_hijset(double efrm, std::string frame, std::string proj, std::string targ, int iap, int izp, int iat, int izt)
bool declareStableParticles(const std::vector< int > &)
float hijran_(int *)
#define himain1
Definition: HijingWrapper.h:33