CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/GeneratorInterface/HydjetInterface/interface/HydjetHadronizer.h

Go to the documentation of this file.
00001 #ifndef HydjetHadronizer_h
00002 #define HydjetHadronizer_h
00003 
00004 // $Id: HydjetHadronizer.h,v 1.9 2011/02/17 20:54:45 yarba Exp $
00005 
00014 #define PYCOMP pycomp_
00015 
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00018 
00019 #include <map>
00020 #include <string>
00021 #include <vector>
00022 #include <math.h>
00023 
00024 namespace HepMC {
00025   class GenEvent;
00026   class GenParticle;
00027   class GenVertex;
00028 }
00029 
00030 namespace gen
00031 {
00032    class Pythia6Service;
00033 
00034   class HydjetHadronizer : public BaseHadronizer {
00035   public:
00036     HydjetHadronizer(const edm::ParameterSet &);
00037     virtual ~HydjetHadronizer();
00038 
00039     bool generatePartonsAndHadronize();
00040     bool hadronize();
00041     bool decay();
00042     bool residualDecay();
00043     bool readSettings( int );
00044     bool initializeForExternalPartons();
00045     bool initializeForInternalPartons();
00046     bool declareStableParticles( const std::vector<int> );
00047     bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
00048 
00049     void finalizeEvent();
00050     void statistics();
00051     const char* classname() const;
00052 
00053   private:
00054     void                                        add_heavy_ion_rec(HepMC::GenEvent *evt);
00055     HepMC::GenParticle*                         build_hyjet( int index, int barcode );  
00056     HepMC::GenVertex*                           build_hyjet_vertex(int i, int id);
00057     bool                                        get_particles(HepMC::GenEvent* evt);
00058     bool                                        call_hyinit(double energy, double a, int ifb, double bmin,
00059                                                             double bmax,double bfix,int nh);
00060     bool                                        hydjet_init(const edm::ParameterSet &pset);
00061     inline double                               nuclear_radius() const;
00062     void                                        rotateEvtPlane();
00063 
00064     HepMC::GenEvent   *evt;
00065     edm::ParameterSet  pset_;
00066     double             abeamtarget_;           // beam/target atomic mass number 
00067     double             bfixed_;                // fixed impact param (fm); valid only if cflag_=0
00068     double             bmax_;                  // max impact param; 
00069                                                // units of nucl radius
00070     double             bmin_;                  // min impact param; 
00071                                                // units of nucl radius
00072     int                cflag_;                 // centrality flag 
00073                                                // =  0 fixed impact param, 
00074                                                // <> 0 between bmin and bmax
00075     bool             embedding_;               // Switch for embedding mode
00076     double             comenergy;              // collision energy   
00077     bool               doradiativeenloss_;     
00078     bool               docollisionalenloss_;   
00079     double             fracsoftmult_;          // fraction of soft hydro induced hadronic multiplicity
00080                                                // proportional to no of nucleon participants
00081                                                // (1-fracsoftmult_)--- fraction of soft 
00082                                                // multiplicity proportional to the numebr 
00083                                                // of nucleon-nucleon binary collisions
00084                                                // DEFAULT=1., allowed range [0.01,1]
00085     double             hadfreeztemp_;          // hadron freez-out temperature
00086                                                // DEFAULT=0.14MeV, allowed ranges [0.08,0.2]MeV
00087     std::string        hymode_;                // Hydjet running mode
00088     unsigned int       maxEventsToPrint_;      // Events to print if verbosity  
00089     double             maxlongy_;              // max longitudinal collective rapidity: 
00090                                                // controls width of eta-spectra
00091                                                // DEFAULT=4, allowed range [0.01,7.0]
00092     double             maxtrany_;              // max transverse collective rapidity: 
00093                                                // controls slope of low-pt spectra
00094                                                // DEFAULT=1.5, allowed range [0.01,3.0]
00095     int                nsub_;                  // number of sub-events
00096     int                nhard_;                 // multiplicity of PYTHIA(+PYQUEN)-induced particles in event
00097     int                nmultiplicity_;         // mean soft multiplicity in central PbPb
00098                                                // automatically calculated for other centralitie and beams         
00099     int                nsoft_;                 // multiplicity of HYDRO-induced particles in event 
00100     unsigned int       nquarkflavor_;          
00101 
00102     unsigned int       pythiaPylistVerbosity_; // pythia verbosity; def=1 
00103     double             qgpt0_;                 // initial temperature of QGP
00104                                                // DEFAULT = 1GeV; allowed range [0.2,2.0]GeV; 
00105     double             qgptau0_;               // proper time of QGP formation
00106                                                // DEFAULT = 0.1 fm/c; allowed range [0.01,10.0]fm/
00107 
00108     double             phi0_;                  // Event plane angle
00109     double             sinphi0_;
00110     double             cosphi0_;
00111     bool               rotate_;                // Switch to rotate event plane
00112 
00113     unsigned int       shadowingswitch_;       // shadowing switcher
00114                                                // 1-ON, 0-OFF
00115     double             signn_;                 // inelastic nucleon nucleon cross section [mb]
00116                                                // DEFAULT= 58 mb
00117 
00118     Pythia6Service* pythia6Service_;
00119     edm::InputTag   src_;
00120   };
00121 
00122   double HydjetHadronizer::nuclear_radius() const
00123   {
00124     // Return the nuclear radius derived from the 
00125     // beam/target atomic mass number.
00126 
00127     return 1.15 * pow((double)abeamtarget_, 1./3.);
00128   }
00129 } /*end namespace*/
00130 #endif