CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/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/UserHooks.h"
00016 
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00020 
00021 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00022 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00023 
00024 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00025 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00026 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00027 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
00028 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00029 
00030 #include "GeneratorInterface/Pythia8Interface/interface/LHAupLesHouches.h"
00031 
00032 #include "HepPID/ParticleIDTranslations.hh"
00033 
00034 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00035 
00036 using namespace gen;
00037 using namespace Pythia8;
00038 
00039 class Pythia8Hadronizer : public BaseHadronizer {
00040     public:
00041         Pythia8Hadronizer(const edm::ParameterSet &params);
00042         ~Pythia8Hadronizer();
00043  
00044         bool initializeForInternalPartons();
00045         bool initializeForExternalPartons();
00046         
00047         bool declareStableParticles(const std::vector<int> &pdgIds);
00048         bool declareSpecialSettings( const std::vector<std::string> );
00049 
00050         void statistics();
00051 
00052         bool generatePartonsAndHadronize();
00053         bool hadronize();
00054         bool decay();
00055         bool residualDecay();
00056         void finalizeEvent();
00057 
00058         const char *classname() const { return "Pythia8Hadronizer"; }
00059 
00060     private:
00061         ParameterCollector      parameters;
00062 
00064         double                  comEnergy;
00066         unsigned int            pythiaPylistVerbosity;
00068         bool                    pythiaHepMCVerbosity;
00070         unsigned int            maxEventsToPrint;
00071 
00072     string LHEInputFileName;
00073 
00075     bool            useUserHook;
00076 
00077     std::auto_ptr<LHAupLesHouches>      lhaUP;
00078 
00079         std::auto_ptr<Pythia>   pythia;
00080         Event                   *pythiaEvent;
00081         HepMC::I_Pythia8        toHepMC;   
00082 
00083 };
00084 
00085 
00086 Pythia8Hadronizer::Pythia8Hadronizer(const edm::ParameterSet &params) :
00087         BaseHadronizer(params),
00088         parameters(params.getParameter<edm::ParameterSet>("PythiaParameters")),
00089         comEnergy(params.getParameter<double>("comEnergy")),
00090         pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
00091         pythiaHepMCVerbosity(params.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
00092         maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00093     LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
00094     useUserHook(false)
00095 {
00096     if( params.exists( "useUserHook" ) )
00097       useUserHook = params.getParameter<bool>("useUserHook");
00098     randomEngine = &getEngineReference();
00099 }
00100 
00101 Pythia8Hadronizer::~Pythia8Hadronizer()
00102 {
00103 }
00104 
00105 bool Pythia8Hadronizer::initializeForInternalPartons()
00106 {
00107         //Old code that used Pythia8 own random engine
00108         //edm::Service<edm::RandomNumberGenerator> rng;
00109         //uint32_t seed = rng->mySeed();
00110         //Pythia8::Rndm::init(seed);
00111 
00112     RandomP8* RP8 = new RandomP8();
00113 
00114         pythia.reset(new Pythia);
00115 
00116     pythia->setRndmEnginePtr(RP8);
00117     if(useUserHook) pythia->setUserHooksPtr(new PtHatReweightUserHook());
00118 
00119         pythiaEvent = &pythia->event;
00120 
00121         for(ParameterCollector::const_iterator line = parameters.begin();
00122             line != parameters.end(); ++line) {
00123                 if (line->find("Random:") != std::string::npos)
00124                         throw cms::Exception("PythiaError")
00125                                 << "Attempted to set random number "
00126                                    "using Pythia commands.  Please use "
00127                                    "the RandomNumberGeneratorService."
00128                                 << std::endl;
00129 
00130                 if (!pythia->readString(*line))
00131                         throw cms::Exception("PythiaError")
00132                                 << "Pythia 8 did not accept \""
00133                                 << *line << "\"." << std::endl;
00134         }
00135 
00136     if(pythiaPylistVerbosity > 10) {
00137       if(pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13)
00138         pythia->settings.listAll();
00139       if(pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13)
00140         pythia->particleData.listAll();
00141     }
00142 
00143         pythia->init(2212, 2212, comEnergy);
00144 
00145         pythia->settings.listChanged();
00146 
00147         return true;
00148 }
00149 
00150 
00151 bool Pythia8Hadronizer::initializeForExternalPartons()
00152 {
00153 
00154     std::cout << "Initializing for external partons" << std::endl;
00155 
00156     RandomP8* RP8 = new RandomP8();
00157 
00158     pythia.reset(new Pythia);
00159 
00160     pythia->setRndmEnginePtr(RP8);
00161 
00162     pythiaEvent = &pythia->event;
00163 
00164     for(ParameterCollector::const_iterator line = parameters.begin();
00165         line != parameters.end(); ++line) {
00166         if (line->find("Random:") != std::string::npos)
00167             throw cms::Exception("PythiaError")
00168                 << "Attempted to set random number "
00169                    "using Pythia commands.  Please use "
00170                    "the RandomNumberGeneratorService."
00171                 << std::endl;
00172 
00173         if (!pythia->readString(*line))
00174             throw cms::Exception("PythiaError")
00175                 << "Pythia 8 did not accept \""
00176                 << *line << "\"." << std::endl;
00177     }
00178 
00179     if(pythiaPylistVerbosity > 10) {
00180       if(pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13)
00181         pythia->settings.listAll();
00182       if(pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13)
00183         pythia->particleData.listAll();
00184     }
00185 
00186     if(LHEInputFileName != string()) {
00187 
00188       cout << endl;
00189       cout << "LHE Input File Name = " << LHEInputFileName << endl;
00190       cout << endl;
00191       pythia->init(LHEInputFileName);
00192 
00193     } else {
00194 
00195       lhaUP.reset(new LHAupLesHouches());
00196       lhaUP->loadRunInfo(lheRunInfo());
00197       pythia->init(lhaUP.get());
00198 
00199     }
00200 
00201     return true;
00202 }
00203 
00204 
00205 #if 0
00206 // naive Pythia8 HepMC status fixup
00207 static int getStatus(const HepMC::GenParticle *p)
00208 {
00209         int status = p->status();
00210         if (status > 0)
00211                 return status;
00212         else if (status > -30 && status < 0)
00213                 return 3;
00214         else
00215                 return 2;
00216 }
00217 #endif
00218 
00219 bool Pythia8Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00220 {
00221 
00222    for ( size_t i=0; i<pdgIds.size(); i++ )
00223    {
00224       // FIXME: need to double check if PID's are the same in Py6 & Py8,
00225       //        because the HepPDT translation tool is actually for **Py6** 
00226       // 
00227       int PyID = HepPID::translatePDTtoPythia( pdgIds[i] ); 
00228       std::ostringstream pyCard ;
00229       pyCard << PyID <<":mayDecay=false";
00230       pythia->readString( pyCard.str() );
00231       // alternative:
00232       // set the 2nd input argument warn=false 
00233       // - this way Py8 will NOT print warnings about unknown particle code(s)
00234       // pythia->readString( pyCard.str(), false )
00235    }
00236    
00237    return true;
00238 
00239 }
00240 bool Pythia8Hadronizer::declareSpecialSettings( const std::vector<std::string> )
00241 {
00242    return true;
00243 }
00244 
00245 
00246 void Pythia8Hadronizer::statistics()
00247 {
00248         pythia->statistics();
00249 
00250         double xsec = pythia->info.sigmaGen(); // cross section in mb
00251     xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
00252         runInfo().setInternalXSec(xsec);
00253 }
00254 
00255 bool Pythia8Hadronizer::generatePartonsAndHadronize()
00256 {
00257         if (!pythia->next())
00258                 return false;
00259 
00260         event().reset(new HepMC::GenEvent);
00261         toHepMC.fill_next_event(*pythiaEvent, event().get());
00262 
00263         return true;
00264 }
00265 
00266 bool Pythia8Hadronizer::hadronize()
00267 {
00268     if(LHEInputFileName == string()) {
00269       //cout << "start loading event" << endl;
00270       lhaUP->loadEvent(lheEvent());
00271       //cout << "finish loading event" << endl;
00272     }
00273 
00274     if (!pythia->next())
00275         return false;
00276 
00277     // update LHE matching statistics
00278     //
00279     lheEvent()->count( lhef::LHERunInfo::kAccepted );
00280 
00281     event().reset(new HepMC::GenEvent);
00282     toHepMC.fill_next_event(*pythiaEvent, event().get());
00283 
00284     return true;
00285 }
00286 
00287 bool Pythia8Hadronizer::decay()
00288 {
00289         return true;
00290 }
00291 
00292 bool Pythia8Hadronizer::residualDecay()
00293 {
00294         return true;
00295 }
00296 
00297 void Pythia8Hadronizer::finalizeEvent()
00298 {
00299   bool lhe = lheEvent() != 0;
00300 
00301   event()->set_signal_process_id(pythia->info.code());
00302   event()->set_event_scale(pythia->info.pTHat());       //FIXME
00303 
00304   int id1 = pythia->info.id1();
00305   int id2 = pythia->info.id2();
00306   if (id1 == 21) id1 = 0;
00307   if (id2 == 21) id2 = 0;
00308   double x1 = pythia->info.x1();
00309   double x2 = pythia->info.x2();
00310   double Q = pythia->info.QRen();
00311   double pdf1 = pythia->info.pdf1() / pythia->info.x1();
00312   double pdf2 = pythia->info.pdf2() / pythia->info.x2();
00313   event()->set_pdf_info(HepMC::PdfInfo(id1,id2,x1,x2,Q,pdf1,pdf2));
00314 
00315   event()->weights().push_back(pythia->info.weight());
00316 
00317   // now create the GenEventInfo product from the GenEvent and fill
00318   // the missing pieces
00319   eventInfo().reset( new GenEventInfoProduct( event().get() ) );
00320 
00321   // in pythia pthat is used to subdivide samples into different bins
00322   // in LHE mode the binning is done by the external ME generator
00323   // which is likely not pthat, so only filling it for Py6 internal mode
00324   if (!lhe) {
00325     eventInfo()->setBinningValues(std::vector<double>(1, pythia->info.pTHat()));
00326   }
00327 
00328   //******** Verbosity ********
00329 
00330   if (maxEventsToPrint > 0 &&
00331       (pythiaPylistVerbosity || pythiaHepMCVerbosity)) {
00332     maxEventsToPrint--;
00333     if (pythiaPylistVerbosity) {
00334       pythia->info.list(std::cout); 
00335       pythia->event.list(std::cout);
00336     } 
00337 
00338         if (pythiaHepMCVerbosity) {
00339                         std::cout << "Event process = "
00340                                   << pythia->info.code() << "\n"
00341                                   << "----------------------" << std::endl;
00342                         event()->print();
00343                 }
00344         }
00345 }
00346 
00347 typedef edm::GeneratorFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8GeneratorFilter;
00348 DEFINE_FWK_MODULE(Pythia8GeneratorFilter);
00349 
00350 typedef edm::HadronizerFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8HadronizerFilter;
00351 DEFINE_FWK_MODULE(Pythia8HadronizerFilter);