CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
gen::Py8JetGun Class Reference
Inheritance diagram for gen::Py8JetGun:
gen::Py8GunBase gen::Py8InterfaceBase

Public Member Functions

const char * classname () const
 
bool generatePartonsAndHadronize ()
 
 Py8JetGun (edm::ParameterSet const &)
 
 ~Py8JetGun ()
 
- Public Member Functions inherited from gen::Py8GunBase
void finalizeEvent ()
 
edm::EventgetEDMEvent () const
 
HepMC::GenEvent * getGenEvent ()
 
GenEventInfoProductgetGenEventInfo ()
 
GenRunInfoProductgetGenRunInfo ()
 
bool initializeForInternalPartons ()
 
 Py8GunBase (edm::ParameterSet const &ps)
 
void resetEvent (HepMC::GenEvent *event)
 
void resetEventInfo (GenEventInfoProduct *eventInfo)
 
virtual bool residualDecay ()
 
virtual bool select (HepMC::GenEvent *) const
 
void setEDMEvent (edm::Event &event)
 
void statistics ()
 
 ~Py8GunBase ()
 
- Public Member Functions inherited from gen::Py8InterfaceBase
bool decay ()
 
bool declareSpecialSettings (const std::vector< std::string > &)
 
bool declareStableParticles (const std::vector< int > &)
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
bool readSettings (int)
 
 ~Py8InterfaceBase ()
 

Private Attributes

double fMaxE
 
double fMaxEta
 
double fMaxP
 
double fMinE
 
double fMinEta
 
double fMinP
 

Additional Inherited Members

- Protected Member Functions inherited from gen::Py8GunBase
std::auto_ptr< HepMC::GenEvent > & event ()
 
std::auto_ptr
< GenEventInfoProduct > & 
eventInfo ()
 
GenRunInfoProductrunInfo ()
 
- Protected Attributes inherited from gen::Py8GunBase
double fMaxPhi
 
double fMinPhi
 
std::vector< int > fPartIDs
 
- Protected Attributes inherited from gen::Py8InterfaceBase
std::auto_ptr< Pythia8::Pythia > fDecayer
 
std::auto_ptr< Pythia8::Pythia > fMasterGen
 
ParameterCollector fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
unsigned int pythiaPylistVerbosity
 
HepMC::Pythia8ToHepMC toHepMC
 

Detailed Description

Definition at line 10 of file Py8JetGun.cc.

Constructor & Destructor Documentation

gen::Py8JetGun::Py8JetGun ( edm::ParameterSet const &  ps)

Definition at line 34 of file Py8JetGun.cc.

References fMaxE, fMaxEta, fMaxP, fMinE, fMinEta, fMinP, and edm::ParameterSet::getParameter().

35  : Py8GunBase(ps)
36 {
37 
38  // ParameterSet defpset ;
39  edm::ParameterSet pgun_params =
40  ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
41  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
42  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
43  fMinP = pgun_params.getParameter<double>("MinP"); // , 0.);
44  fMaxP = pgun_params.getParameter<double>("MaxP"); // , 0.);
45  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
46  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 0.);
47 
48 }
T getParameter(std::string const &) const
double fMinEta
Definition: Py8JetGun.cc:23
double fMinE
Definition: Py8JetGun.cc:27
double fMaxP
Definition: Py8JetGun.cc:26
double fMaxE
Definition: Py8JetGun.cc:28
Py8GunBase(edm::ParameterSet const &ps)
Definition: Py8GunBase.cc:7
double fMinP
Definition: Py8JetGun.cc:25
double fMaxEta
Definition: Py8JetGun.cc:24
gen::Py8JetGun::~Py8JetGun ( )
inline

Definition at line 15 of file Py8JetGun.cc.

15 {}

Member Function Documentation

const char * gen::Py8JetGun::classname ( ) const
virtual

Implements gen::Py8InterfaceBase.

Definition at line 129 of file Py8JetGun.cc.

130 {
131  return "Py8JetGun";
132 }
bool gen::Py8JetGun::generatePartonsAndHadronize ( )
virtual

Implements gen::Py8InterfaceBase.

Definition at line 50 of file Py8JetGun.cc.

References python.multivaluedict::append(), funct::cos(), eta(), gen::Py8GunBase::event(), create_public_lumi_plots::exp, gen::Py8InterfaceBase::fMasterGen, fMaxE, fMaxEta, fMaxP, gen::Py8GunBase::fMaxPhi, fMinE, fMinEta, fMinP, gen::Py8GunBase::fMinPhi, gen::Py8GunBase::fPartIDs, i, M_PI, scaleCards::mass, phi, createTree::pp, randomEngine, funct::sin(), mathSSE::sqrt(), and gen::Py8InterfaceBase::toHepMC.

51 {
52 
53  fMasterGen->event.reset();
54 
55  double totPx = 0.;
56  double totPy = 0.;
57  double totPz = 0.;
58  double totE = 0.;
59  double totM = 0.;
60  double phi, eta, the, ee, pp;
61 
62  for ( size_t i=0; i<fPartIDs.size(); i++ )
63  {
64 
65  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
66 
67  // FIXME !!!
68  // Ouch, it's using bare randomEngine pointer - that's NOT safe.
69  // Need to hold a pointer somewhere properly !!!
70  //
71  phi = 2. * M_PI * randomEngine->flat() ;
72  the = acos( -1. + 2.*randomEngine->flat() );
73 
74  // from input
75  //
76  ee = (fMaxE-fMinE)*randomEngine->flat() + fMinE;
77 
78  double mass = (fMasterGen->particleData).m0( particleID );
79 
80  pp = sqrt( ee*ee - mass*mass );
81 
82  double px = pp * sin(the) * cos(phi);
83  double py = pp * sin(the) * sin(phi);
84  double pz = pp * cos(the);
85 
86  if ( !((fMasterGen->particleData).isParticle( particleID )) )
87  {
88  particleID = std::fabs(particleID) ;
89  }
90  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
91 
92  // values for computing total mass
93  //
94  totPx += px;
95  totPy += py;
96  totPz += pz;
97  totE += ee;
98 
99  }
100 
101  totM = sqrt( totE*totE - (totPx*totPx+totPy*totPy+totPz*totPz) );
102 
103  //now the boost (from input params)
104  //
105  pp = (fMaxP-fMinP)*randomEngine->flat() + fMinP;
106  ee = sqrt( totM*totM + pp*pp );
107 
108  //the boost direction (from input params)
109  //
110  phi = (fMaxPhi-fMinPhi)*randomEngine->flat() + fMinPhi;
111  eta = (fMaxEta-fMinEta)*randomEngine->flat() + fMinEta;
112  the = 2.*atan(exp(-eta));
113 
114  double betaX = pp/ee * std::sin(the) * std::cos(phi);
115  double betaY = pp/ee * std::sin(the) * std::sin(phi);
116  double betaZ = pp/ee * std::cos(the);
117 
118  // boost all particles
119  //
120  (fMasterGen->event).bst( betaX, betaY, betaZ );
121 
122  if ( !fMasterGen->next() ) return false;
123 
124  event().reset(new HepMC::GenEvent);
125  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
126 
127 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
double fMinEta
Definition: Py8JetGun.cc:23
double fMinE
Definition: Py8JetGun.cc:27
tuple pp
Definition: createTree.py:15
CLHEP::HepRandomEngine * randomEngine
Definition: PYR.cc:4
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:66
double fMaxP
Definition: Py8JetGun.cc:26
T eta() const
double fMaxE
Definition: Py8JetGun.cc:28
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< int > fPartIDs
Definition: Py8GunBase.h:65
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:60
#define M_PI
Definition: BFit3D.cc:3
double fMaxPhi
Definition: Py8GunBase.h:67
tuple mass
Definition: scaleCards.py:27
HepMC::Pythia8ToHepMC toHepMC
double fMinP
Definition: Py8JetGun.cc:25
double fMaxEta
Definition: Py8JetGun.cc:24
Definition: DDAxes.h:10

Member Data Documentation

double gen::Py8JetGun::fMaxE
private

Definition at line 28 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMaxEta
private

Definition at line 24 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMaxP
private

Definition at line 26 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinE
private

Definition at line 27 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinEta
private

Definition at line 23 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinP
private

Definition at line 25 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().