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.
1 
3 // #include "GeneratorInterface/Pythia8Interface/interface/RandomP8.h"
4 // #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
6 
7 using namespace Pythia8;
8 
10 
11 namespace gen {
12 
13 Py8GunBase::Py8GunBase( edm::ParameterSet const& ps )
14  : Py8InterfaceBase(ps)
15 {
16 
18  ps.getUntrackedParameter<double>("filterEfficiency", -1.) );
20  GenRunInfoProduct::XSec(ps.getUntrackedParameter<double>("crossSection", -1.)) );
22  GenRunInfoProduct::XSec(ps.getUntrackedParameter<double>("crossSectionNLO", -1.)) );
23 
24 
25  // PGun specs
26  //
27  edm::ParameterSet pgun_params =
28  ps.getParameter<edm::ParameterSet>("PGunParameters");
29 
30  // although there's the method ParameterSet::empty(),
31  // it looks like it's NOT even necessary to check if it is,
32  // before trying to extract parameters - if it is empty,
33  // the default values seem to be taken
34  //
35  fPartIDs = pgun_params.getParameter< std::vector<int> >("ParticleID");
36  fMinPhi = pgun_params.getParameter<double>("MinPhi"); // ,-3.14159265358979323846);
37  fMaxPhi = pgun_params.getParameter<double>("MaxPhi"); // , 3.14159265358979323846);
38 
39 }
40 
41 // specific to Py8GunHad !!!
42 //
44 {
45 
46  // NO MATTER what was this setting below, override it before init
47  // - this is essencial for the PGun mode
48 
49  // Key requirement: switch off ProcessLevel, and thereby also PartonLevel.
50  fMasterGen->readString("ProcessLevel:all = off");
51  fMasterGen->readString("Standalone:allowResDec=on");
52  // pythia->readString("ProcessLevel::resonanceDecays=on");
53  fMasterGen->init();
54 
55  // init decayer
56  fDecayer->readString("ProcessLevel:all = off"); // Same trick!
57  fDecayer->readString("Standalone:allowResDec=on");
58  // pythia->readString("ProcessLevel::resonanceDecays=on");
59  fDecayer->init();
60 
61  return true;
62 
63 }
64 
66 {
67 
68  Event* pythiaEvent = &(fMasterGen->event);
69 
70  int NPartsBeforeDecays = pythiaEvent->size()-1; // do NOT count the very 1st "system" particle
71  // in Pythia8::Event record; it does NOT even
72  // get translated by the HepMCInterface to the
73  // HepMC::GenEvent record !!!
74  int NPartsAfterDecays = event().get()->particles_size();
75  int NewBarcode = NPartsAfterDecays;
76 
77  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
78  {
79 
80  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
81 
82  if ( part->status() == 1 )
83  {
84  fDecayer->event.reset();
85  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
86  part->momentum().x(),
87  part->momentum().y(),
88  part->momentum().z(),
89  part->momentum().t(),
90  part->generated_mass() );
91  HepMC::GenVertex* ProdVtx = part->production_vertex();
92  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
93  ProdVtx->position().z(), ProdVtx->position().t() );
94  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
95  fDecayer->event.append( py8part );
96  int nentries = fDecayer->event.size();
97  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
98  fDecayer->next();
99  int nentries1 = fDecayer->event.size();
100  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
101 
102  part->set_status(2);
103 
104  Particle& py8daughter = fDecayer->event[nentries]; // the 1st daughter
105  HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
106  py8daughter.yProd(),
107  py8daughter.zProd(),
108  py8daughter.tProd()) );
109 
110  DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
111  // I presume (vtx) barcode will be given automatically
112 
113  HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
114 
115  HepMC::GenParticle* daughter =
116  new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
117 
118  NewBarcode++;
119  daughter->suggest_barcode( NewBarcode );
120  DecVtx->add_particle_out( daughter );
121 
122  for ( ipart=nentries+1; ipart<nentries1; ipart++ )
123  {
124  py8daughter = fDecayer->event[ipart];
125  HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
126  HepMC::GenParticle* daughterN =
127  new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
128  NewBarcode++;
129  daughterN->suggest_barcode( NewBarcode );
130  DecVtx->add_particle_out( daughterN );
131  }
132 
133  event().get()->add_vertex( DecVtx );
134 
135  }
136  }
137 
138  return true;
139 
140 }
141 
143 {
144 
145 
146 // FIXME !!!
147 
148  //******** Verbosity ********
149 
150  if (maxEventsToPrint > 0 &&
152  {
155  {
156  fMasterGen->info.list(std::cout);
157  fMasterGen->event.list(std::cout);
158  }
159 
161  {
162  std::cout << "Event process = "
163  << fMasterGen->info.code() << "\n"
164  << "----------------------" << std::endl;
165  event()->print();
166  }
167  }
168 
169  return;
170 }
171 
173 {
174 
175  fMasterGen->statistics();
176 
177  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
178  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
179  runInfo().setInternalXSec(xsec);
180  return;
181 
182 }
183 
184 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::auto_ptr< Pythia8::Pythia > fMasterGen
void setFilterEfficiency(double effic)
double fMinPhi
Definition: Py8GunBase.h:75
bool initializeForInternalPartons()
Definition: Py8GunBase.cc:43
void setInternalXSec(const XSec &xsec)
virtual bool residualDecay()
Definition: Py8GunBase.cc:65
void finalizeEvent()
Definition: Py8GunBase.cc:142
static const std::vector< std::string > p8SharedResources
Definition: Py8GunBase.h:85
unsigned int pythiaPylistVerbosity
std::vector< int > fPartIDs
Definition: Py8GunBase.h:74
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:69
part
Definition: HCALResponse.h:20
unsigned int maxEventsToPrint
double fMaxPhi
Definition: Py8GunBase.h:76
std::auto_ptr< Pythia8::Pythia > fDecayer
tuple cout
Definition: gather_cfg.py:121
static const std::string kPythia8
void setExternalXSecNLO(const XSec &xsec)
GenRunInfoProduct & runInfo()
Definition: Py8GunBase.h:68
void setExternalXSecLO(const XSec &xsec)