CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/GeneratorInterface/PyquenInterface/src/PyquenHadronizer.cc

Go to the documentation of this file.
00001 /*
00002  *
00003  * Generates PYQUEN HepMC events
00004  *
00005  * $Id: PyquenHadronizer.cc,v 1.15 2013/01/20 23:55:21 mnguyen Exp $
00006 */
00007 
00008 #include <iostream>
00009 #include "time.h"
00010 
00011 #include "GeneratorInterface/PyquenInterface/interface/PyquenHadronizer.h"
00012 #include "GeneratorInterface/PyquenInterface/interface/PyquenWrapper.h"
00013 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Declarations.h"
00014 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00015 
00016 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00017 
00018 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00023 
00024 #include "GeneratorInterface/HiGenCommon/interface/HiGenEvtSelectorFactory.h"
00025 
00026 #include "HepMC/IO_HEPEVT.h"
00027 #include "HepMC/PythiaWrapper.h"
00028 
00029 using namespace gen;
00030 using namespace edm;
00031 using namespace std;
00032 
00033 HepMC::IO_HEPEVT hepevtio;
00034 
00035 PyquenHadronizer :: PyquenHadronizer(const ParameterSet & pset):
00036    BaseHadronizer(pset),
00037    pset_(pset),
00038 abeamtarget_(pset.getParameter<double>("aBeamTarget")),
00039 angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
00040 bmin_(pset.getParameter<double>("bMin")),
00041 bmax_(pset.getParameter<double>("bMax")),
00042 bfixed_(pset.getParameter<double>("bFixed")),
00043 cflag_(pset.getParameter<int>("cFlag")),
00044 comenergy(pset.getParameter<double>("comEnergy")),
00045 doquench_(pset.getParameter<bool>("doQuench")),
00046 doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
00047 docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
00048 doIsospin_(pset.getParameter<bool>("doIsospin")),
00049    protonSide_(pset.getUntrackedParameter<int>("protonSide",0)),
00050 embedding_(pset.getParameter<bool>("embeddingMode")),
00051 evtPlane_(0),
00052 nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
00053 qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
00054 qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
00055 maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
00056 pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
00057 pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
00058 pythia6Service_(new Pythia6Service(pset)),
00059 filterType_(pset.getUntrackedParameter<string>("filterType","None"))
00060 {
00061   // Verbosity Level
00062   // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
00063   LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl;
00064   // HepMC event verbosity Level
00065   pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false);
00066   LogDebug("HepMCverbosity")  << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl; 
00067 
00068   //Max number of events printed on verbosity level 
00069   maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
00070   LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl;
00071 
00072   if(embedding_){ 
00073      cflag_ = 0;
00074      src_ = pset.getParameter<InputTag>("backgroundLabel");
00075   }
00076   selector_ = HiGenEvtSelectorFactory::get(filterType_,pset);
00077 
00078 }
00079 
00080 
00081 //_____________________________________________________________________
00082 PyquenHadronizer::~PyquenHadronizer()
00083 {
00084   // distructor
00085   call_pystat(1);
00086 
00087   delete pythia6Service_;
00088 
00089 }
00090 
00091 
00092 //_____________________________________________________________________
00093 void PyquenHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
00094 {
00095   HepMC::HeavyIon *hi = new HepMC::HeavyIon(
00096     1,                                 // Ncoll_hard
00097     -1,                                 // Npart_proj
00098     -1,                                 // Npart_targ
00099     1,                                 // Ncoll
00100     -1,                                 // spectator_neutrons
00101     -1,                                 // spectator_protons
00102     -1,                                 // N_Nwounded_collisions
00103     -1,                                 // Nwounded_N_collisions
00104     -1,                                 // Nwounded_Nwounded_collisions
00105     plfpar.bgen,                        // impact_parameter in [fm]
00106     evtPlane_,                                  // event_plane_angle
00107     0,                                  // eccentricity
00108     -1                                  // sigma_inel_NN
00109   );
00110 
00111   evt->set_heavy_ion(*hi);
00112 
00113   delete hi;
00114 }
00115 
00116 //_____________________________________________________________________
00117 bool PyquenHadronizer::generatePartonsAndHadronize()
00118 {
00119    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00120 
00121    // Not possible to retrieve impact paramter and event plane info
00122    // at this part, need to overwrite filter() in 
00123    // PyquenGeneratorFilter 
00124 
00125    const edm::Event& e = getEDMEvent();
00126 
00127    if(embedding_){
00128       Handle<HepMCProduct> input;
00129       e.getByLabel(src_,input);
00130       const HepMC::GenEvent * inev = input->GetEvent();
00131       const HepMC::HeavyIon* hi = inev->heavy_ion();
00132       if(hi){
00133          bfixed_ = hi->impact_parameter();
00134          evtPlane_ = hi->event_plane_angle();
00135       }else{
00136          LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
00137       }
00138    }
00139 
00140    // Generate PYQUEN event
00141    // generate single partonic PYTHIA jet event
00142    
00143    // Take into account whether it's a nn or pp or pn interaction
00144    if(doIsospin_){
00145      string projN = "p";
00146      string targN = "p";
00147      if(protonSide_ != 1) projN = nucleon();
00148      if(protonSide_ != 2) targN = nucleon();
00149      call_pyinit("CMS", projN.data(), targN.data(), comenergy);
00150    }
00151    call_pyevnt();
00152    
00153    // call PYQUEN to apply parton rescattering and energy loss 
00154    // if doQuench=FALSE, it is pure PYTHIA
00155    if( doquench_ ){
00156      PYQUEN(abeamtarget_,cflag_,bfixed_,bmin_,bmax_);
00157      edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
00158    } else {
00159      edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
00160    }
00161    
00162    // call PYTHIA to finish the hadronization
00163    pyexec_();
00164    
00165    // fill the HEPEVT with the PYJETS event record
00166    call_pyhepc(1);
00167    
00168    // event information
00169     HepMC::GenEvent* evt = hepevtio.read_next_event();
00170 
00171    evt->set_signal_process_id(pypars.msti[0]);      // type of the process
00172    evt->set_event_scale(pypars.pari[16]);           // Q^2
00173    
00174    if(embedding_) rotateEvtPlane(evt,evtPlane_);
00175    add_heavy_ion_rec(evt);
00176    
00177    event().reset(evt);
00178 
00179   return true;
00180 }
00181 
00182 bool PyquenHadronizer::readSettings( int )
00183 {
00184 
00185    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00186    pythia6Service_->setGeneralParams();
00187    pythia6Service_->setCSAParams();
00188 
00189    //Proton to Nucleon fraction
00190    pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
00191 
00192    //initialize pythia         
00193    pyqpythia_init(pset_);
00194 
00195    //initilize pyquen          
00196    pyquen_init(pset_);
00197 
00198    return true;
00199 
00200 }
00201 
00202 bool PyquenHadronizer::initializeForInternalPartons(){
00203 
00204 
00205    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00206 
00207    // Call PYTHIA              
00208    call_pyinit("CMS", "p", "p", comenergy);
00209 
00210    return true;
00211 }
00212 
00213 
00214 //_____________________________________________________________________
00215 bool PyquenHadronizer::pyqpythia_init(const ParameterSet & pset)
00216 {
00217 
00218   //Turn Hadronization Off whether or not there is quenching  
00219   // PYEXEC is called later anyway
00220   string sHadOff("MSTP(111)=0");
00221   gen::call_pygive(sHadOff);
00222 
00223   return true;
00224 }
00225 
00226 
00227 //_________________________________________________________________
00228 bool PyquenHadronizer::pyquen_init(const ParameterSet &pset)
00229 {
00230   // PYQUEN initialization
00231 
00232   // angular emitted gluon  spectrum selection 
00233   pyqpar.ianglu = angularspecselector_;
00234 
00235   // type of medium induced partonic energy loss
00236   if( doradiativeenloss_ && docollisionalenloss_ ){
00237     edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
00238     pyqpar.ienglu = 0; 
00239   } else if ( doradiativeenloss_ ) {
00240     edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
00241     pyqpar.ienglu = 1; 
00242   } else if ( docollisionalenloss_ ) {
00243     edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
00244     pyqpar.ienglu = 2; 
00245   } else {
00246     edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
00247     pyqpar.ienglu = 0; 
00248   }
00249 
00250   // number of active quark flavors in qgp
00251   pyqpar.nfu    = nquarkflavor_;
00252   // initial temperature of QGP
00253   pyqpar.T0u    = qgpt0_;
00254   // proper time of QGP formation
00255   pyqpar.tau0u  = qgptau0_;
00256 
00257   return true;
00258 }
00259 
00260 const char* PyquenHadronizer::nucleon(){
00261   int* dummy = 0;
00262   double random = gen::pyr_(dummy);
00263   const char* nuc = 0;
00264   if(random > pfrac_) nuc = "n";
00265   else nuc = "p";
00266   
00267   return nuc;
00268 }
00269 
00270 void PyquenHadronizer::rotateEvtPlane(HepMC::GenEvent* evt, double angle){
00271 
00272    double sinphi0 = sin(angle);
00273    double cosphi0 = cos(angle);
00274 
00275    for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
00276          vt!=evt->vertices_end(); ++vt )
00277       {
00278          
00279          double x0 = (*vt)->position().x();
00280          double y0 = (*vt)->position().y();
00281          double z = (*vt)->position().z();
00282          double t = (*vt)->position().t();
00283 
00284          double x = x0*cosphi0-y0*sinphi0;
00285          double y = y0*cosphi0+x0*sinphi0;
00286 
00287          (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;      
00288       }
00289 
00290    for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
00291          vt!=evt->particles_end(); ++vt )
00292       {
00293 
00294          double x0 = (*vt)->momentum().x();
00295          double y0 = (*vt)->momentum().y();
00296          double z = (*vt)->momentum().z();
00297          double t = (*vt)->momentum().t();
00298 
00299          double x = x0*cosphi0-y0*sinphi0;
00300          double y = y0*cosphi0+x0*sinphi0;
00301 
00302          (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
00303       }
00304 }
00305 
00306 bool PyquenHadronizer::declareStableParticles( std::vector<int> pdg )
00307 {
00308    for ( size_t i=0; i < pdg.size(); i++ )
00309       {
00310          int pyCode = pycomp_( pdg[i] );
00311          std::ostringstream pyCard ;
00312          pyCard << "MDCY(" << pyCode << ",1)=0";
00313          std::cout << pyCard.str() << std::endl;
00314          call_pygive( pyCard.str() );
00315       }
00316 
00317    return true;
00318 
00319 }
00320 
00321 
00322 
00323 //____________________________________________________________________
00324 
00325 bool PyquenHadronizer::hadronize()
00326 {
00327    return false;
00328 }
00329 
00330 bool PyquenHadronizer::decay()
00331 {
00332    return true;
00333 }
00334 
00335 bool PyquenHadronizer::residualDecay()
00336 {
00337    return true;
00338 }
00339 
00340 void PyquenHadronizer::finalizeEvent(){
00341 
00342 }
00343 
00344 void PyquenHadronizer::statistics(){
00345 }
00346 
00347 const char* PyquenHadronizer::classname() const
00348 {
00349    return "gen::PyquenHadronizer";
00350 }
00351 
00352 
00353