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