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("Standalone:allowResDec=on");
46  // pythia->readString("ProcessLevel::resonanceDecays=on");
47  fMasterGen->init();
48 
49  // init decayer
50  fDecayer->readString("ProcessLevel:all = off"); // Same trick!
51  fDecayer->readString("Standalone:allowResDec=on");
52  // pythia->readString("ProcessLevel::resonanceDecays=on");
53  fDecayer->init();
54 
55  return true;
56 
57 }
58 
60 {
61 
62  Event* pythiaEvent = &(fMasterGen->event);
63 
64  int NPartsBeforeDecays = pythiaEvent->size()-1; // do NOT count the very 1st "system" particle
65  // in Pythia8::Event record; it does NOT even
66  // get translated by the HepMCInterface to the
67  // HepMC::GenEvent record !!!
68  int NPartsAfterDecays = event().get()->particles_size();
69  int NewBarcode = NPartsAfterDecays;
70 
71  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
72  {
73 
74  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
75 
76  if ( part->status() == 1 )
77  {
78  fDecayer->event.reset();
79  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
80  part->momentum().x(),
81  part->momentum().y(),
82  part->momentum().z(),
83  part->momentum().t(),
84  part->generated_mass() );
85  HepMC::GenVertex* ProdVtx = part->production_vertex();
86  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
87  ProdVtx->position().z(), ProdVtx->position().t() );
88  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
89  fDecayer->event.append( py8part );
90  int nentries = fDecayer->event.size();
91  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
92  fDecayer->next();
93  int nentries1 = fDecayer->event.size();
94  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
95 
96  part->set_status(2);
97 
98  Particle& py8daughter = fDecayer->event[nentries]; // the 1st daughter
99  HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
100  py8daughter.yProd(),
101  py8daughter.zProd(),
102  py8daughter.tProd()) );
103 
104  DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
105  // I presume (vtx) barcode will be given automatically
106 
107  HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
108 
109  HepMC::GenParticle* daughter =
110  new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
111 
112  NewBarcode++;
113  daughter->suggest_barcode( NewBarcode );
114  DecVtx->add_particle_out( daughter );
115 
116  for ( int ipart1=nentries+1; ipart1<nentries1; ipart1++ )
117  {
118  py8daughter = fDecayer->event[ipart1];
119  HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
120  HepMC::GenParticle* daughterN =
121  new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
122  NewBarcode++;
123  daughterN->suggest_barcode( NewBarcode );
124  DecVtx->add_particle_out( daughterN );
125  }
126 
127  event().get()->add_vertex( DecVtx );
128 
129  }
130  }
131 
132  return true;
133 
134 }
135 
137 {
138 
139 
140 // FIXME !!!
141 
142  //******** Verbosity ********
143 
144  if (maxEventsToPrint > 0 &&
146  {
149  {
150  fMasterGen->info.list(std::cout);
151  fMasterGen->event.list(std::cout);
152  }
153 
155  {
156  std::cout << "Event process = "
157  << fMasterGen->info.code() << "\n"
158  << "----------------------" << std::endl;
159  event()->print();
160  }
161  }
162 
163  return;
164 }
165 
167 {
168 
169  fMasterGen->statistics();
170 
171  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
172  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
173  runInfo().setInternalXSec(xsec);
174  return;
175 
176 }
177 
178 }
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:37
void setInternalXSec(const XSec &xsec)
virtual bool residualDecay()
Definition: Py8GunBase.cc:59
void finalizeEvent()
Definition: Py8GunBase.cc:136
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)