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