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
00018
00019 #include "GeneratorInterface/Pythia8Interface/interface/JetMatchingHook.h"
00020
00021
00022
00023 #include "GeneratorInterface/Pythia8Interface/interface/EmissionVetoHook.h"
00024 #include "GeneratorInterface/Pythia8Interface/interface/EmissionVetoHook1.h"
00025
00026 #include "FWCore/ServiceRegistry/interface/Service.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00029 #include "FWCore/ParameterSet/interface/FileInPath.h"
00030
00031 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00032 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00033
00034 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00035 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00036 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00037 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
00038 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00039
00040 #include "GeneratorInterface/Pythia8Interface/interface/LHAupLesHouches.h"
00041
00042 #include "HepPID/ParticleIDTranslations.hh"
00043
00044 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00045
00046 using namespace gen;
00047 using namespace Pythia8;
00048
00049 class Pythia8Hadronizer : public BaseHadronizer {
00050 public:
00051 Pythia8Hadronizer(const edm::ParameterSet ¶ms);
00052 ~Pythia8Hadronizer();
00053
00054 bool readSettings( int );
00055 bool initializeForInternalPartons();
00056 bool initializeForExternalPartons();
00057
00058 bool declareStableParticles(const std::vector<int> &pdgIds);
00059 bool declareSpecialSettings( const std::vector<std::string> );
00060
00061 void statistics();
00062
00063 bool generatePartonsAndHadronize();
00064 bool hadronize();
00065 bool decay();
00066 bool residualDecay();
00067 void finalizeEvent();
00068
00069 const char *classname() const { return "Pythia8Hadronizer"; }
00070
00071 private:
00072 ParameterCollector parameters;
00073
00075 double comEnergy;
00077 unsigned int pythiaPylistVerbosity;
00079 bool pythiaHepMCVerbosity;
00081 unsigned int maxEventsToPrint;
00082
00083 string LHEInputFileName;
00084
00085 std::auto_ptr<LHAupLesHouches> lhaUP;
00086
00087 std::auto_ptr<Pythia> pythia;
00088 std::auto_ptr<Pythia> decayer;
00089 Event* pythiaEvent;
00090 HepMC::I_Pythia8 toHepMC;
00091
00092 enum { PP, PPbar, ElectronPositron };
00093 int fInitialState ;
00094
00095 double fBeam1PZ;
00096 double fBeam2PZ;
00097
00098
00099
00100 UserHooks* fReweightUserHook;
00101
00102
00103
00104 JetMatchingHook* fJetMatchingHook;
00105
00106
00107
00108 EmissionVetoHook* fEmissionVetoHook;
00109 EmissionVetoHook1* fEmissionVetoHook1;
00110
00111 int EV1_nFinal;
00112 bool EV1_vetoOn;
00113 int EV1_maxVetoCount;
00114 int EV1_pThardMode;
00115 int EV1_pTempMode;
00116 int EV1_emittedMode;
00117 int EV1_pTdefMode;
00118 bool EV1_MPIvetoOn;
00119
00120 };
00121
00122
00123 Pythia8Hadronizer::Pythia8Hadronizer(const edm::ParameterSet ¶ms) :
00124 BaseHadronizer(params),
00125 parameters(params.getParameter<edm::ParameterSet>("PythiaParameters")),
00126 comEnergy(params.getParameter<double>("comEnergy")),
00127 pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
00128 pythiaHepMCVerbosity(params.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
00129 maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00130 LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
00131 fInitialState(PP),
00132 fReweightUserHook(0),
00133 fJetMatchingHook(0),
00134 fEmissionVetoHook(0),fEmissionVetoHook1(0)
00135 {
00136 randomEngine = &getEngineReference();
00137
00138
00139
00140
00141
00142
00143 RandomP8* RP8 = new RandomP8();
00144
00145
00146
00147 if ( params.exists( "PPbarInitialState" ) )
00148 {
00149 if ( fInitialState == PP )
00150 {
00151 fInitialState = PPbar;
00152 edm::LogInfo("GeneratorInterface|Pythia6Interface")
00153 << "Pythia6 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
00154 << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00155 std::cout << "Pythia6 will be initialized for PROTON-ANTIPROTON INITIAL STATE." << std::endl;
00156 std::cout << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00157 }
00158 else
00159 {
00160
00161 }
00162 }
00163 else if ( params.exists( "ElectronPositronInitialState" ) )
00164 {
00165 if ( fInitialState == PP )
00166 {
00167 fInitialState = ElectronPositron;
00168 edm::LogInfo("GeneratorInterface|Pythia6Interface")
00169 << "Pythia6 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
00170 << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00171 std::cout << "Pythia6 will be initialized for ELECTRON-POSITRON INITIAL STATE." << std::endl;
00172 std::cout << "This is a user-request change from the DEFAULT PROTON-PROTON initial state." << std::endl;
00173 }
00174 else
00175 {
00176
00177 }
00178 }
00179 else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
00180 {
00181
00182 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00183 <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
00184 }
00185
00186 pythia.reset(new Pythia);
00187 decayer.reset(new Pythia);
00188
00189 pythia->setRndmEnginePtr(RP8);
00190 decayer->setRndmEnginePtr(RP8);
00191
00192 if( params.exists( "SLHAFileForPythia8" ) ) {
00193 std::string slhafilenameshort = params.getParameter<string>("SLHAFileForPythia8");
00194 edm::FileInPath f1( slhafilenameshort );
00195 std::string slhafilename = f1.fullPath();
00196 std::string pythiacommandslha = std::string("SLHA:file = ") + slhafilename;
00197 pythia->readString(pythiacommandslha);
00198 for ( ParameterCollector::const_iterator line = parameters.begin();
00199 line != parameters.end(); ++line ) {
00200 if (line->find("SLHA:file") != std::string::npos)
00201 throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
00202 << "using Pythia8 card SLHA:file and Pythia8Interface card SLHAFileForPythia8"
00203 << std::endl;
00204 }
00205 }
00206
00207
00208
00209 if( params.exists( "reweightGen" ) )
00210 fReweightUserHook = new PtHatReweightUserHook();
00211
00212
00213
00214 if ( params.exists("jetMatching") )
00215 {
00216 edm::ParameterSet jmParams =
00217 params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
00218 fJetMatchingHook = new JetMatchingHook( jmParams, &pythia->info );
00219 }
00220
00221
00222
00223 if ( params.exists("emissionVeto") )
00224 {
00225 fEmissionVetoHook = new EmissionVetoHook(0);
00226 }
00227
00228 if ( params.exists("emissionVeto1") )
00229 {
00230 EV1_nFinal = -1;
00231 if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
00232 EV1_vetoOn = true;
00233 if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
00234 EV1_maxVetoCount = 10;
00235 if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
00236 EV1_pThardMode = 1;
00237 if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
00238 EV1_pTempMode = 0;
00239 if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
00240 if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
00241 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00242 <<" Wrong value for EV1_pTempMode code\n";
00243 EV1_emittedMode = 0;
00244 if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
00245 EV1_pTdefMode = 1;
00246 if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
00247 EV1_MPIvetoOn = false;
00248 if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
00249 fEmissionVetoHook1 = new EmissionVetoHook1(EV1_nFinal, EV1_vetoOn,
00250 EV1_maxVetoCount, EV1_pThardMode, EV1_pTempMode,
00251 EV1_emittedMode, EV1_pTdefMode, EV1_MPIvetoOn, 0);
00252 }
00253
00254 int NHooks=0;
00255 if(fReweightUserHook) NHooks++;
00256 if(fJetMatchingHook) NHooks++;
00257 if(fEmissionVetoHook) NHooks++;
00258 if(fEmissionVetoHook1) NHooks++;
00259 if(NHooks > 1)
00260 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00261 <<" Too many User Hooks. \n Please choose one from: reweightGen, jetMatching, emissionVeto \n";
00262
00263 if(fReweightUserHook) pythia->setUserHooksPtr(fReweightUserHook);
00264 if(fJetMatchingHook) pythia->setUserHooksPtr(fJetMatchingHook);
00265 if(fEmissionVetoHook || fEmissionVetoHook1) {
00266 cout << "Turning on Emission Veto Hook";
00267 if(fEmissionVetoHook1) cout << " 1";
00268 cout << endl;
00269 int nversion = (int)(1000.*(pythia->settings.parm("Pythia:versionNumber") - 8.));
00270 if(nversion < 157) {
00271 cout << "obsolete pythia8 version for this Emission Veto code" << endl;
00272 cout << "Please update pythia8 version using the instructions here:" << endl;
00273 cout << "https://twiki.cern.ch/twiki/bin/view/CMS/Pythia8Interface" << endl;
00274 cout << "or try to use tag V00-01-28 of this interface" << endl;
00275 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00276 <<" Obsolete pythia8 version for this Emission Veto code\n";
00277 }
00278 if(fEmissionVetoHook) pythia->setUserHooksPtr(fEmissionVetoHook);
00279 if(fEmissionVetoHook1) pythia->setUserHooksPtr(fEmissionVetoHook1);
00280 }
00281 }
00282
00283
00284 Pythia8Hadronizer::~Pythia8Hadronizer()
00285 {
00286
00287
00288 if(fEmissionVetoHook) {delete fEmissionVetoHook; fEmissionVetoHook=0;}
00289 if(fEmissionVetoHook1) {delete fEmissionVetoHook1; fEmissionVetoHook1=0;}
00290 }
00291
00292
00293 bool Pythia8Hadronizer::readSettings( int )
00294 {
00295 for ( ParameterCollector::const_iterator line = parameters.begin();
00296 line != parameters.end(); ++line ) {
00297 if (line->find("Random:") != std::string::npos)
00298 throw cms::Exception("PythiaError") << "Attempted to set random number "
00299 "using Pythia commands. Please use " "the RandomNumberGeneratorService."
00300 << std::endl;
00301
00302 if (!pythia->readString(*line)) throw cms::Exception("PythiaError")
00303 << "Pythia 8 did not accept \""
00304 << *line << "\"." << std::endl;
00305 }
00306
00307 if ( pythiaPylistVerbosity > 10 ) {
00308 if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
00309 pythia->settings.listAll();
00310 if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
00311 pythia->particleData.listAll();
00312 }
00313
00314 return true;
00315 }
00316
00317
00318 bool Pythia8Hadronizer::initializeForInternalPartons()
00319 {
00320
00321 pythiaEvent = &pythia->event;
00322
00323 if ( fInitialState == PP )
00324 {
00325 pythia->init(2212, 2212, comEnergy);
00326 }
00327 else if ( fInitialState == PPbar )
00328 {
00329 pythia->init(2212, -2212, comEnergy);
00330 }
00331 else if ( fInitialState == ElectronPositron )
00332 {
00333 pythia->init(11, -11, comEnergy);
00334 }
00335 else
00336 {
00337
00338 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00339 <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
00340 }
00341
00342 pythia->settings.listChanged();
00343
00344 return true;
00345 }
00346
00347
00348 bool Pythia8Hadronizer::initializeForExternalPartons()
00349 {
00350
00351 std::cout << "Initializing for external partons" << std::endl;
00352
00353 pythiaEvent = &pythia->event;
00354
00355 if(LHEInputFileName != string()) {
00356
00357 cout << endl;
00358 cout << "LHE Input File Name = " << LHEInputFileName << endl;
00359 cout << endl;
00360 pythia->init(LHEInputFileName);
00361
00362 } else {
00363
00364 lhaUP.reset(new LHAupLesHouches());
00365 lhaUP->loadRunInfo(lheRunInfo());
00366 pythia->init(lhaUP.get());
00367
00368 }
00369
00370
00371
00372 if ( fJetMatchingHook )
00373 {
00374
00375
00376 fJetMatchingHook->init ( lheRunInfo() );
00377 }
00378
00379 return true;
00380 }
00381
00382
00383 #if 0
00384
00385 static int getStatus(const HepMC::GenParticle *p)
00386 {
00387 int status = p->status();
00388 if (status > 0)
00389 return status;
00390 else if (status > -30 && status < 0)
00391 return 3;
00392 else
00393 return 2;
00394 }
00395 #endif
00396
00397
00398 bool Pythia8Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00399 {
00400 for ( size_t i=0; i<pdgIds.size(); i++ )
00401 {
00402
00403
00404
00405
00406
00407
00408 int PyID = pdgIds[i];
00409 std::ostringstream pyCard ;
00410 pyCard << PyID <<":mayDecay=false";
00411 pythia->readString( pyCard.str() );
00412
00413
00414
00415
00416 }
00417
00418
00419 decayer->readString("ProcessLevel:all = off");
00420 decayer->init();
00421
00422 return true;
00423 }
00424
00425
00426 bool Pythia8Hadronizer::declareSpecialSettings( const std::vector<std::string> settings )
00427 {
00428 for ( unsigned int iss=0; iss<settings.size(); iss++ )
00429 {
00430 if ( settings[iss].find("QED-brem-off") == std::string::npos ) continue;
00431 pythia->readString( "TimeShower:QEDshowerByL=off" );
00432 }
00433
00434 return true;
00435 }
00436
00437
00438 void Pythia8Hadronizer::statistics()
00439 {
00440 pythia->statistics();
00441
00442 double xsec = pythia->info.sigmaGen();
00443 xsec *= 1.0e9;
00444 runInfo().setInternalXSec(xsec);
00445 }
00446
00447
00448 bool Pythia8Hadronizer::generatePartonsAndHadronize()
00449 {
00450 if (!pythia->next()) return false;
00451
00452 event().reset(new HepMC::GenEvent);
00453 toHepMC.fill_next_event(*pythiaEvent, event().get());
00454
00455 return true;
00456 }
00457
00458
00459 bool Pythia8Hadronizer::hadronize()
00460 {
00461 if(LHEInputFileName == string()) lhaUP->loadEvent(lheEvent());
00462
00463 if ( fJetMatchingHook )
00464 {
00465 fJetMatchingHook->resetMatchingStatus();
00466 fJetMatchingHook->beforeHadronization( lheEvent() );
00467 }
00468
00469 bool py8next = pythia->next();
00470
00471 if (!py8next)
00472 {
00473 lheEvent()->count( lhef::LHERunInfo::kSelected );
00474 event().reset();
00475 return false;
00476 }
00477
00478
00479
00480 lheEvent()->count( lhef::LHERunInfo::kAccepted );
00481
00482 event().reset(new HepMC::GenEvent);
00483 toHepMC.fill_next_event(*pythiaEvent, event().get());
00484
00485 return true;
00486 }
00487
00488
00489 bool Pythia8Hadronizer::decay()
00490 {
00491 return true;
00492 }
00493
00494
00495 bool Pythia8Hadronizer::residualDecay()
00496 {
00497
00498 int NPartsBeforeDecays = pythiaEvent->size();
00499 int NPartsAfterDecays = event().get()->particles_size();
00500 int NewBarcode = NPartsAfterDecays;
00501
00502 for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
00503 {
00504
00505 HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
00506
00507 if ( part->status() == 1 )
00508 {
00509 decayer->event.reset();
00510 Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
00511 part->momentum().x(),
00512 part->momentum().y(),
00513 part->momentum().z(),
00514 part->momentum().t(),
00515 part->generated_mass() );
00516 HepMC::GenVertex* ProdVtx = part->production_vertex();
00517 py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
00518 ProdVtx->position().z(), ProdVtx->position().t() );
00519 py8part.tau( (decayer->particleData).tau0( part->pdg_id() ) );
00520 decayer->event.append( py8part );
00521 int nentries = decayer->event.size();
00522 if ( !decayer->event[nentries-1].mayDecay() ) continue;
00523 decayer->next();
00524 int nentries1 = decayer->event.size();
00525
00526 if ( nentries1 <= nentries ) continue;
00527
00528 part->set_status(2);
00529
00530 Particle& py8daughter = decayer->event[nentries];
00531 HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
00532 py8daughter.yProd(),
00533 py8daughter.zProd(),
00534 py8daughter.tProd()) );
00535
00536 DecVtx->add_particle_in( part );
00537
00538
00539 HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
00540
00541 HepMC::GenParticle* daughter =
00542 new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
00543
00544 NewBarcode++;
00545 daughter->suggest_barcode( NewBarcode );
00546 DecVtx->add_particle_out( daughter );
00547
00548 for ( ipart=nentries+1; ipart<nentries1; ipart++ )
00549 {
00550 py8daughter = decayer->event[ipart];
00551 HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
00552 HepMC::GenParticle* daughterN =
00553 new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
00554 NewBarcode++;
00555 daughterN->suggest_barcode( NewBarcode );
00556 DecVtx->add_particle_out( daughterN );
00557 }
00558
00559 event().get()->add_vertex( DecVtx );
00560
00561 }
00562 }
00563
00564 return true;
00565 }
00566
00567
00568 void Pythia8Hadronizer::finalizeEvent()
00569 {
00570 bool lhe = lheEvent() != 0;
00571
00572 event()->set_signal_process_id(pythia->info.code());
00573 event()->set_event_scale(pythia->info.pTHat());
00574
00575
00576
00577
00578
00579
00580 if (event()->alphaQED() <= 0)
00581 event()->set_alphaQED( pythia->info.alphaEM() );
00582 if (event()->alphaQCD() <= 0)
00583 event()->set_alphaQCD( pythia->info.alphaS() );
00584
00585 HepMC::GenCrossSection xsec;
00586 xsec.set_cross_section( pythia->info.sigmaGen() * 1e9,
00587 pythia->info.sigmaErr() * 1e9);
00588 event()->set_cross_section(xsec);
00589
00590
00591
00592
00593
00594
00595
00596
00597 int id1 = pythia->info.id1();
00598 int id2 = pythia->info.id2();
00599 if (id1 == 21) id1 = 0;
00600 if (id2 == 21) id2 = 0;
00601 double x1 = pythia->info.x1();
00602 double x2 = pythia->info.x2();
00603
00604 double Q = pythia->info.QFac();
00605 double pdf1 = pythia->info.pdf1() / pythia->info.x1();
00606 double pdf2 = pythia->info.pdf2() / pythia->info.x2();
00607 event()->set_pdf_info(HepMC::PdfInfo(id1,id2,x1,x2,Q,pdf1,pdf2));
00608
00609
00610
00611 if (lhe && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
00612
00613 event()->weights().push_back( pythia->info.weight() * 1.0e9 );
00614 else
00615 event()->weights().push_back( pythia->info.weight() );
00616
00617
00618
00619 eventInfo().reset( new GenEventInfoProduct( event().get() ) );
00620
00621
00622
00623
00624 if (!lhe) {
00625 eventInfo()->setBinningValues(std::vector<double>(1, pythia->info.pTHat()));
00626 }
00627
00628
00629
00630 if (maxEventsToPrint > 0 &&
00631 (pythiaPylistVerbosity || pythiaHepMCVerbosity)) {
00632 maxEventsToPrint--;
00633 if (pythiaPylistVerbosity) {
00634 pythia->info.list(std::cout);
00635 pythia->event.list(std::cout);
00636 }
00637
00638 if (pythiaHepMCVerbosity) {
00639 std::cout << "Event process = "
00640 << pythia->info.code() << "\n"
00641 << "----------------------" << std::endl;
00642 event()->print();
00643 }
00644 }
00645 }
00646
00647
00648 typedef edm::GeneratorFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8GeneratorFilter;
00649 DEFINE_FWK_MODULE(Pythia8GeneratorFilter);
00650
00651
00652 typedef edm::HadronizerFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8HadronizerFilter;
00653 DEFINE_FWK_MODULE(Pythia8HadronizerFilter);