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