CMS 3D CMS Logo

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

Public Member Functions

const char * classname () const override
 
bool generatePartonsAndHadronize () override
 
 Py8EGun (edm::ParameterSet const &)
 
 ~Py8EGun () 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 makeTmpSLHA (const std::string &)
 
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 () noexcept(false)
 

Private Attributes

bool fAddAntiParticle
 
double fMaxE
 
double fMaxEta
 
double fMinE
 
double fMinEta
 

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< Pythia8::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 11 of file Py8EGun.cc.

Constructor & Destructor Documentation

◆ Py8EGun()

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

Definition at line 30 of file Py8EGun.cc.

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

30  : Py8GunBase(ps) {
31  // ParameterSet defpset ;
32  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
33  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
34  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
35  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
36  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 0.);
37  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
38  }
bool fAddAntiParticle
Definition: Py8EGun.cc:25
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double fMinE
Definition: Py8EGun.cc:23
double fMinEta
Definition: Py8EGun.cc:21
double fMaxEta
Definition: Py8EGun.cc:22
double fMaxE
Definition: Py8EGun.cc:24
Py8GunBase(edm::ParameterSet const &ps)
Definition: Py8GunBase.cc:15

◆ ~Py8EGun()

gen::Py8EGun::~Py8EGun ( )
inlineoverride

Definition at line 14 of file Py8EGun.cc.

14 {}

Member Function Documentation

◆ classname()

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

Implements gen::Py8InterfaceBase.

Definition at line 127 of file Py8EGun.cc.

127 { return "Py8EGun"; }

◆ generatePartonsAndHadronize()

bool gen::Py8EGun::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 40 of file Py8EGun.cc.

References funct::abs(), mps_setup::append, funct::cos(), PVValHelper::eta, gen::BaseHadronizer::event(), gen::Py8GunBase::evtGenDecay(), Exception, JetChargeProducer_cfi::exp, fAddAntiParticle, gen::P8RndmEngine::flat(), gen::Py8InterfaceBase::fMasterGen, fMaxE, fMaxEta, gen::Py8GunBase::fMaxPhi, fMinE, fMinEta, gen::Py8GunBase::fMinPhi, gen::Py8GunBase::fPartIDs, mps_fire::i, dqm-mbProfile::log, EgHLTOffHistBins_cfi::mass, EgammaObjectsElectrons_cfi::particleID, createTree::pp, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, gen::Py8InterfaceBase::randomEngine(), funct::sin(), findQualityFiles::size, mathSSE::sqrt(), metsig::tau, and gen::Py8InterfaceBase::toHepMC.

40  {
41  fMasterGen->event.reset();
42 
43  int NTotParticles = fPartIDs.size();
44  if (fAddAntiParticle)
45  NTotParticles *= 2;
46 
47  // energy below is dummy, it is not used
48  (fMasterGen->event).append(990, -11, 0, 0, 2, 1 + NTotParticles, 0, 0, 0., 0., 0., 15000., 15000.);
49 
50  int colorindex = 101;
51 
52  for (size_t i = 0; i < fPartIDs.size(); i++) {
53  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
54  if ((std::abs(particleID) <= 6 || particleID == 21) && !(fAddAntiParticle)) {
55  throw cms::Exception("PythiaError") << "Attempting to generate quarks or gluons without setting "
56  "AddAntiParticle to true. This will not handle color properly."
57  << std::endl;
58  }
59 
60  double phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
61  double ee = (fMaxE - fMinE) * randomEngine().flat() + fMinE;
62  double eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
63  double the = 2. * atan(exp(-eta));
64 
65  double mass = (fMasterGen->particleData).m0(particleID);
66 
67  double pp = sqrt(ee * ee - mass * mass);
68  double px = pp * sin(the) * cos(phi);
69  double py = pp * sin(the) * sin(phi);
70  double pz = pp * cos(the);
71 
72  if (!((fMasterGen->particleData).isParticle(particleID))) {
73  particleID = std::fabs(particleID);
74  }
75  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
76  (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, colorindex, 0, px, py, pz, ee, mass);
77  if (!fAddAntiParticle)
78  colorindex += 1;
79  } else if (std::abs(particleID) == 21) { // gluons
80  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, colorindex, colorindex + 1, px, py, pz, ee, mass);
81  if (!fAddAntiParticle) {
82  colorindex += 2;
83  }
84  }
85  // other
86  else {
87  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
88  int eventSize = (fMasterGen->event).size() - 1;
89  // -log(flat) = exponential distribution
90  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
91  (fMasterGen->event)[eventSize].tau(tauTmp);
92  }
93 
94  // Here also need to add anti-particle (if any)
95  // otherwise just add a 2nd particle of the same type
96  // (for example, gamma)
97  //
98  if (fAddAntiParticle) {
99  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
100  (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, colorindex, -px, -py, -pz, ee, mass);
101  colorindex += 1;
102  } else if (std::abs(particleID) == 21) { // gluons
103  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, colorindex + 1, colorindex, -px, -py, -pz, ee, mass);
104  colorindex += 2;
105  } else {
106  if ((fMasterGen->particleData).isParticle(-particleID)) {
107  (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
108  } else {
109  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
110  }
111  int eventSize = (fMasterGen->event).size() - 1;
112  // -log(flat) = exponential distribution
113  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
114  (fMasterGen->event)[eventSize].tau(tauTmp);
115  }
116  }
117  }
118 
119  if (!fMasterGen->next())
120  return false;
121  evtGenDecay();
122 
123  event() = std::make_unique<HepMC::GenEvent>();
124  return toHepMC.fill_next_event(fMasterGen->event, event().get());
125  }
size
Write out results.
bool fAddAntiParticle
Definition: Py8EGun.cc:25
double fMinE
Definition: Py8EGun.cc:23
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double fMinEta
Definition: Py8EGun.cc:21
double flat() override
Definition: P8RndmEngine.cc:7
void evtGenDecay()
Definition: Py8GunBase.cc:144
T sqrt(T t)
Definition: SSEVec.h:19
double fMaxEta
Definition: Py8EGun.cc:22
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:57
P8RndmEngine & randomEngine()
double fMaxE
Definition: Py8EGun.cc:24
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: event.py:1

Member Data Documentation

◆ fAddAntiParticle

bool gen::Py8EGun::fAddAntiParticle
private

Definition at line 25 of file Py8EGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8EGun().

◆ fMaxE

double gen::Py8EGun::fMaxE
private

Definition at line 24 of file Py8EGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8EGun().

◆ fMaxEta

double gen::Py8EGun::fMaxEta
private

Definition at line 22 of file Py8EGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8EGun().

◆ fMinE

double gen::Py8EGun::fMinE
private

Definition at line 23 of file Py8EGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8EGun().

◆ fMinEta

double gen::Py8EGun::fMinEta
private

Definition at line 21 of file Py8EGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8EGun().