CMS 3D CMS Logo

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