CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HydjetHadronizer.h
Go to the documentation of this file.
1 #ifndef HydjetHadronizer_h
2 #define HydjetHadronizer_h
3 
11 #define PYCOMP pycomp_
12 
15 
21 
22 #include <map>
23 #include <string>
24 #include <vector>
25 #include <cmath>
26 
27 namespace CLHEP {
28  class HepRandomEngine;
29 }
30 
31 namespace HepMC {
32  class GenEvent;
33  class GenParticle;
34  class GenVertex;
35  class FourVector;
36 } // namespace HepMC
37 
38 namespace gen {
39  class Pythia6Service;
40 
42  public:
44  ~HydjetHadronizer() override;
45 
47  bool hadronize();
48  bool decay();
49  bool residualDecay();
50  bool readSettings(int);
53  bool declareStableParticles(const std::vector<int>&);
54  bool declareSpecialSettings(const std::vector<std::string>&) { return true; }
55 
56  void finalizeEvent();
57  void statistics();
58  const char* classname() const;
59 
60  private:
61  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override;
62  std::vector<std::string> const& doSharedResources() const override { return theSharedResources; }
63 
64  static const std::vector<std::string> theSharedResources;
65 
67  HepMC::GenParticle* build_hyjet(int index, int barcode);
68  HepMC::GenVertex* build_hyjet_vertex(int i, int id);
70  bool call_hyinit(double energy, double a, int ifb, double bmin, double bmax, double bfix, int nh);
71  bool hydjet_init(const edm::ParameterSet& pset);
72  inline double nuclear_radius() const;
73  void rotateEvtPlane();
74  int convertStatus(int st);
75 
78  double abeamtarget_;
80  double bfixed_;
81  double bmax_;
82  double bmin_;
84  int cflag_;
86  bool embedding_;
89  double comenergy;
92  double fracsoftmult_;
93  double hadfreeztemp_;
101  unsigned int maxEventsToPrint_;
102  double maxlongy_;
103  double maxtrany_;
106  int nsub_;
109  int nhard_;
111  int nsoft_;
113  unsigned int nquarkflavor_;
114  unsigned int pythiaPylistVerbosity_;
116  double qgpt0_;
117  double qgptau0_;
119 
121  double phi0_;
122  double sinphi0_;
123  double cosphi0_;
124  bool rotate_;
125 
126  unsigned int shadowingswitch_;
127  double signn_;
129  HepMC::FourVector* fVertex_;
131  std::vector<double> signalVtx_;
132 
135  };
136 
141  double HydjetHadronizer::nuclear_radius() const { return 1.15 * pow((double)abeamtarget_, 1. / 3.); }
142 } // namespace gen
143 #endif
double bfixed_
fixed impact param (fm); valid only if cflag_=0
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
HepMC::GenEvent * evt
bool docollisionalenloss_
DEFAULT = true.
bool hydjet_init(const edm::ParameterSet &pset)
void add_heavy_ion_rec(HepMC::GenEvent *evt)
double comenergy
collision energy
double v[5][pyjets_maxn]
bool embedding_
Switch for embedding mode.
std::string hymode_
Hydjet running mode.
int nsoft_
multiplicity of HYDRO-induced particles in event
bool initializeForExternalPartons()
const char * classname() const
unsigned int maxEventsToPrint_
Events to print if verbosity.
uint32_t nh
HepMC::GenParticle * build_hyjet(int index, int barcode)
Pythia6Service * pythia6Service_
bool declareSpecialSettings(const std::vector< std::string > &)
int nsub_
number of sub-events
std::vector< double > signalVtx_
Pset double vector to set event signal vertex.
int nhard_
multiplicity of PYTHIA(+PYQUEN)-induced particles in event
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_
pythia verbosity; def=1
double abeamtarget_
beam/target atomic mass number
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > src_
double a
Definition: hdecay.h:119
edm::ParameterSet pset_
unsigned int shadowingswitch_
bool doradiativeenloss_
DEFAULT = true.
HydjetHadronizer(const edm::ParameterSet &, edm::ConsumesCollector &&)
int angularspecselector_
angular emitted gluon spectrum selection
bool declareStableParticles(const std::vector< int > &)
bool rotate_
Switch to rotate event plane.
double phi0_
Event plane angle.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
HepMC::FourVector * fVertex_
Event signal vertex.
static const std::vector< std::string > theSharedResources
bool get_particles(HepMC::GenEvent *evt)
std::vector< std::string > const & doSharedResources() const override