CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HydjetHadronizer.h
Go to the documentation of this file.
1 #ifndef HydjetHadronizer_h
2 #define HydjetHadronizer_h
3 
4 
13 #define PYCOMP pycomp_
14 
17 
18 #include <map>
19 #include <string>
20 #include <vector>
21 #include <math.h>
22 
23 namespace HepMC {
24  class GenEvent;
25  class GenParticle;
26  class GenVertex;
27 }
28 
29 namespace gen
30 {
31  class Pythia6Service;
32 
34  public:
36  virtual ~HydjetHadronizer();
37 
39  bool hadronize();
40  bool decay();
41  bool residualDecay();
42  bool readSettings( int );
45  bool declareStableParticles( const std::vector<int>& );
46  bool declareSpecialSettings( const std::vector<std::string>& ) { return true; }
47 
48  void finalizeEvent();
49  void statistics();
50  const char* classname() const;
51 
52  private:
53  void add_heavy_ion_rec(HepMC::GenEvent *evt);
54  HepMC::GenParticle* build_hyjet( int index, int barcode );
55  HepMC::GenVertex* build_hyjet_vertex(int i, int id);
56  bool get_particles(HepMC::GenEvent* evt);
57  bool call_hyinit(double energy, double a, int ifb, double bmin,
58  double bmax,double bfix,int nh);
59  bool hydjet_init(const edm::ParameterSet &pset);
60  inline double nuclear_radius() const;
61  void rotateEvtPlane();
62 
63  HepMC::GenEvent *evt;
65  double abeamtarget_; // beam/target atomic mass number
66  double bfixed_; // fixed impact param (fm); valid only if cflag_=0
67  double bmax_; // max impact param;
68  // units of nucl radius
69  double bmin_; // min impact param;
70  // units of nucl radius
71  int cflag_; // centrality flag
72  // = 0 fixed impact param,
73  // <> 0 between bmin and bmax
74  bool embedding_; // Switch for embedding mode
75  double comenergy; // collision energy
78  double fracsoftmult_; // fraction of soft hydro induced hadronic multiplicity
79  // proportional to no of nucleon participants
80  // (1-fracsoftmult_)--- fraction of soft
81  // multiplicity proportional to the numebr
82  // of nucleon-nucleon binary collisions
83  // DEFAULT=1., allowed range [0.01,1]
84  double hadfreeztemp_; // hadron freez-out temperature
85  // DEFAULT=0.14MeV, allowed ranges [0.08,0.2]MeV
86  std::string hymode_; // Hydjet running mode
87  unsigned int maxEventsToPrint_; // Events to print if verbosity
88  double maxlongy_; // max longitudinal collective rapidity:
89  // controls width of eta-spectra
90  // DEFAULT=4, allowed range [0.01,7.0]
91  double maxtrany_; // max transverse collective rapidity:
92  // controls slope of low-pt spectra
93  // DEFAULT=1.5, allowed range [0.01,3.0]
94  int nsub_; // number of sub-events
95  int nhard_; // multiplicity of PYTHIA(+PYQUEN)-induced particles in event
96  int nmultiplicity_; // mean soft multiplicity in central PbPb
97  // automatically calculated for other centralitie and beams
98  int nsoft_; // multiplicity of HYDRO-induced particles in event
99  unsigned int nquarkflavor_;
100  unsigned int pythiaPylistVerbosity_; // pythia verbosity; def=1
102  double qgpt0_; // initial temperature of QGP
103  // DEFAULT = 1GeV; allowed range [0.2,2.0]GeV;
104  double qgptau0_; // proper time of QGP formation
105  // DEFAULT = 0.1 fm/c; allowed range [0.01,10.0]fm/
106 
107  double phi0_; // Event plane angle
108  double sinphi0_;
109  double cosphi0_;
110  bool rotate_; // Switch to rotate event plane
111 
112  unsigned int shadowingswitch_; // shadowing switcher
113  // 1-ON, 0-OFF
114  double signn_; // inelastic nucleon nucleon cross section [mb]
115  // DEFAULT= 58 mb
116 
119  };
120 
122  {
123  // Return the nuclear radius derived from the
124  // beam/target atomic mass number.
125 
126  return 1.15 * pow((double)abeamtarget_, 1./3.);
127  }
128 } /*end namespace*/
129 #endif
int i
Definition: DBlmapReader.cc:9
HepMC::GenEvent * evt
bool docollisionalenloss_
DEFAULT = true.
HydjetHadronizer(const edm::ParameterSet &)
bool hydjet_init(const edm::ParameterSet &pset)
void add_heavy_ion_rec(HepMC::GenEvent *evt)
bool initializeForExternalPartons()
const char * classname() const
unsigned int maxEventsToPrint_
HepMC::GenParticle * build_hyjet(int index, int barcode)
Pythia6Service * pythia6Service_
bool declareSpecialSettings(const std::vector< std::string > &)
double fracsoftmult_
DEFAULT = true.
HepMC::GenVertex * build_hyjet_vertex(int i, int id)
double nuclear_radius() const
bool call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh)
unsigned int pythiaPylistVerbosity_
double a
Definition: hdecay.h:121
edm::ParameterSet pset_
unsigned int shadowingswitch_
unsigned int nquarkflavor_
bool declareStableParticles(const std::vector< int > &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
bool get_particles(HepMC::GenEvent *evt)