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