CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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.11 2010/02/24 05:05:03 wmtan 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 embedding_(pset.getParameter<bool>("embeddingMode")),
00050 evtPlane_(0),
00051 filterType_(pset.getUntrackedParameter<string>("filterType","None")),
00052 maxTries_(pset.getUntrackedParameter<int>("maxTries",1000)),
00053 nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
00054 qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
00055 qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
00056 maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint",1)),
00057 pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
00058 pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
00059 pythia6Service_(new Pythia6Service(pset))
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 
00077   selector_ = HiGenEvtSelectorFactory::get(filterType_,pset);
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    int counter = 0;
00141    bool pass = false;
00142 
00143    HepMC::GenEvent* evt = 0;
00144    while(!pass){
00145       if(counter == maxTries_) throw edm::Exception(edm::errors::Configuration,"InfiniteLoop")<<"Pyquen tried "<<counter<<" times to generate event with "<<filterType_.data()<<" ."<<endl;
00146       
00147       // Generate PYQUEN event
00148       // generate single partonic PYTHIA jet event
00149 
00150       // Take into account whether it's a nn or pp or pn interaction
00151       if(doIsospin_) call_pyinit("CMS", nucleon(), nucleon(), comenergy);
00152       call_pyevnt();
00153 
00154       // call PYQUEN to apply parton rescattering and energy loss 
00155       // if doQuench=FALSE, it is pure PYTHIA
00156       if( doquench_ ){
00157          PYQUEN(abeamtarget_,cflag_,bfixed_,bmin_,bmax_);
00158          edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN("<<abeamtarget_<<","<<cflag_<<","<<bfixed_<<") ####";
00159       } else {
00160          edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
00161       }
00162       
00163       // call PYTHIA to finish the hadronization
00164       pyexec_();
00165 
00166       // fill the HEPEVT with the PYJETS event record
00167       call_pyhepc(1);
00168       
00169       // event information
00170       evt = hepevtio.read_next_event();
00171       pass = selector_->filter(evt);
00172       counter++;
00173    }
00174 
00175    evt->set_signal_process_id(pypars.msti[0]);      // type of the process
00176    evt->set_event_scale(pypars.pari[16]);           // Q^2
00177    
00178    if(embedding_) rotateEvtPlane(evt,evtPlane_);
00179    add_heavy_ion_rec(evt);
00180    
00181    event().reset(evt);
00182 
00183   return true;
00184 }
00185 
00186 bool PyquenHadronizer::initializeForInternalPartons(){
00187 
00188    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00189    pythia6Service_->setGeneralParams();
00190 
00191    //Proton to Nucleon fraction
00192    pfrac_ = 1./(1.98+0.015*pow(abeamtarget_,2./3));
00193 
00194    //initialize pythia         
00195    pyqpythia_init(pset_);
00196 
00197    //initilize pyquen          
00198    pyquen_init(pset_);
00199 
00200    // Call PYTHIA              
00201    call_pyinit("CMS", "p", "p", comenergy);
00202 
00203    return true;
00204 }
00205 
00206 
00207 //_____________________________________________________________________
00208 bool PyquenHadronizer::pyqpythia_init(const ParameterSet & pset)
00209 {
00210 
00211    //Turn Hadronization Off if there is quenching  
00212    if(doquench_){
00213       string sHadOff("MSTP(111)=0");
00214       gen::call_pygive(sHadOff);
00215    }
00216   return true;
00217 }
00218 
00219 
00220 //_________________________________________________________________
00221 bool PyquenHadronizer::pyquen_init(const ParameterSet &pset)
00222 {
00223   // PYQUEN initialization
00224 
00225   // angular emitted gluon  spectrum selection 
00226   pyqpar.ianglu = angularspecselector_;
00227 
00228   // type of medium induced partonic energy loss
00229   if( doradiativeenloss_ && docollisionalenloss_ ){
00230     edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
00231     pyqpar.ienglu = 0; 
00232   } else if ( doradiativeenloss_ ) {
00233     edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
00234     pyqpar.ienglu = 1; 
00235   } else if ( docollisionalenloss_ ) {
00236     edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
00237     pyqpar.ienglu = 2; 
00238   } else {
00239     edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
00240     pyqpar.ienglu = 0; 
00241   }
00242 
00243   // number of active quark flavors in qgp
00244   pyqpar.nfu    = nquarkflavor_;
00245   // initial temperature of QGP
00246   pyqpar.T0u    = qgpt0_;
00247   // proper time of QGP formation
00248   pyqpar.tau0u  = qgptau0_;
00249 
00250   return true;
00251 }
00252 
00253 const char* PyquenHadronizer::nucleon(){
00254   int* dummy = 0;
00255   double random = gen::pyr_(dummy);
00256   const char* nuc = 0;
00257   if(random > pfrac_) nuc = "n";
00258   else nuc = "p";
00259   
00260   return nuc;
00261 }
00262 
00263 void PyquenHadronizer::rotateEvtPlane(HepMC::GenEvent* evt, double angle){
00264 
00265    double sinphi0 = sin(angle);
00266    double cosphi0 = cos(angle);
00267 
00268    for ( HepMC::GenEvent::vertex_iterator vt=evt->vertices_begin();
00269          vt!=evt->vertices_end(); ++vt )
00270       {
00271          
00272          double x0 = (*vt)->position().x();
00273          double y0 = (*vt)->position().y();
00274          double z = (*vt)->position().z();
00275          double t = (*vt)->position().t();
00276 
00277          double x = x0*cosphi0-y0*sinphi0;
00278          double y = y0*cosphi0+x0*sinphi0;
00279 
00280          (*vt)->set_position( HepMC::FourVector(x,y,z,t) ) ;      
00281       }
00282 
00283    for ( HepMC::GenEvent::particle_iterator vt=evt->particles_begin();
00284          vt!=evt->particles_end(); ++vt )
00285       {
00286 
00287          double x0 = (*vt)->momentum().x();
00288          double y0 = (*vt)->momentum().y();
00289          double z = (*vt)->momentum().z();
00290          double t = (*vt)->momentum().t();
00291 
00292          double x = x0*cosphi0-y0*sinphi0;
00293          double y = y0*cosphi0+x0*sinphi0;
00294 
00295          (*vt)->set_momentum( HepMC::FourVector(x,y,z,t) ) ;
00296       }
00297 }
00298 
00299 bool PyquenHadronizer::declareStableParticles( std::vector<int> pdg )
00300 {
00301    for ( size_t i=0; i < pdg.size(); i++ )
00302       {
00303          int pyCode = pycomp_( pdg[i] );
00304          std::ostringstream pyCard ;
00305          pyCard << "MDCY(" << pyCode << ",1)=0";
00306          std::cout << pyCard.str() << std::endl;
00307          call_pygive( pyCard.str() );
00308       }
00309 
00310    return true;
00311 
00312 }
00313 
00314 
00315 
00316 //____________________________________________________________________
00317 
00318 bool PyquenHadronizer::hadronize()
00319 {
00320    return false;
00321 }
00322 
00323 bool PyquenHadronizer::decay()
00324 {
00325    return true;
00326 }
00327 
00328 bool PyquenHadronizer::residualDecay()
00329 {
00330    return true;
00331 }
00332 
00333 void PyquenHadronizer::finalizeEvent(){
00334 
00335 }
00336 
00337 void PyquenHadronizer::statistics(){
00338 }
00339 
00340 const char* PyquenHadronizer::classname() const
00341 {
00342    return "gen::PyquenHadronizer";
00343 }
00344 
00345 
00346