CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/GeneratorInterface/LHEInterface/plugins/Pythia8Hadronisation.cc

Go to the documentation of this file.
00001 #include <algorithm>
00002 #include <iostream>
00003 #include <iterator>
00004 #include <sstream>
00005 #include <string>
00006 #include <memory>
00007 #include <assert.h>
00008 
00009 #include <boost/shared_ptr.hpp>
00010 
00011 #include <HepMC/GenEvent.h>
00012 #include <HepMC/GenParticle.h>
00013 
00014 #include <Pythia.h>
00015 #include <LesHouches.h>
00016 #include <HepMCInterface.h>
00017 
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "FWCore/ServiceRegistry/interface/Service.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00022 
00023 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
00024 #include "GeneratorInterface/LHEInterface/interface/LHERunInfo.h"
00025 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
00026 #include "GeneratorInterface/LHEInterface/interface/Hadronisation.h"
00027 
00028 using namespace Pythia8;
00029 
00030 namespace lhef {
00031 
00032 class Pythia8Hadronisation : public Hadronisation {
00033     public:
00034         Pythia8Hadronisation(const edm::ParameterSet &params);
00035         ~Pythia8Hadronisation();
00036 
00037     private:
00038         void doInit();
00039         std::auto_ptr<HepMC::GenEvent> doHadronisation();
00040         void newRunInfo(const boost::shared_ptr<LHERunInfo> &runInfo);
00041 
00042         const int                               pythiaPylistVerbosity;
00043         int                                     maxEventsToPrint;
00044         std::vector<std::string>                paramLines;
00045 
00046         class LHAupLesHouches;
00047 
00048         std::auto_ptr<Pythia>                   pythia;
00049         std::auto_ptr<LHAupLesHouches>          lhaUP;
00050         std::auto_ptr<HepMC::I_Pythia8>         conv;
00051 };
00052 
00053 class Pythia8Hadronisation::LHAupLesHouches : public LHAup {
00054     public:
00055         LHAupLesHouches(Hadronisation *hadronisation) :
00056                                         hadronisation(hadronisation) {}
00057 
00058         void loadRunInfo(const boost::shared_ptr<LHERunInfo> &runInfo)
00059         { this->runInfo = runInfo; }
00060 
00061         void loadEvent(const boost::shared_ptr<LHEEvent> &event)
00062         { this->event = event; }
00063 
00064     private:
00065 
00066         bool setInit();
00067         bool setEvent(int idProcIn);
00068 
00069         Hadronisation                   *hadronisation;
00070         boost::shared_ptr<LHERunInfo>   runInfo;
00071         boost::shared_ptr<LHEEvent>     event;
00072 };
00073 
00074 bool Pythia8Hadronisation::LHAupLesHouches::setInit()
00075 {
00076         if (!runInfo)
00077                 return false;
00078         const HEPRUP &heprup = *runInfo->getHEPRUP();
00079 
00080         setBeamA(heprup.IDBMUP.first, heprup.EBMUP.first,
00081                  heprup.PDFGUP.first, heprup.PDFSUP.first);
00082         setBeamB(heprup.IDBMUP.second, heprup.EBMUP.second,
00083                  heprup.PDFGUP.second, heprup.PDFSUP.second);
00084         setStrategy(heprup.IDWTUP);
00085 
00086         for(int i = 0; i < heprup.NPRUP; i++)
00087                 addProcess(heprup.LPRUP[i], heprup.XSECUP[i],
00088                            heprup.XERRUP[i], heprup.XMAXUP[i]);
00089 
00090         hadronisation->onInit().emit();
00091 
00092         runInfo.reset();
00093         return true;
00094 }
00095 
00096 bool Pythia8Hadronisation::LHAupLesHouches::setEvent(int inProcId)
00097 {
00098         if (!event)
00099                 return false;
00100         const HEPEUP &hepeup = *event->getHEPEUP();
00101 
00102         setProcess(hepeup.IDPRUP, hepeup.XWGTUP, hepeup.SCALUP,
00103                    hepeup.AQEDUP, hepeup.AQCDUP);
00104 
00105         for(int i = 0; i < hepeup.NUP; i++)
00106                 addParticle(hepeup.IDUP[i], hepeup.ISTUP[i],
00107                             hepeup.MOTHUP[i].first, hepeup.MOTHUP[i].second,
00108                             hepeup.ICOLUP[i].first, hepeup.ICOLUP[i].second,
00109                             hepeup.PUP[i][0], hepeup.PUP[i][1],
00110                             hepeup.PUP[i][2], hepeup.PUP[i][3],
00111                             hepeup.PUP[i][4], hepeup.VTIMUP[i],
00112                             hepeup.SPINUP[i]);
00113 
00114         const LHEEvent::PDF *pdf = event->getPDF();
00115         if (pdf)
00116                 this->setPdf(pdf->id.first, pdf->id.second,
00117                              pdf->x.first, pdf->x.second,
00118                              pdf->scalePDF,
00119                              pdf->xPDF.first, pdf->xPDF.second, true);
00120 
00121         hadronisation->onBeforeHadronisation().emit();
00122 
00123         event.reset();
00124         return true;
00125 }
00126 
00127 Pythia8Hadronisation::Pythia8Hadronisation(const edm::ParameterSet &params) :
00128         Hadronisation(params),
00129         pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
00130         maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0))
00131 {
00132         std::vector<std::string> setNames =
00133                 params.getParameter<std::vector<std::string> >("parameterSets");
00134 
00135         for(std::vector<std::string>::const_iterator iter = setNames.begin();
00136             iter != setNames.end(); ++iter) {
00137                 std::vector<std::string> lines =
00138                         params.getParameter< std::vector<std::string> >(*iter);
00139 
00140                 for(std::vector<std::string>::const_iterator line = lines.begin();
00141                     line != lines.end(); ++line )
00142                         if (line->substr(0, 14) == "Random:setSeed" ||
00143                             line->substr(0, 11) == "Random:seed")
00144                                 throw cms::Exception("PythiaError")
00145                                         << "Attempted to set random number"
00146                                            " using Pythia command 'MRPY(1)'."
00147                                            " Please use the"
00148                                            " RandomNumberGeneratorService."
00149                                         << std::endl;
00150 
00151                 std::copy(lines.begin(), lines.end(),
00152                           std::back_inserter(paramLines));
00153         }
00154 }
00155 
00156 Pythia8Hadronisation::~Pythia8Hadronisation()
00157 {
00158 }
00159 
00160 void Pythia8Hadronisation::doInit()
00161 {
00162         pythia.reset(new Pythia);
00163         lhaUP.reset(new LHAupLesHouches(this));
00164         conv.reset(new HepMC::I_Pythia8);
00165 
00166         for(std::vector<std::string>::const_iterator iter = paramLines.begin();
00167             iter != paramLines.end(); ++iter)
00168                 if (!pythia->readString(*iter))
00169                         throw cms::Exception("PythiaError")
00170                                 << "Pythia did not accept \""
00171                                 << *iter << "\"." << std::endl;
00172 
00173         edm::Service<edm::RandomNumberGenerator> rng;
00174         std::ostringstream ss;
00175         ss << "Random:seed = " << rng->mySeed();
00176         pythia->readString(ss.str());
00177         pythia->readString("Random:setSeed = on");
00178 }
00179 
00180 // naive Pythia8 HepMC status fixup
00181 static int getStatus(const HepMC::GenParticle *p)
00182 {
00183         int status = p->status();
00184         if (status > 0)
00185                 return status;
00186         else if (status > -30 && status < 0)
00187                 return 3;
00188         else
00189                 return 2;
00190 }
00191 
00192 std::auto_ptr<HepMC::GenEvent> Pythia8Hadronisation::doHadronisation()
00193 {
00194         lhaUP->loadEvent(getRawEvent());
00195         if (!pythia->next())
00196                 throw cms::Exception("PythiaError")
00197                         << "Pythia did not want to process event."
00198                         << std::endl;
00199 
00200         std::auto_ptr<HepMC::GenEvent> event(new HepMC::GenEvent);
00201         conv->fill_next_event(pythia->event, event.get());
00202 
00203         for(HepMC::GenEvent::particle_iterator iter = event->particles_begin();
00204             iter != event->particles_end(); iter++)
00205                 (*iter)->set_status(getStatus(*iter));
00206 
00207         event->set_signal_process_id(pythia->info.code());
00208         event->set_event_scale(pythia->info.pTHat());
00209 
00210         if (maxEventsToPrint > 0) {
00211                 maxEventsToPrint--;
00212                 if (pythiaPylistVerbosity)
00213                         pythia->event.list(std::cout);
00214         }
00215 
00216         return event;
00217 }
00218 
00219 void Pythia8Hadronisation::newRunInfo(
00220                                 const boost::shared_ptr<LHERunInfo> &runInfo)
00221 {
00222         lhaUP->loadRunInfo(runInfo);
00223         pythia->init(lhaUP.get());
00224 }
00225 
00226 DEFINE_LHE_HADRONISATION_PLUGIN(Pythia8Hadronisation);
00227 
00228 } // namespace lhef