CMS 3D CMS Logo

PyquenHadronizer.h
Go to the documentation of this file.
1 #ifndef Pyquen_Hadronizer_h
2 #define Pyquen_Hadronizer_h
3 
15 
21 
22 #include <map>
23 #include <string>
24 #include <vector>
25 
26 #include "HepMC/GenEvent.h"
27 
28 namespace CLHEP {
29  class HepRandomEngine;
30 }
31 
32 namespace gen {
33  class Pythia6Service;
34 
36  public:
38  ~PyquenHadronizer() override;
39 
41  bool hadronize();
42  bool decay();
43  bool residualDecay();
44  bool readSettings(int);
47  bool declareStableParticles(const std::vector<int>&);
48  bool declareSpecialSettings(const std::vector<std::string>&) { return true; }
49  bool select(HepMC::GenEvent* evtTry) const override { return selector_->filter(evtTry); }
50  void finalizeEvent();
51  void statistics();
52  const char* classname() const;
53 
54  private:
55  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override;
56  std::vector<std::string> const& doSharedResources() const override { return theSharedResources; }
57 
58  static const std::vector<std::string> theSharedResources;
59 
61 
63  bool pyquen_init(const edm::ParameterSet& pset);
64  const char* nucleon();
65  void rotateEvtPlane(HepMC::GenEvent* evt, double angle);
66 
68  double abeamtarget_;
69  unsigned int angularspecselector_;
70  double bmin_;
74  double bmax_;
75  double bfixed_;
76  int cflag_;
77  double comenergy;
78  bool doquench_;
81  bool doIsospin_;
83  bool embedding_;
84  double evtPlane_;
85  double pfrac_;
86 
87  unsigned int nquarkflavor_;
88  double qgpt0_;
90  double qgptau0_;
92  unsigned int maxEventsToPrint_;
94 
95  HepMC::FourVector* fVertex_;
96  std::vector<double> signalVtx_;
97 
99  unsigned int pythiaPylistVerbosity_;
100 
101  // CLHEP::HepRandomEngine* fRandomEngine;
106  };
107 } // namespace gen
108 
109 #endif
const char * classname() const
unsigned int maxEventsToPrint_
Events to print if verbosity.
bool pyquen_init(const edm::ParameterSet &pset)
bool declareStableParticles(const std::vector< int > &)
bool doradiativeenloss_
DEFAULT = true.
bool docollisionalenloss_
DEFAULT = true.
double v[5][pyjets_maxn]
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
bool initializeForExternalPartons()
double bmax_
max impact param (fm); valid only if cflag_!=0
double abeamtarget_
beam/target atomic mass number
bool declareSpecialSettings(const std::vector< std::string > &)
int cflag_
centrality flag =0 fixed impact param, <>0 minbias
unsigned int pythiaPylistVerbosity_
Pythia PYLIST Verbosity flag.
bool doIsospin_
Run n&p with proper ratios; if false, only p+p collisions.
unsigned int angularspecselector_
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > src_
virtual bool filter(HepMC::GenEvent *)
HepMC::FourVector * fVertex_
Event signal vertex.
void add_heavy_ion_rec(HepMC::GenEvent *evt)
Pythia6Service * pythia6Service_
PyquenHadronizer(const edm::ParameterSet &, edm::ConsumesCollector &&)
edm::ParameterSet pset_
bool pyqpythia_init(const edm::ParameterSet &pset)
unsigned int nquarkflavor_
bool doquench_
if true perform quenching (default = true)
std::vector< double > signalVtx_
Pset double vector to set event signal vertex.
double comenergy
collision energy
bool select(HepMC::GenEvent *evtTry) const override
bool pythiaHepMCVerbosity_
HepMC verbosity flag.
static const std::vector< std::string > theSharedResources
BaseHiGenEvtSelector * selector_
double bmin_
min impact param (fm); valid only if cflag_!=0
std::vector< std::string > const & doSharedResources() const override
double bfixed_
fixed impact param (fm); valid only if cflag_=0
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
double pfrac_
Proton fraction in the nucleus.
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11