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 override
 
bool generatePartonsAndHadronize () override
 
 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 setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
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 > &)
 
void p8SetRandomEngine (CLHEP::HepRandomEngine *v)
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
P8RndmEnginerandomEngine ()
 
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
HepMC::IO_AsciiParticles * ascii_io
 
std::auto_ptr< Pythia8::Pythia > fDecayer
 
std::auto_ptr< Pythia8::Pythia > fMasterGen
 
ParameterCollector fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
bool pythiaHepMCVerbosityParticles
 
unsigned int pythiaPylistVerbosity
 
HepMC::Pythia8ToHepMC toHepMC
 

Detailed Description

Definition at line 9 of file Py8JetGun.cc.

Constructor & Destructor Documentation

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

Definition at line 33 of file Py8JetGun.cc.

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

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

Definition at line 14 of file Py8JetGun.cc.

14 {}

Member Function Documentation

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

Implements gen::Py8InterfaceBase.

Definition at line 124 of file Py8JetGun.cc.

125 {
126  return "Py8JetGun";
127 }
bool gen::Py8JetGun::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 49 of file Py8JetGun.cc.

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

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

Member Data Documentation

double gen::Py8JetGun::fMaxE
private

Definition at line 27 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMaxEta
private

Definition at line 23 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMaxP
private

Definition at line 25 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinE
private

Definition at line 26 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinEta
private

Definition at line 22 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().

double gen::Py8JetGun::fMinP
private

Definition at line 24 of file Py8JetGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8JetGun().