CMS 3D CMS Logo

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

Public Member Functions

const char * classname () const override
 
bool generatePartonsAndHadronize () override
 
 Py8MassGun (edm::ParameterSet const &)
 
 ~Py8MassGun () 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< HepMC3::GenEvent > getGenEvent3 ()
 
std::unique_ptr< GenEventInfoProductgetGenEventInfo ()
 
std::unique_ptr< GenEventInfoProduct3getGenEventInfo3 ()
 
virtual std::unique_ptr< GenLumiInfoHeadergetGenLumiInfoHeader () const
 
GenRunInfoProductgetGenRunInfo ()
 
std::unique_ptr< lhef::LHEEventgetLHEEvent ()
 
const std::shared_ptr< lhef::LHERunInfo > & getLHERunInfo () const
 
unsigned int getVHepMC ()
 
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 resetEvent3 (std::unique_ptr< HepMC3::GenEvent > event3)
 
void resetEventInfo (std::unique_ptr< GenEventInfoProduct > eventInfo)
 
void resetEventInfo3 (std::unique_ptr< GenEventInfoProduct3 > 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

double fMaxEta
 
double fMaxM
 
double fMaxP
 
double fMaxPt
 
double fMinEta
 
double fMinM
 
double fMinP
 
double fMinPt
 
int fMomMode
 

Additional Inherited Members

- Protected Member Functions inherited from gen::BaseHadronizer
std::unique_ptr< HepMC::GenEvent > & event ()
 
std::unique_ptr< HepMC3::GenEvent > & event3 ()
 
std::unique_ptr< GenEventInfoProduct > & eventInfo ()
 
std::unique_ptr< GenEventInfoProduct3 > & eventInfo3 ()
 
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
unsigned int ivhepmc = 2
 
std::string lheFile_
 
int randomIndex_
 

Detailed Description

Definition at line 11 of file Py8MassGun.cc.

Constructor & Destructor Documentation

◆ Py8MassGun()

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

Definition at line 34 of file Py8MassGun.cc.

References fMaxEta, fMaxM, fMaxP, fMaxPt, fMinEta, fMinM, fMinP, fMinPt, fMomMode, and edm::ParameterSet::getParameter().

34  : Py8GunBase(ps) {
35  // ParameterSet defpset ;
36  edm::ParameterSet pgun_params = 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  fMinP = pgun_params.getParameter<double>("MinP"); // , 0.);
40  fMaxP = pgun_params.getParameter<double>("MaxP"); // , 0.);
41  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
42  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 0.);
43  fMinM = pgun_params.getParameter<double>("MinM"); // , 0.);
44  fMaxM = pgun_params.getParameter<double>("MaxM"); // , 0.);
45  fMomMode = pgun_params.getParameter<int>("MomMode"); // , 1);
46  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Py8GunBase(edm::ParameterSet const &ps)
Definition: Py8GunBase.cc:15

◆ ~Py8MassGun()

gen::Py8MassGun::~Py8MassGun ( )
inlineoverride

Definition at line 14 of file Py8MassGun.cc.

14 {}

Member Function Documentation

◆ classname()

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

Implements gen::Py8InterfaceBase.

Definition at line 120 of file Py8MassGun.cc.

120 { return "Py8MassGun"; }

◆ generatePartonsAndHadronize()

bool gen::Py8MassGun::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 48 of file Py8MassGun.cc.

References mps_setup::append, funct::cos(), StorageManager_cfg::e1, PVValHelper::eta, gen::BaseHadronizer::event(), JetChargeProducer_cfi::exp, gen::P8RndmEngine::flat(), gen::Py8InterfaceBase::fMasterGen, fMaxEta, fMaxM, fMaxP, gen::Py8GunBase::fMaxPhi, fMaxPt, fMinEta, fMinM, fMinP, gen::Py8GunBase::fMinPhi, fMinPt, fMomMode, gen::Py8GunBase::fPartIDs, callgraph::m2, M_PI, pv::pT, gen::Py8InterfaceBase::randomEngine(), funct::sin(), mathSSE::sqrt(), and gen::Py8InterfaceBase::toHepMC.

48  {
49  fMasterGen->event.reset();
50  size_t pSize = fPartIDs.size();
51  if (pSize > 2)
52  return false;
53 
54  int NTotParticles = fPartIDs.size();
55 
56  // energy below is dummy, it is not used
57  (fMasterGen->event).append(990, -11, 0, 0, 2, 1 + NTotParticles, 0, 0, 0., 0., 0., 15000., 15000.);
58 
59  // Pick a flat mass range
60  double phi, eta, the, ee, pp;
61  double m0 = (fMaxM - fMinM) * randomEngine().flat() + fMinM;
62  // Global eta
64 
65  if (pSize == 2) {
66  // Masses.
67  double m1 = fMasterGen->particleData.m0(fPartIDs[0]);
68  double m2 = fMasterGen->particleData.m0(fPartIDs[1]);
69 
70  // Energies and absolute momentum in the rest frame.
71  if (m1 + m2 > m0)
72  return false;
73  double e1 = 0.5 * (m0 * m0 + m1 * m1 - m2 * m2) / m0;
74  double e2 = 0.5 * (m0 * m0 + m2 * m2 - m1 * m1) / m0;
75  double pAbs = 0.5 * sqrt((m0 - m1 - m2) * (m0 + m1 + m2) * (m0 + m1 - m2) * (m0 - m1 + m2)) / m0;
76  // Isotropic angles in rest frame give three-momentum.
77  double cosTheta = 2. * randomEngine().flat() - 1.;
78  double sinTheta = sqrt(1. - cosTheta * cosTheta);
79  phi = 2. * M_PI * randomEngine().flat();
80 
81  double pX = pAbs * sinTheta * cos(phi);
82  double pY = pAbs * sinTheta * sin(phi);
83  double pZ = pAbs * cosTheta;
84 
85  (fMasterGen->event).append(fPartIDs[0], 1, 1, 0, 0, 0, 0, 0, pX, pY, pZ, e1, m1);
86  (fMasterGen->event).append(fPartIDs[1], 1, 1, 0, 0, 0, 0, 0, -pX, -pY, -pZ, e2, m2);
87  } else {
88  (fMasterGen->event).append(fPartIDs[0], 1, 1, 0, 0, 0, 0, 0, 0.0, 0.0, 0.0, m0, m0);
89  }
90 
91  //now the boost (from input params)
92  if (fMomMode == 0) {
93  pp = (fMaxP - fMinP) * randomEngine().flat() + fMinP;
94  } else {
95  double pT = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
96  pp = pT * cosh(eta);
97  }
98  ee = sqrt(m0 * m0 + pp * pp);
99 
100  //the boost direction (from input params)
101  //
102  phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
103  the = 2. * atan(exp(-eta));
104 
105  double betaX = pp / ee * std::sin(the) * std::cos(phi);
106  double betaY = pp / ee * std::sin(the) * std::sin(phi);
107  double betaZ = pp / ee * std::cos(the);
108 
109  // boost all particles
110  //
111  (fMasterGen->event).bst(betaX, betaY, betaZ);
112 
113  if (!fMasterGen->next())
114  return false;
115 
116  event() = std::make_unique<HepMC::GenEvent>();
117  return toHepMC.fill_next_event(fMasterGen->event, event().get());
118  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double flat() override
Definition: P8RndmEngine.cc:7
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::unique_ptr< HepMC::GenEvent > & event()
std::vector< int > fPartIDs
Definition: Py8GunBase.h:57
#define M_PI
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: event.py:1

Member Data Documentation

◆ fMaxEta

double gen::Py8MassGun::fMaxEta
private

Definition at line 22 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMaxM

double gen::Py8MassGun::fMaxM
private

Definition at line 28 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMaxP

double gen::Py8MassGun::fMaxP
private

Definition at line 24 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMaxPt

double gen::Py8MassGun::fMaxPt
private

Definition at line 26 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMinEta

double gen::Py8MassGun::fMinEta
private

Definition at line 21 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMinM

double gen::Py8MassGun::fMinM
private

Definition at line 27 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMinP

double gen::Py8MassGun::fMinP
private

Definition at line 23 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMinPt

double gen::Py8MassGun::fMinPt
private

Definition at line 25 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().

◆ fMomMode

int gen::Py8MassGun::fMomMode
private

Definition at line 29 of file Py8MassGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8MassGun().