CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/GeneratorInterface/HydjetInterface/src/HydjetHadronizer.cc

Go to the documentation of this file.
00001 /*
00002  * $Id: HydjetHadronizer.cc,v 1.13 2012/03/10 09:27:19 mnguyen Exp $
00003  *
00004  * Interface to the HYDJET generator, produces HepMC events
00005  *
00006  * Original Author: Camelia Mironov
00007  */
00008 
00009 #include <iostream>
00010 #include <cmath>
00011 
00012 #include "boost/lexical_cast.hpp"
00013 
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/Run.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00020 #include "FWCore/Utilities/interface/EDMException.h"
00021 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00022 
00023 #include "GeneratorInterface/HydjetInterface/interface/HydjetHadronizer.h"
00024 #include "GeneratorInterface/HydjetInterface/interface/HydjetWrapper.h"
00025 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Declarations.h"
00026 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00027 
00028 #include "HepMC/PythiaWrapper6_4.h"
00029 #include "HepMC/GenEvent.h"
00030 #include "HepMC/HeavyIon.h"
00031 #include "HepMC/SimpleVector.h"
00032 
00033 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00034 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00035 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00036 #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
00037 
00038 using namespace edm;
00039 using namespace std;
00040 using namespace gen;
00041 
00042 namespace {
00043    int convertStatus(int st){      
00044       if(st<= 0) return 0;
00045       if(st<=10) return 1;
00046       if(st<=20) return 2;
00047       if(st<=30) return 3;      
00048       else return st;
00049    }
00050 }
00051 
00052 
00053 //_____________________________________________________________________
00054 HydjetHadronizer::HydjetHadronizer(const ParameterSet &pset) :
00055     BaseHadronizer(pset),
00056     evt(0), 
00057     pset_(pset),
00058     abeamtarget_(pset.getParameter<double>("aBeamTarget")),
00059     bfixed_(pset.getParameter<double>("bFixed")),
00060     bmax_(pset.getParameter<double>("bMax")),
00061     bmin_(pset.getParameter<double>("bMin")),
00062     cflag_(pset.getParameter<int>("cFlag")),
00063     embedding_(pset.getParameter<bool>("embeddingMode")),
00064     comenergy(pset.getParameter<double>("comEnergy")),
00065     doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
00066     docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
00067     fracsoftmult_(pset.getParameter<double>("fracSoftMultiplicity")),
00068     hadfreeztemp_(pset.getParameter<double>("hadronFreezoutTemperature")),
00069     hymode_(pset.getParameter<string>("hydjetMode")),
00070     maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)),
00071     maxlongy_(pset.getParameter<double>("maxLongitudinalRapidity")),
00072     maxtrany_(pset.getParameter<double>("maxTransverseRapidity")),
00073     nsub_(0),
00074     nhard_(0),
00075     nmultiplicity_(pset.getParameter<int>("nMultiplicity")),
00076     nsoft_(0),
00077     nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
00078     pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
00079     qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
00080     qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
00081     phi0_(0.),
00082     sinphi0_(0.),
00083     cosphi0_(1.),
00084     rotate_(pset.getParameter<bool>("rotateEventPlane")),
00085     shadowingswitch_(pset.getParameter<int>("shadowingSwitch")),
00086     signn_(pset.getParameter<double>("sigmaInelNN")),
00087     pythia6Service_(new Pythia6Service(pset))
00088 {
00089   // Default constructor
00090 
00091   // PYLIST Verbosity Level
00092   // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
00093   pythiaPylistVerbosity_ = pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0);
00094   LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_;
00095 
00096   //Max number of events printed on verbosity level 
00097   maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint",0);
00098   LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_;
00099 
00100   if(embedding_) src_ = pset.getParameter<edm::InputTag>("backgroundLabel");
00101   
00102 }
00103 
00104 
00105 //_____________________________________________________________________
00106 HydjetHadronizer::~HydjetHadronizer()
00107 {
00108   // destructor
00109   call_pystat(1);
00110   delete pythia6Service_;
00111 }
00112 
00113 
00114 //_____________________________________________________________________
00115 void HydjetHadronizer::add_heavy_ion_rec(HepMC::GenEvent *evt)
00116 {
00117   // heavy ion record in the final CMSSW Event
00118    double npart = hyfpar.npart;
00119    int nproj = static_cast<int>(npart / 2);
00120    int ntarg = static_cast<int>(npart - nproj);
00121 
00122   HepMC::HeavyIon* hi = new HepMC::HeavyIon(
00123     nsub_,                               // Ncoll_hard/N of SubEvents
00124     nproj,                               // Npart_proj
00125     ntarg,                               // Npart_targ
00126     static_cast<int>(hyfpar.nbcol),      // Ncoll
00127     0,                                   // spectator_neutrons
00128     0,                                   // spectator_protons
00129     0,                                   // N_Nwounded_collisions
00130     0,                                   // Nwounded_N_collisions
00131     0,                                   // Nwounded_Nwounded_collisions
00132     hyfpar.bgen * nuclear_radius(),      // impact_parameter in [fm]
00133     phi0_,                                // event_plane_angle
00134     0,                                   // eccentricity
00135     hyjpar.sigin                         // sigma_inel_NN
00136   );
00137 
00138   evt->set_heavy_ion(*hi);
00139   delete hi;
00140 }
00141 
00142 //___________________________________________________________________     
00143 HepMC::GenParticle* HydjetHadronizer::build_hyjet(int index, int barcode)
00144 {
00145    // Build particle object corresponding to index in hyjets (soft+hard)  
00146 
00147    double x0 = hyjets.phj[0][index];
00148    double y0 = hyjets.phj[1][index];
00149 
00150    double x = x0*cosphi0_-y0*sinphi0_;
00151    double y = y0*cosphi0_+x0*sinphi0_;
00152 
00153    HepMC::GenParticle* p = new HepMC::GenParticle(
00154      HepMC::FourVector(x,                                 // px
00155                        y,                                 // py
00156                        hyjets.phj[2][index],              // pz
00157                        hyjets.phj[3][index]),             // E
00158                        hyjets.khj[1][index],              // id
00159                        convertStatus(hyjets.khj[0][index] // status
00160                       )
00161    );
00162 
00163    p->suggest_barcode(barcode);
00164    return p;
00165 }
00166 
00167 //___________________________________________________________________     
00168 HepMC::GenVertex* HydjetHadronizer::build_hyjet_vertex(int i,int id)
00169 {
00170    // build verteces for the hyjets stored events                        
00171 
00172    double x0=hyjets.vhj[0][i];
00173    double y0=hyjets.vhj[1][i];
00174    double x = x0*cosphi0_-y0*sinphi0_;
00175    double y = y0*cosphi0_+x0*sinphi0_;
00176    double z=hyjets.vhj[2][i];
00177    double t=hyjets.vhj[4][i];
00178 
00179    HepMC::GenVertex* vertex = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t),id);
00180    return vertex;
00181 }
00182 
00183 //___________________________________________________________________
00184 
00185 bool HydjetHadronizer::generatePartonsAndHadronize()
00186 {
00187    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00188    
00189    // generate single event
00190    if(embedding_){
00191      cflag_ = 0;
00192      const edm::Event& e = getEDMEvent();
00193      Handle<HepMCProduct> input;
00194      e.getByLabel(src_,input);
00195      const HepMC::GenEvent * inev = input->GetEvent();
00196      const HepMC::HeavyIon* hi = inev->heavy_ion();
00197     if(hi){
00198        bfixed_ = hi->impact_parameter();
00199        phi0_ = hi->event_plane_angle();
00200        sinphi0_ = sin(phi0_);
00201        cosphi0_ = cos(phi0_);
00202      }else{
00203        LogWarning("EventEmbedding")<<"Background event does not have heavy ion record!";
00204      }
00205    }else if(rotate_) rotateEvtPlane();
00206 
00207    nsoft_    = 0;
00208    nhard_    = 0;
00209 
00210    edm::LogInfo("HYDJETmode") << "##### HYDJET  nhsel = " << hyjpar.nhsel;
00211    edm::LogInfo("HYDJETfpart") << "##### HYDJET fpart = " << hyflow.fpart;
00212    edm::LogInfo("HYDJETtf") << "##### HYDJET hadron freez-out temp, Tf = " << hyflow.Tf;
00213    edm::LogInfo("HYDJETinTemp") << "##### HYDJET: QGP init temperature, T0 ="<<pyqpar.T0u;
00214    edm::LogInfo("HYDJETinTau") << "##### HYDJET: QGP formation time,tau0 ="<<pyqpar.tau0u;
00215 
00216    // generate a HYDJET event
00217    int ntry = 0;
00218    while(nsoft_ == 0 && nhard_ == 0){
00219       if(ntry > 100){
00220          edm::LogError("HydjetEmptyEvent") << "##### HYDJET: No Particles generated, Number of tries ="<<ntry;
00221 
00222          // Throw an exception.  Use the EventCorruption exception since it maps onto SkipEvent
00223          // which is what we want to do here.
00224 
00225          std::ostringstream sstr;
00226          sstr << "HydjetHadronizerProducer: No particles generated after " << ntry << " tries.\n";
00227          edm::Exception except(edm::errors::EventCorruption, sstr.str());
00228          throw except;
00229       } else {
00230          HYEVNT();
00231          nsoft_    = hyfpar.nhyd;
00232          nsub_     = hyjpar.njet;
00233          nhard_    = hyfpar.npyt;
00234          ++ntry;
00235       }
00236    }
00237 
00238    if(hyjpar.nhsel < 3) nsub_++;
00239 
00240    // event information
00241    HepMC::GenEvent *evt = new HepMC::GenEvent();
00242 
00243    if(nhard_>0 || nsoft_>0) get_particles(evt); 
00244 
00245    evt->set_signal_process_id(pypars.msti[0]);      // type of the process
00246    evt->set_event_scale(pypars.pari[16]);           // Q^2
00247    add_heavy_ion_rec(evt);
00248 
00249    event().reset(evt);
00250    return true;
00251 }
00252 
00253 
00254 //_____________________________________________________________________  
00255 bool HydjetHadronizer::get_particles(HepMC::GenEvent *evt )
00256 {
00257    // Hard particles. The first nhard_ lines from hyjets array.                
00258    // Pythia/Pyquen sub-events (sub-collisions) for a given event              
00259    // Return T/F if success/failure
00260    // Create particles from lujet entries, assign them into vertices and
00261    // put the vertices in the GenEvent, for each SubEvent
00262    // The SubEvent information is kept by storing indeces of main vertices
00263    // of subevents as a vector in GenHIEvent.
00264 
00265    LogDebug("SubEvent")<< "Number of sub events "<<nsub_;
00266    LogDebug("Hydjet")<<"Number of hard events "<<hyjpar.njet;
00267    LogDebug("Hydjet")<<"Number of hard particles "<<nhard_;
00268    LogDebug("Hydjet")<<"Number of soft particles "<<nsoft_;
00269 
00270    vector<HepMC::GenVertex*>  sub_vertices(nsub_);
00271 
00272    int ihy  = 0;
00273    for(int isub=0;isub<nsub_;isub++){
00274       LogDebug("SubEvent") <<"Sub Event ID : "<<isub;
00275 
00276       int sub_up = (isub+1)*50000; // Upper limit in mother index, determining the range of Sub-Event
00277       vector<HepMC::GenParticle*> particles;
00278       vector<int>                 mother_ids;
00279       vector<HepMC::GenVertex*>   prods;
00280 
00281       sub_vertices[isub] = new HepMC::GenVertex(HepMC::FourVector(0,0,0,0),isub);
00282       evt->add_vertex(sub_vertices[isub]);
00283       if(!evt->signal_process_vertex()) evt->set_signal_process_vertex(sub_vertices[isub]);
00284 
00285       while(ihy<nhard_+nsoft_ && (hyjets.khj[2][ihy] < sub_up || ihy > nhard_ )){
00286          particles.push_back(build_hyjet(ihy,ihy+1));
00287          prods.push_back(build_hyjet_vertex(ihy,isub));
00288          mother_ids.push_back(hyjets.khj[2][ihy]);
00289          LogDebug("DecayChain")<<"Mother index : "<<hyjets.khj[2][ihy];
00290 
00291          ihy++;
00292       }
00293 
00294       //Produce Vertices and add them to the GenEvent. Remember that GenParticles are adopted by
00295       //GenVertex and GenVertex is adopted by GenEvent.
00296 
00297       LogDebug("Hydjet")<<"Number of particles in vector "<<particles.size();
00298 
00299       for (unsigned int i = 0; i<particles.size(); i++) {
00300          HepMC::GenParticle* part = particles[i];
00301 
00302          //The Fortran code is modified to preserve mother id info, by seperating the beginning
00303          //mother indices of successive subevents by 5000
00304          int mid = mother_ids[i]-isub*50000-1;
00305          LogDebug("DecayChain")<<"Particle "<<i;
00306          LogDebug("DecayChain")<<"Mother's ID "<<mid;
00307          LogDebug("DecayChain")<<"Particle's PDG ID "<<part->pdg_id();
00308 
00309          if(mid <= 0){
00310             sub_vertices[isub]->add_particle_out(part);
00311             continue;
00312          }
00313 
00314          if(mid > 0){
00315             HepMC::GenParticle* mother = particles[mid];
00316             LogDebug("DecayChain")<<"Mother's PDG ID "<<mother->pdg_id();
00317 
00318             HepMC::GenVertex* prod_vertex = mother->end_vertex();
00319             if(!prod_vertex){
00320                prod_vertex = prods[i];
00321                prod_vertex->add_particle_in(mother);
00322                evt->add_vertex(prod_vertex);
00323                prods[i]=0; // mark to protect deletion
00324             }
00325             prod_vertex->add_particle_out(part);
00326          }
00327       }
00328       // cleanup vertices not assigned to evt
00329       for (unsigned int i = 0; i<prods.size(); i++) {
00330          if(prods[i]) delete prods[i];
00331       }
00332    }
00333    return true;
00334 }
00335 
00336 
00337 //______________________________________________________________
00338 bool HydjetHadronizer::call_hyinit(double energy,double a, int ifb, double bmin,
00339                                    double bmax,double bfix,int nh)
00340 {
00341   // initialize hydjet  
00342 
00343    pydatr.mrpy[2]=1;
00344    HYINIT(energy,a,ifb,bmin,bmax,bfix,nh);
00345    return true;
00346 }
00347 
00348 
00349 //______________________________________________________________
00350 bool HydjetHadronizer::hydjet_init(const ParameterSet &pset)
00351 {
00352   // set hydjet options
00353 
00354   // hydjet running mode mode
00355   // kHydroOnly --- nhsel=0 jet production off (pure HYDRO event), nhsel=0
00356   // kHydroJets --- nhsle=1 jet production on, jet quenching off (HYDRO+njet*PYTHIA events)
00357   // kHydroQJet --- nhsel=2 jet production & jet quenching on (HYDRO+njet*PYQUEN events)
00358   // kJetsOnly  --- nhsel=3 jet production on, jet quenching off, HYDRO off (njet*PYTHIA events)
00359   // kQJetsOnly --- nhsel=4 jet production & jet quenching on, HYDRO off (njet*PYQUEN events)
00360 
00361   if(hymode_ == "kHydroOnly") hyjpar.nhsel=0;
00362   else if ( hymode_ == "kHydroJets") hyjpar.nhsel=1;
00363   else if ( hymode_ == "kHydroQJets") hyjpar.nhsel=2;
00364   else if ( hymode_ == "kJetsOnly") hyjpar.nhsel=3;
00365   else if ( hymode_ == "kQJetsOnly") hyjpar.nhsel=4;
00366   else  hyjpar.nhsel=2;
00367 
00368   // fraction of soft hydro induced multiplicity 
00369   hyflow.fpart =  fracsoftmult_; 
00370 
00371   // hadron freez-out temperature
00372   hyflow.Tf   = hadfreeztemp_;
00373 
00374   // maximum longitudinal collective rapidity
00375   hyflow.ylfl = maxlongy_;
00376   
00377   // maximum transverse collective rapidity
00378   hyflow.ytfl = maxtrany_;  
00379 
00380   // shadowing on=1, off=0
00381   hyjpar.ishad  = shadowingswitch_;
00382 
00383   // set inelastic nucleon-nucleon cross section
00384   hyjpar.sigin  = signn_;
00385 
00386   // number of active quark flavors in qgp
00387   pyqpar.nfu    = nquarkflavor_;
00388 
00389   // initial temperature of QGP
00390   pyqpar.T0u    = qgpt0_;
00391 
00392   // proper time of QGP formation
00393   pyqpar.tau0u  = qgptau0_;
00394 
00395   // type of medium induced partonic energy loss
00396   if( doradiativeenloss_ && docollisionalenloss_ ){
00397     edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
00398     pyqpar.ienglu = 0; 
00399   } else if ( doradiativeenloss_ ) {
00400     edm::LogInfo("HydjetenLoss") << "##### Only RADIATIVE partonic energy loss ON ####";
00401     pyqpar.ienglu = 1; 
00402   } else if ( docollisionalenloss_ ) {
00403     edm::LogInfo("HydjetEnLoss") << "##### Only COLLISIONAL partonic energy loss ON ####";
00404     pyqpar.ienglu = 2; 
00405   } else {
00406     edm::LogInfo("HydjetEnLoss") << "##### Radiative AND Collisional partonic energy loss ON ####";
00407     pyqpar.ienglu = 0; 
00408   }
00409   return true;
00410 }
00411 
00412 //_____________________________________________________________________
00413 
00414 bool HydjetHadronizer::readSettings( int ) {
00415 
00416    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00417    pythia6Service_->setGeneralParams();
00418 
00419    return true;
00420 
00421 }
00422 
00423 //_____________________________________________________________________
00424 
00425 bool HydjetHadronizer::initializeForInternalPartons(){
00426 
00427    Pythia6Service::InstanceWrapper guard(pythia6Service_);
00428    // pythia6Service_->setGeneralParams();
00429 
00430    // the input impact parameter (bxx_) is in [fm]; transform in [fm/RA] for hydjet usage
00431    const float ra = nuclear_radius();
00432    LogInfo("RAScaling")<<"Nuclear radius(RA) =  "<<ra;
00433    bmin_     /= ra;
00434    bmax_     /= ra;
00435    bfixed_   /= ra;
00436 
00437    // hydjet running options 
00438    hydjet_init(pset_);
00439    // initialize hydjet
00440    LogInfo("HYDJETinAction") << "##### Calling HYINIT("<<comenergy<<","<<abeamtarget_<<","
00441                              <<cflag_<<","<<bmin_<<","<<bmax_<<","<<bfixed_<<","<<nmultiplicity_<<") ####";
00442    call_hyinit(comenergy,abeamtarget_,cflag_,bmin_,bmax_,bfixed_,nmultiplicity_);
00443    return true;
00444 }
00445 
00446 bool HydjetHadronizer::declareStableParticles( std::vector<int> pdg )
00447 {
00448   for ( size_t i=0; i < pdg.size(); i++ ) {
00449     int pyCode = pycomp_( pdg[i] );
00450     std::ostringstream pyCard ;
00451     pyCard << "MDCY(" << pyCode << ",1)=0";
00452     std::cout << pyCard.str() << std::endl;
00453     call_pygive( pyCard.str() );
00454   }
00455   return true;
00456 }
00457 
00458 //________________________________________________________________
00459 void HydjetHadronizer::rotateEvtPlane()
00460 {
00461   const double pi = 3.14159265358979;
00462   phi0_ = 2.*pi*gen::pyr_(0) - pi;
00463   sinphi0_ = sin(phi0_);
00464   cosphi0_ = cos(phi0_);
00465 }
00466 
00467 
00468 //________________________________________________________________
00469 bool HydjetHadronizer::hadronize()
00470 {
00471   return false;
00472 }
00473 
00474 bool HydjetHadronizer::decay()
00475 {
00476   return true;
00477 }
00478 
00479 bool HydjetHadronizer::residualDecay()
00480 {
00481   return true;
00482 }
00483 
00484 void HydjetHadronizer::finalizeEvent()
00485 {
00486 }
00487 
00488 void HydjetHadronizer::statistics()
00489 {
00490 }
00491 
00492 const char* HydjetHadronizer::classname() const
00493 {
00494    return "gen::HydjetHadronizer";
00495 }