CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/GeneratorInterface/Pythia8Interface/src/Pythia8Hadronizer.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <string>
00004 #include <memory>
00005 #include <stdint.h>
00006 
00007 #include <HepMC/GenEvent.h>
00008 #include <HepMC/GenParticle.h>
00009 
00010 #include <Pythia.h>
00011 #include <HepMCInterface.h>
00012 
00013 #include "GeneratorInterface/Pythia8Interface/interface/RandomP8.h"
00014 
00015 #include "GeneratorInterface/Pythia8Interface/interface/ReweightUserHooks.h"
00016 
00017 // PS matchning prototype
00018 //
00019 #include "GeneratorInterface/Pythia8Interface/interface/JetMatchingHook.h"
00020 
00021 // Emission Veto Hook
00022 //
00023 #include "GeneratorInterface/Pythia8Interface/interface/EmissionVetoHook.h"
00024 
00025 #include "FWCore/ServiceRegistry/interface/Service.h"
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00028 
00029 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00030 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00031 
00032 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00033 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00034 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00035 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
00036 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00037 
00038 #include "GeneratorInterface/Pythia8Interface/interface/LHAupLesHouches.h"
00039 
00040 #include "HepPID/ParticleIDTranslations.hh"
00041 
00042 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00043 
00044 using namespace gen;
00045 using namespace Pythia8;
00046 
00047 class Pythia8Hadronizer : public BaseHadronizer {
00048   public:
00049     Pythia8Hadronizer(const edm::ParameterSet &params);
00050    ~Pythia8Hadronizer();
00051  
00052     bool readSettings( int );
00053     bool initializeForInternalPartons();
00054     bool initializeForExternalPartons();
00055         
00056     bool declareStableParticles(const std::vector<int> &pdgIds);
00057     bool declareSpecialSettings( const std::vector<std::string> );
00058 
00059     void statistics();
00060 
00061     bool generatePartonsAndHadronize();
00062     bool hadronize();
00063     bool decay();
00064     bool residualDecay();
00065     void finalizeEvent();
00066 
00067     const char *classname() const { return "Pythia8Hadronizer"; }
00068 
00069   private:
00070     ParameterCollector  parameters;
00071 
00073     double       comEnergy;
00075     unsigned int pythiaPylistVerbosity;
00077     bool         pythiaHepMCVerbosity;
00079     unsigned int maxEventsToPrint;
00080 
00081     string LHEInputFileName;
00082 
00083     std::auto_ptr<LHAupLesHouches>  lhaUP;
00084 
00085     std::auto_ptr<Pythia>   pythia;
00086     std::auto_ptr<Pythia>   decayer;
00087     Event*                  pythiaEvent;
00088     HepMC::I_Pythia8        toHepMC;
00089 
00090     enum { PP, PPbar, ElectronPositron };
00091     int  fInitialState ; // pp, ppbar, or e-e+
00092 
00093     double fBeam1PZ;
00094     double fBeam2PZ;
00095 
00096     // Reweight user hook
00097     //
00098     UserHooks* fReweightUserHook;
00099         
00100     // PS matching protot6ype
00101     //
00102     JetMatchingHook* fJetMatchingHook;
00103         
00104     // Emission Veto Hook
00105     //
00106     EmissionVetoHook* fEmissionVetoHook;
00107 
00108 };
00109 
00110 
00111 Pythia8Hadronizer::Pythia8Hadronizer(const edm::ParameterSet &params) :
00112   BaseHadronizer(params),
00113   parameters(params.getParameter<edm::ParameterSet>("PythiaParameters")),
00114   comEnergy(params.getParameter<double>("comEnergy")),
00115   pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
00116   pythiaHepMCVerbosity(params.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
00117   maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00118   LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
00119   fInitialState(PP),
00120   fReweightUserHook(0),
00121   fJetMatchingHook(0),
00122   fEmissionVetoHook(0)
00123 {
00124   randomEngine = &getEngineReference();
00125 
00126   //Old code that used Pythia8 own random engine
00127   //edm::Service<edm::RandomNumberGenerator> rng;
00128   //uint32_t seed = rng->mySeed();
00129   //Pythia8::Rndm::init(seed);
00130 
00131   RandomP8* RP8 = new RandomP8();
00132 
00133   // J.Y.: the following 3 parameters are hacked "for a reason"
00134   //
00135   if ( params.exists( "PPbarInitialState" ) )
00136   {
00137     if ( fInitialState == PP )
00138     {
00139       fInitialState = PPbar;
00140       edm::LogInfo("GeneratorInterface|Pythia6Interface")
00141       << "Pythia6 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
00142       << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00143       std::cout << "Pythia6 will be initialized for PROTON-ANTIPROTON INITIAL STATE." << std::endl;
00144       std::cout << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00145     }
00146     else
00147     {   
00148       // probably need to throw on attempt to override ?
00149     }
00150   }   
00151   else if ( params.exists( "ElectronPositronInitialState" ) )
00152   {
00153     if ( fInitialState == PP )
00154     {
00155       fInitialState = ElectronPositron;
00156       edm::LogInfo("GeneratorInterface|Pythia6Interface")
00157       << "Pythia6 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
00158       << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00159       std::cout << "Pythia6 will be initialized for ELECTRON-POSITRON INITIAL STATE." << std::endl; 
00160       std::cout << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00161     }
00162     else
00163     {   
00164        // probably need to throw on attempt to override ?
00165     }
00166   }
00167   else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
00168   {
00169     // throw on unknown initial state !
00170     throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00171       <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
00172   }
00173 
00174   pythia.reset(new Pythia);
00175   decayer.reset(new Pythia);
00176 
00177   pythia->setRndmEnginePtr(RP8);
00178   decayer->setRndmEnginePtr(RP8);
00179     
00180   if( params.exists( "SLHAFileForPythia8" ) ) {
00181     std::string slhafilenameshort = params.getParameter<string>("SLHAFileForPythia8");
00182     edm::FileInPath f1( slhafilenameshort );
00183     std::string slhafilename = f1.fullPath();
00184     std::string pythiacommandslha = std::string("SLHA:file = ") + slhafilename;
00185     pythia->readString(pythiacommandslha);
00186     for ( ParameterCollector::const_iterator line = parameters.begin();
00187           line != parameters.end(); ++line ) {
00188       if (line->find("SLHA:file") != std::string::npos)
00189         throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
00190         << "using Pythia8 card SLHA:file and Pythia8Interface card SLHAFileForPythia8"
00191         << std::endl;
00192      }
00193   } 
00194 
00195   // Reweight user hook
00196   //
00197   if( params.exists( "reweightGen" ) )
00198     fReweightUserHook = new PtHatReweightUserHook();
00199 
00200   // PS matching prototype
00201   //
00202   if ( params.exists("jetMatching") )
00203   {
00204     edm::ParameterSet jmParams =
00205       params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
00206     fJetMatchingHook = new JetMatchingHook( jmParams, &pythia->info );
00207   }
00208 
00209   // Emission veto
00210   //
00211   if ( params.exists("emissionVeto") )
00212   {   
00213     fEmissionVetoHook = new EmissionVetoHook(0);
00214     pythia->setUserHooksPtr( fEmissionVetoHook );
00215   }  
00216 
00217   int NHooks=0;
00218   if(fReweightUserHook) NHooks++;
00219   if(fJetMatchingHook) NHooks++;
00220   if(fEmissionVetoHook) NHooks++;
00221   if(NHooks > 1)
00222     throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00223       <<" Too many User Hooks. \n Please choose one from: reweightGen, jetMatching, emissionVeto \n";
00224 
00225   if(fReweightUserHook) pythia->setUserHooksPtr(fReweightUserHook);
00226   if(fJetMatchingHook) pythia->setUserHooksPtr(fJetMatchingHook);
00227   if(fEmissionVetoHook) pythia->setUserHooksPtr(fEmissionVetoHook);
00228 }
00229 
00230 
00231 Pythia8Hadronizer::~Pythia8Hadronizer()
00232 {
00233 // do we need to delete UserHooks/JetMatchingHook here ???
00234 
00235   if(fEmissionVetoHook) {delete fEmissionVetoHook; fEmissionVetoHook=0;}
00236 }
00237 
00238 
00239 bool Pythia8Hadronizer::readSettings( int )
00240 {
00241   for ( ParameterCollector::const_iterator line = parameters.begin();
00242         line != parameters.end(); ++line ) {
00243     if (line->find("Random:") != std::string::npos)
00244       throw cms::Exception("PythiaError") << "Attempted to set random number "
00245         "using Pythia commands. Please use " "the RandomNumberGeneratorService."
00246         << std::endl;
00247 
00248     if (!pythia->readString(*line)) throw cms::Exception("PythiaError")
00249                                       << "Pythia 8 did not accept \""
00250                                       << *line << "\"." << std::endl;
00251   }
00252 
00253   if ( pythiaPylistVerbosity > 10 ) {
00254     if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
00255       pythia->settings.listAll();
00256     if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
00257       pythia->particleData.listAll();
00258   }
00259 
00260   return true;
00261 }
00262 
00263 
00264 bool Pythia8Hadronizer::initializeForInternalPartons()
00265 {
00266 
00267   pythiaEvent = &pythia->event;
00268 
00269   if ( fInitialState == PP ) // default
00270   {
00271     pythia->init(2212, 2212, comEnergy);
00272   }
00273   else if ( fInitialState == PPbar )
00274   {
00275     pythia->init(2212, -2212, comEnergy);
00276   }
00277   else if ( fInitialState == ElectronPositron )
00278   {
00279     pythia->init(11, -11, comEnergy);
00280   }    
00281   else 
00282   {
00283     // throw on unknown initial state !
00284     throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00285       <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
00286   }
00287 
00288   pythia->settings.listChanged();
00289 
00290   return true;
00291 }
00292 
00293 
00294 bool Pythia8Hadronizer::initializeForExternalPartons()
00295 {
00296 
00297   std::cout << "Initializing for external partons" << std::endl;
00298 
00299   pythiaEvent = &pythia->event;
00300     
00301   if(LHEInputFileName != string()) {
00302 
00303     cout << endl;
00304     cout << "LHE Input File Name = " << LHEInputFileName << endl;
00305     cout << endl;
00306     pythia->init(LHEInputFileName);
00307 
00308   } else {
00309 
00310     lhaUP.reset(new LHAupLesHouches());
00311     lhaUP->loadRunInfo(lheRunInfo());
00312     pythia->init(lhaUP.get());
00313 
00314   }
00315 
00316   // PS matching prototype
00317   //
00318   if ( fJetMatchingHook ) 
00319   {
00320     // matcher will be init as well, inside init(...)
00321     //
00322     fJetMatchingHook->init ( lheRunInfo() );
00323   }
00324 
00325     return true;
00326 }
00327 
00328 
00329 #if 0
00330 // naive Pythia8 HepMC status fixup
00331 static int getStatus(const HepMC::GenParticle *p)
00332 {
00333   int status = p->status();
00334   if (status > 0)
00335     return status;
00336   else if (status > -30 && status < 0)
00337     return 3;
00338   else
00339     return 2;
00340 }
00341 #endif
00342 
00343 
00344 bool Pythia8Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00345 {
00346   for ( size_t i=0; i<pdgIds.size(); i++ )
00347   {
00348     // FIXME: need to double check if PID's are the same in Py6 & Py8,
00349     //        because the HepPDT translation tool is actually for **Py6** 
00350     // 
00351     // well, actually it looks like Py8 operates in PDT id's rather than Py6's
00352     //
00353     // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] ); 
00354     int PyID = pdgIds[i]; 
00355     std::ostringstream pyCard ;
00356     pyCard << PyID <<":mayDecay=false";
00357     pythia->readString( pyCard.str() );
00358     // alternative:
00359     // set the 2nd input argument warn=false 
00360     // - this way Py8 will NOT print warnings about unknown particle code(s)
00361     // pythia->readString( pyCard.str(), false )
00362   }
00363       
00364   // init decayer
00365   decayer->readString("ProcessLevel:all = off"); // The trick!
00366   decayer->init();
00367    
00368   return true;
00369 }
00370 
00371 
00372 bool Pythia8Hadronizer::declareSpecialSettings( const std::vector<std::string> settings )
00373 {
00374   for ( unsigned int iss=0; iss<settings.size(); iss++ )
00375   {
00376     if ( settings[iss].find("QED-brem-off") == std::string::npos ) continue;
00377     pythia->readString( "TimeShower:QEDshowerByL=off" );
00378   }
00379 
00380   return true;
00381 }
00382 
00383 
00384 void Pythia8Hadronizer::statistics()
00385 {
00386   pythia->statistics();
00387 
00388   double xsec = pythia->info.sigmaGen(); // cross section in mb
00389   xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
00390   runInfo().setInternalXSec(xsec);
00391 }
00392 
00393 
00394 bool Pythia8Hadronizer::generatePartonsAndHadronize()
00395 {
00396   if (!pythia->next()) return false;
00397 
00398   event().reset(new HepMC::GenEvent);
00399   toHepMC.fill_next_event(*pythiaEvent, event().get());
00400 
00401   return true;
00402 }
00403 
00404 
00405 bool Pythia8Hadronizer::hadronize()
00406 {
00407   if(LHEInputFileName == string()) lhaUP->loadEvent(lheEvent());
00408 
00409   if ( fJetMatchingHook ) 
00410   {
00411     fJetMatchingHook->resetMatchingStatus(); 
00412     fJetMatchingHook->beforeHadronization( lheEvent() );
00413   }
00414 
00415   bool py8next = pythia->next();
00416   // if (!pythia->next())
00417   if (!py8next)
00418   {
00419     lheEvent()->count( lhef::LHERunInfo::kSelected );
00420     event().reset();
00421     return false;
00422   }
00423 
00424   // update LHE matching statistics
00425   //
00426   lheEvent()->count( lhef::LHERunInfo::kAccepted );
00427 
00428   event().reset(new HepMC::GenEvent);
00429   toHepMC.fill_next_event(*pythiaEvent, event().get());
00430 
00431   return true;
00432 }
00433 
00434 
00435 bool Pythia8Hadronizer::decay()
00436 {
00437    return true;
00438 }
00439 
00440 
00441 bool Pythia8Hadronizer::residualDecay()
00442 {
00443 
00444   int NPartsBeforeDecays = pythiaEvent->size();
00445   int NPartsAfterDecays = event().get()->particles_size();
00446   int NewBarcode = NPartsAfterDecays;
00447    
00448   for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
00449   {
00450 
00451     HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
00452 
00453     if ( part->status() == 1 )
00454     {
00455       decayer->event.reset();
00456       Particle py8part(  part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
00457                          part->momentum().x(),
00458                          part->momentum().y(),
00459                          part->momentum().z(),
00460                          part->momentum().t(),
00461                          part->generated_mass() );
00462       HepMC::GenVertex* ProdVtx = part->production_vertex();
00463       py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(), 
00464                      ProdVtx->position().z(), ProdVtx->position().t() );
00465       py8part.tau( (decayer->particleData).tau0( part->pdg_id() ) );
00466       decayer->event.append( py8part );
00467       int nentries = decayer->event.size();
00468       if ( !decayer->event[nentries-1].mayDecay() ) continue;
00469       decayer->next();
00470       int nentries1 = decayer->event.size();
00471       // --> decayer->event.list(std::cout);
00472       if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
00473             
00474       part->set_status(2);
00475             
00476       Particle& py8daughter = decayer->event[nentries]; // the 1st daughter
00477       HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
00478                                                        py8daughter.yProd(),
00479                                                        py8daughter.zProd(),
00480                                                        py8daughter.tProd()) );
00481 
00482       DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
00483                                        // I presume (vtx) barcode will be given automatically
00484             
00485       HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
00486             
00487       HepMC::GenParticle* daughter =
00488                         new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
00489             
00490       NewBarcode++;
00491       daughter->suggest_barcode( NewBarcode );
00492       DecVtx->add_particle_out( daughter );
00493                     
00494       for ( ipart=nentries+1; ipart<nentries1; ipart++ )
00495       {
00496         py8daughter = decayer->event[ipart];
00497         HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );           
00498         HepMC::GenParticle* daughterN =
00499                         new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
00500         NewBarcode++;
00501         daughterN->suggest_barcode( NewBarcode );
00502         DecVtx->add_particle_out( daughterN );
00503       }
00504             
00505       event().get()->add_vertex( DecVtx );
00506 
00507     }
00508  } 
00509    
00510  return true;
00511 }
00512 
00513 
00514 void Pythia8Hadronizer::finalizeEvent()
00515 {
00516   bool lhe = lheEvent() != 0;
00517 
00518   event()->set_signal_process_id(pythia->info.code());
00519   event()->set_event_scale(pythia->info.pTHat());       //FIXME
00520 
00521   //cout.precision(10);
00522   //cout << " pt = " << pythia->info.pTHat() << " weights = "
00523   //     << pythia->info.weight() << " "
00524   //     << fReweightUserHook->biasedSelectionWeight() << endl;
00525 
00526   if (event()->alphaQED() <= 0)
00527     event()->set_alphaQED( pythia->info.alphaEM() );
00528   if (event()->alphaQCD() <= 0)
00529     event()->set_alphaQCD( pythia->info.alphaS() );
00530 
00531   HepMC::GenCrossSection xsec;
00532   xsec.set_cross_section( pythia->info.sigmaGen() * 1e9,
00533                           pythia->info.sigmaErr() * 1e9);
00534   event()->set_cross_section(xsec);
00535 
00536   // Putting pdf info into the HepMC record
00537   // There is the overloaded pythia8 HepMCInterface method fill_next_event
00538   // that does this, but CMSSW GeneratorInterface does not fill HepMC
00539   // record according to the HepMC convention (stores f(x) instead of x*f(x)
00540   // and converts gluon PDG ID to zero). For this reason we use the
00541   // method fill_next_event (above) that does NOT this and fill pdf info here
00542   //
00543   int id1 = pythia->info.id1();
00544   int id2 = pythia->info.id2();
00545   if (id1 == 21) id1 = 0;
00546   if (id2 == 21) id2 = 0;
00547   double x1 = pythia->info.x1();
00548   double x2 = pythia->info.x2();
00549   //double Q = pythia->info.QRen();
00550   double Q = pythia->info.QFac();
00551   double pdf1 = pythia->info.pdf1() / pythia->info.x1();
00552   double pdf2 = pythia->info.pdf2() / pythia->info.x2();
00553   event()->set_pdf_info(HepMC::PdfInfo(id1,id2,x1,x2,Q,pdf1,pdf2));
00554 
00555   // Storing weights. Will be moved to pythia8 HepMCInterface
00556   //
00557   if (lhe && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
00558     // translate mb to pb (CMS/Gen "convention" as of May 2009)
00559     event()->weights().push_back( pythia->info.weight() * 1.0e9 );
00560   else
00561     event()->weights().push_back( pythia->info.weight() );
00562 
00563   // now create the GenEventInfo product from the GenEvent and fill
00564   // the missing pieces
00565   eventInfo().reset( new GenEventInfoProduct( event().get() ) );
00566 
00567   // in pythia pthat is used to subdivide samples into different bins
00568   // in LHE mode the binning is done by the external ME generator
00569   // which is likely not pthat, so only filling it for Py6 internal mode
00570   if (!lhe) {
00571     eventInfo()->setBinningValues(std::vector<double>(1, pythia->info.pTHat()));
00572   }
00573 
00574   //******** Verbosity ********
00575 
00576   if (maxEventsToPrint > 0 &&
00577       (pythiaPylistVerbosity || pythiaHepMCVerbosity)) {
00578     maxEventsToPrint--;
00579     if (pythiaPylistVerbosity) {
00580       pythia->info.list(std::cout); 
00581       pythia->event.list(std::cout);
00582     } 
00583 
00584     if (pythiaHepMCVerbosity) {
00585       std::cout << "Event process = "
00586                 << pythia->info.code() << "\n"
00587                 << "----------------------" << std::endl;
00588       event()->print();
00589     }
00590   }
00591 }
00592 
00593 
00594 typedef edm::GeneratorFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8GeneratorFilter;
00595 DEFINE_FWK_MODULE(Pythia8GeneratorFilter);
00596 
00597 
00598 typedef edm::HadronizerFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8HadronizerFilter;
00599 DEFINE_FWK_MODULE(Pythia8HadronizerFilter);