CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Py8GunBase.cc
Go to the documentation of this file.
4 
5 // EvtGen plugin
6 //
7 #include "Pythia8Plugins/EvtGen.h"
8 
9 using namespace Pythia8;
10 
12 
13 namespace gen {
14 
15 Py8GunBase::Py8GunBase( edm::ParameterSet const& ps )
16  : Py8InterfaceBase(ps)
17 {
18  // PGun specs
19  //
20  edm::ParameterSet pgun_params =
21  ps.getParameter<edm::ParameterSet>("PGunParameters");
22 
23  // although there's the method ParameterSet::empty(),
24  // it looks like it's NOT even necessary to check if it is,
25  // before trying to extract parameters - if it is empty,
26  // the default values seem to be taken
27  //
28  fPartIDs = pgun_params.getParameter< std::vector<int> >("ParticleID");
29  fMinPhi = pgun_params.getParameter<double>("MinPhi"); // ,-3.14159265358979323846);
30  fMaxPhi = pgun_params.getParameter<double>("MaxPhi"); // , 3.14159265358979323846);
31 
32 }
33 
34 // specific to Py8GunHad !!!
35 //
37 {
38 
39  // NO MATTER what was this setting below, override it before init
40  // - this is essencial for the PGun mode
41 
42  // Key requirement: switch off ProcessLevel, and thereby also PartonLevel.
43  fMasterGen->readString("ProcessLevel:all = off");
44  fMasterGen->readString("ProcessLevel::resonanceDecays=on");
45  fMasterGen->init();
46 
47  // init decayer
48  fDecayer->readString("ProcessLevel:all = off"); // Same trick!
49  fDecayer->readString("ProcessLevel::resonanceDecays=on");
50  fDecayer->init();
51 
52  if (useEvtGen) {
53  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
54 
55  evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile.c_str(), evtgenPdlFile.c_str()));
56 
57  for (unsigned int i=0; i<evtgenUserFiles.size(); i++) {
58  edm::FileInPath evtgenUserFile(evtgenUserFiles.at(i));
59  evtgenDecays->readDecayFile(evtgenUserFile.fullPath().c_str());
60  }
61  }
62 
63  return true;
64 
65 }
66 
68 {
69 
70  Event* pythiaEvent = &(fMasterGen->event);
71 
72  int NPartsBeforeDecays = pythiaEvent->size()-1; // do NOT count the very 1st "system" particle
73  // in Pythia8::Event record; it does NOT even
74  // get translated by the HepMCInterface to the
75  // HepMC::GenEvent record !!!
76  int NPartsAfterDecays = event().get()->particles_size();
77 
78  if(NPartsAfterDecays == NPartsBeforeDecays) return true;
79 
80  bool result = true;
81 
82  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
83  {
84 
85  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
86 
87  if ( part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id()) )
88  {
89  fDecayer->event.reset();
90  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
91  part->momentum().x(),
92  part->momentum().y(),
93  part->momentum().z(),
94  part->momentum().t(),
95  part->generated_mass() );
96  HepMC::GenVertex* ProdVtx = part->production_vertex();
97  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
98  ProdVtx->position().z(), ProdVtx->position().t() );
99  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
100  fDecayer->event.append( py8part );
101  int nentries = fDecayer->event.size();
102  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
103  fDecayer->next();
104  int nentries1 = fDecayer->event.size();
105  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
106 
107  part->set_status(2);
108 
109  result = toHepMC.fill_next_event( *(fDecayer.get()), event().get(), -1, true, part);
110 
111  }
112  }
113 
114  return result;
115 
116 }
117 
119 {
120 
121  //******** Verbosity ********
122 
123  if (maxEventsToPrint > 0 &&
126  {
129  {
130  fMasterGen->info.list();
131  fMasterGen->event.list();
132  }
133 
135  {
136  std::cout << "Event process = "
137  << fMasterGen->info.code() << "\n"
138  << "----------------------" << std::endl;
139  event()->print();
140  }
142  std::cout << "Event process = "
143  << fMasterGen->info.code() << "\n"
144  << "----------------------" << std::endl;
145  ascii_io->write_event(event().get());
146  }
147  }
148 
149  return;
150 }
151 
153 {
154 
155  fMasterGen->stat();
156 
157  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
158  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
159  runInfo().setInternalXSec(xsec);
160  return;
161 
162 }
163 
165 {
166  if (evtgenDecays.get()) evtgenDecays->decay();
167 }
168 
169 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
std::auto_ptr< EvtGenDecays > evtgenDecays
HepMC::IO_AsciiParticles * ascii_io
double fMinPhi
Definition: Py8GunBase.h:60
bool initializeForInternalPartons()
Definition: Py8GunBase.cc:36
void setInternalXSec(const XSec &xsec)
virtual bool residualDecay()
Definition: Py8GunBase.cc:67
std::auto_ptr< HepMC::GenEvent > & event()
std::vector< std::string > evtgenUserFiles
GenRunInfoProduct & runInfo()
void finalizeEvent()
Definition: Py8GunBase.cc:118
void evtGenDecay()
Definition: Py8GunBase.cc:164
tuple result
Definition: query.py:137
static const std::vector< std::string > p8SharedResources
Definition: Py8GunBase.h:64
unsigned int pythiaPylistVerbosity
std::vector< int > fPartIDs
Definition: Py8GunBase.h:59
part
Definition: HCALResponse.h:20
unsigned int maxEventsToPrint
double fMaxPhi
Definition: Py8GunBase.h:61
std::auto_ptr< Pythia8::Pythia > fDecayer
HepMC::Pythia8ToHepMC toHepMC
tuple cout
Definition: gather_cfg.py:121
static const std::string kPythia8