CMS 3D CMS Logo

Hydjet2Hadronizer.h
Go to the documentation of this file.
1 #ifndef Hydjet2Hadronizer_h
2 #define Hydjet2Hadronizer_h
3 
13 
19 
20 #include <cmath>
21 #include <iostream>
22 #include <map>
23 #include <string>
24 #include <vector>
25 
26 #include <Hydjet2.h>
27 #include <InitialParams.h>
28 
29 namespace CLHEP {
30  class HepRandomEngine;
31  class RandFlat;
32  class RandPoisson;
33  class RandGauss;
34 } // namespace CLHEP
35 
36 namespace HepMC {
37  class GenEvent;
38  class GenParticle;
39  class GenVertex;
40  class FourVector;
41 } // namespace HepMC
42 
43 namespace gen {
44  class Pythia6Service;
46  public:
48  ~Hydjet2Hadronizer() override;
49 
50  bool readSettings(int);
51  bool declareSpecialSettings(const std::vector<std::string> &) { return true; }
55  bool declareStableParticles(const std::vector<int> &);
56 
57  bool hadronize(); //0
58  bool decay(); //0
59  bool residualDecay();
60  void finalizeEvent();
61  void statistics();
62  const char *classname() const;
63 
64  private:
65  void doSetRandomEngine(CLHEP::HepRandomEngine *v) override;
66  void rotateEvtPlane();
68  HepMC::GenParticle *build_hyjet2(int index, int barcode);
69  HepMC::GenVertex *build_hyjet2_vertex(int i, int id);
71 
72  std::vector<std::string> const &doSharedResources() const override { return theSharedResources; }
73  static const std::vector<std::string> theSharedResources;
74 
75  inline double nuclear_radius() const;
76 
77  int convertStatusForComponents(int, int);
78  int convertStatus(int);
79 
80  InitialParamsHydjet_t fParams;
81  Hydjet2 *hj2;
82 
83  bool ev = false;
85  bool rotate_; // Switch to rotate event plane
87  int nsub_; // number of sub-events
88  int nhard_; // multiplicity of PYTHIA(+PYQUEN)-induced particles in event
89  int nsoft_; // multiplicity of HYDRO-induced particles in event
90  double phi0_; // Event plane angle
91  double sinphi0_;
92  double cosphi0_;
93 
94  unsigned int pythiaPylistVerbosity_; // pythia verbosity; def=1
95  unsigned int maxEventsToPrint_; // Events to print if verbosity
96 
98  double Sigin, Sigjet;
99 
100  HepMC::FourVector *fVertex_; // Event signal vertex
101 
102  std::vector<double> signalVtx_; // Pset double vector to set event signal vertex
103 
106  };
107 
109  // Return the nuclear radius derived from the
110  // beam/target atomic mass number.
111 
112  return 1.15 * pow((double)fParams.fAw, 1. / 3.);
113  }
114 } // namespace gen
115 #endif
bool declareStableParticles(const std::vector< int > &)
edm::ParameterSet pset
unsigned int maxEventsToPrint_
InitialParamsHydjet_t fParams
void add_heavy_ion_rec(HepMC::GenEvent *evt)
unsigned int pythiaPylistVerbosity_
HepMC::FourVector * fVertex_
std::vector< double > signalVtx_
static const std::vector< std::string > theSharedResources
double v[5][pyjets_maxn]
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > src_
const char * classname() const
HepMC::GenVertex * build_hyjet2_vertex(int i, int id)
double nuclear_radius() const
Pythia6Service * pythia6Service_
std::vector< std::string > const & doSharedResources() const override
bool get_particles(HepMC::GenEvent *evt)
int convertStatusForComponents(int, int)
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
bool declareSpecialSettings(const std::vector< std::string > &)
Hydjet2Hadronizer(const edm::ParameterSet &, edm::ConsumesCollector &&)
HepMC::GenParticle * build_hyjet2(int index, int barcode)
bool initializeForExternalPartons()
HepMC::GenEvent * evt
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29