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