CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
gen::Py8PtGun Class Reference
Inheritance diagram for gen::Py8PtGun:
gen::Py8GunBase gen::Py8InterfaceBase gen::BaseHadronizer

Public Member Functions

const char * classname () const override
 
bool generatePartonsAndHadronize () override
 
 Py8PtGun (edm::ParameterSet const &)
 
 ~Py8PtGun () override
 
- Public Member Functions inherited from gen::Py8GunBase
void evtGenDecay ()
 
void finalizeEvent () override
 
bool initializeForInternalPartons () override
 
 Py8GunBase (edm::ParameterSet const &ps)
 
virtual bool residualDecay ()
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
void statistics () override
 
 ~Py8GunBase () override
 
- 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 () override
 
- Public Member Functions inherited from gen::BaseHadronizer
 BaseHadronizer (edm::ParameterSet const &ps)
 
void cleanLHE ()
 
void generateLHE (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine, unsigned int ncpu)
 
edm::EventgetEDMEvent () const
 
std::unique_ptr< HepMC::GenEventgetGenEvent ()
 
std::unique_ptr< GenEventInfoProductgetGenEventInfo ()
 
virtual std::unique_ptr< GenLumiInfoHeadergetGenLumiInfoHeader () const
 
GenRunInfoProductgetGenRunInfo ()
 
std::unique_ptr< lhef::LHEEventgetLHEEvent ()
 
const std::shared_ptr< lhef::LHERunInfo > & getLHERunInfo () const
 
const std::string & gridpackPath () const
 
int randomIndex () const
 
const std::string & randomInitConfigDescription () const
 
void randomizeIndex (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine)
 
void resetEvent (std::unique_ptr< HepMC::GenEvent > event)
 
void resetEventInfo (std::unique_ptr< GenEventInfoProduct > eventInfo)
 
virtual bool select (HepMC::GenEvent *) const
 
void setEDMEvent (edm::Event &event)
 
void setLHEEvent (std::unique_ptr< lhef::LHEEvent > event)
 
void setLHERunInfo (std::unique_ptr< lhef::LHERunInfo > runInfo)
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
virtual ~BaseHadronizer ()(false)
 

Private Attributes

bool fAddAntiParticle
 
double fMaxEta
 
double fMaxPt
 
double fMinEta
 
double fMinPt
 

Additional Inherited Members

- Protected Member Functions inherited from gen::BaseHadronizer
std::unique_ptr< HepMC::GenEvent > & event ()
 
std::unique_ptr< GenEventInfoProduct > & eventInfo ()
 
lhef::LHEEventlheEvent ()
 
lhef::LHERunInfolheRunInfo ()
 
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::shared_ptr< EvtGenDecays > evtgenDecays
 
std::string evtgenDecFile
 
std::string evtgenPdlFile
 
std::vector< std::string > evtgenUserFiles
 
std::unique_ptr< Pythia8::Pythia > fDecayer
 
std::unique_ptr< Pythia8::Pythia > fMasterGen
 
edm::ParameterSet fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
bool pythiaHepMCVerbosityParticles
 
unsigned int pythiaPylistVerbosity
 
std::string slhafile_
 
HepMC::Pythia8ToHepMC toHepMC
 
bool useEvtGen
 
- Protected Attributes inherited from gen::BaseHadronizer
std::string lheFile_
 
int randomIndex_
 

Detailed Description

Definition at line 8 of file Py8PtGun.cc.

Constructor & Destructor Documentation

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

Definition at line 31 of file Py8PtGun.cc.

References fAddAntiParticle, fMaxEta, fMaxPt, fMinEta, fMinPt, and edm::ParameterSet::getParameter().

32  : Py8GunBase(ps) {
33 
34  // ParameterSet defpset ;
35  edm::ParameterSet pgun_params =
36  ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
37  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
38  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
39  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
40  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 0.);
41  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
42 
43 }
T getParameter(std::string const &) const
bool fAddAntiParticle
Definition: Py8PtGun.cc:25
double fMinPt
Definition: Py8PtGun.cc:23
double fMinEta
Definition: Py8PtGun.cc:21
Py8GunBase(edm::ParameterSet const &ps)
Definition: Py8GunBase.cc:15
double fMaxPt
Definition: Py8PtGun.cc:24
double fMaxEta
Definition: Py8PtGun.cc:22
gen::Py8PtGun::~Py8PtGun ( )
inlineoverride

Definition at line 13 of file Py8PtGun.cc.

References classname(), and generatePartonsAndHadronize().

13 {}

Member Function Documentation

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

Implements gen::Py8InterfaceBase.

Definition at line 116 of file Py8PtGun.cc.

Referenced by ~Py8PtGun().

117 {
118  return "Py8PtGun";
119 }
bool gen::Py8PtGun::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 45 of file Py8PtGun.cc.

References funct::abs(), mps_setup::append, funct::cos(), PVValHelper::eta, gen::BaseHadronizer::event(), gen::Py8GunBase::evtGenDecay(), JetChargeProducer_cfi::exp, fAddAntiParticle, gen::P8RndmEngine::flat(), gen::Py8InterfaceBase::fMasterGen, fMaxEta, gen::Py8GunBase::fMaxPhi, fMaxPt, fMinEta, gen::Py8GunBase::fMinPhi, fMinPt, gen::Py8GunBase::fPartIDs, mps_fire::i, cmsBatch::log, ResonanceBuilder::mass, EgammaObjectsElectrons_cfi::particleID, createTree::pp, EnergyCorrector::pt, gen::Py8InterfaceBase::randomEngine(), funct::sin(), findQualityFiles::size, mathSSE::sqrt(), metsig::tau, and gen::Py8InterfaceBase::toHepMC.

Referenced by ~Py8PtGun().

46 {
47 
48  fMasterGen->event.reset();
49 
50  for ( size_t i=0; i<fPartIDs.size(); i++ ){
51 
52  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
53 
54  double phi = (fMaxPhi-fMinPhi) * randomEngine().flat() + fMinPhi;
55  double eta = (fMaxEta-fMinEta) * randomEngine().flat() + fMinEta;
56  double the = 2.*atan(exp(-eta));
57 
58  double pt = (fMaxPt-fMinPt) * randomEngine().flat() + fMinPt;
59 
60  double mass = (fMasterGen->particleData).m0( particleID );
61 
62  double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
63  double ee = sqrt( pp*pp + mass*mass );
64 
65  double px = pt * cos(phi);
66  double py = pt * sin(phi);
67  double pz = pp * cos(the);
68 
69  if ( !((fMasterGen->particleData).isParticle( particleID )) ){
70  particleID = std::abs(particleID) ;
71  }
72  if( 1<= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
73  (fMasterGen->event).append( particleID, 23, 101, 0, px, py, pz, ee, mass );
74  else if (std::abs(particleID) == 21) // gluons
75  (fMasterGen->event).append( 21, 23, 101, 102, px, py, pz, ee, mass );
76  // other
77  else {
78  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
79  int eventSize = (fMasterGen->event).size()-1;
80  // -log(flat) = exponential distribution
81  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
82  (fMasterGen->event)[eventSize].tau( tauTmp );
83  }
84 
85 // Here also need to add anti-particle (if any)
86 // otherwise just add a 2nd particle of the same type
87 // (for example, gamma)
88 //
89  if ( fAddAntiParticle ) {
90  if( 1 <= std::abs(particleID) && std::abs(particleID) <= 6){ // quarks
91  (fMasterGen->event).append( -particleID, 23, 0, 101, -px, -py, -pz, ee, mass );
92  } else if (std::abs(particleID) == 21){ // gluons
93  (fMasterGen->event).append( 21, 23, 102, 101, -px, -py, -pz, ee, mass );
94  } else {
95  if ( (fMasterGen->particleData).isParticle( -particleID ) ) {
96  (fMasterGen->event).append( -particleID, 1, 0, 0, -px, -py, -pz, ee, mass );
97  } else {
98  (fMasterGen->event).append( particleID, 1, 0, 0, -px, -py, -pz, ee, mass );
99  }
100  int eventSize = (fMasterGen->event).size()-1;
101  // -log(flat) = exponential distribution
102  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
103  (fMasterGen->event)[eventSize].tau( tauTmp );
104  }
105  }
106  }
107 
108  if ( !fMasterGen->next() ) return false;
109  evtGenDecay();
110 
111  event().reset(new HepMC::GenEvent);
112  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
113 
114 }
size
Write out results.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:60
double flat() override
Definition: P8RndmEngine.cc:7
bool fAddAntiParticle
Definition: Py8PtGun.cc:25
void evtGenDecay()
Definition: Py8GunBase.cc:159
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< HepMC::GenEvent > & event()
std::vector< int > fPartIDs
Definition: Py8GunBase.h:59
P8RndmEngine & randomEngine()
double fMinPt
Definition: Py8PtGun.cc:23
double fMinEta
Definition: Py8PtGun.cc:21
double fMaxPhi
Definition: Py8GunBase.h:61
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
double fMaxPt
Definition: Py8PtGun.cc:24
double fMaxEta
Definition: Py8PtGun.cc:22
Definition: event.py:1

Member Data Documentation

bool gen::Py8PtGun::fAddAntiParticle
private

Definition at line 25 of file Py8PtGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtGun().

double gen::Py8PtGun::fMaxEta
private

Definition at line 22 of file Py8PtGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtGun().

double gen::Py8PtGun::fMaxPt
private

Definition at line 24 of file Py8PtGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtGun().

double gen::Py8PtGun::fMinEta
private

Definition at line 21 of file Py8PtGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtGun().

double gen::Py8PtGun::fMinPt
private

Definition at line 23 of file Py8PtGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtGun().