CMS 3D CMS Logo

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

Public Member Functions

const char * classname () const override
 
bool generatePartonsAndHadronize () override
 
 Py8PtAndLxyGun (edm::ParameterSet const &)
 
 ~Py8PtAndLxyGun () 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

bool fAddAntiParticle
 
double fConeH
 
double fConeRadius
 
double fDistanceToAPEX
 
double fDxyMax
 
double fDzMax
 
double fLxyBackFraction
 
double fLxyMax
 
double fLxyMin
 
double fLzMax
 
double fLzOppositeFraction
 
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< 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 Py8PtAndLxyGun.cc.

Constructor & Destructor Documentation

◆ Py8PtAndLxyGun()

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

Definition at line 40 of file Py8PtAndLxyGun.cc.

References fAddAntiParticle, fConeH, fConeRadius, fDistanceToAPEX, fDxyMax, fDzMax, fLxyBackFraction, fLxyMax, fLxyMin, fLzMax, fLzOppositeFraction, fMaxEta, fMaxPt, fMinEta, fMinPt, and edm::ParameterSet::getParameter().

40  : Py8GunBase(ps) {
41  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters");
42  fMinEta = pgun_params.getParameter<double>("MinEta");
43  fMaxEta = pgun_params.getParameter<double>("MaxEta");
44  fMinPt = pgun_params.getParameter<double>("MinPt");
45  fMaxPt = pgun_params.getParameter<double>("MaxPt");
46  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle");
47  fDxyMax = pgun_params.getParameter<double>("dxyMax");
48  fDzMax = pgun_params.getParameter<double>("dzMax");
49  fLxyMin = pgun_params.getParameter<double>("LxyMin");
50  fLxyMax = pgun_params.getParameter<double>("LxyMax");
51  fLzMax = pgun_params.getParameter<double>("LzMax");
52  fConeRadius = pgun_params.getParameter<double>("ConeRadius");
53  fConeH = pgun_params.getParameter<double>("ConeH");
54  fDistanceToAPEX = pgun_params.getParameter<double>("DistanceToAPEX");
55  fLxyBackFraction = std::clamp(pgun_params.getParameter<double>("LxyBackFraction"), 0., 1.);
56  fLzOppositeFraction = std::clamp(pgun_params.getParameter<double>("LzOppositeFraction"), 0., 1.);
57  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Py8GunBase(edm::ParameterSet const &ps)
Definition: Py8GunBase.cc:15

◆ ~Py8PtAndLxyGun()

gen::Py8PtAndLxyGun::~Py8PtAndLxyGun ( )
inlineoverride

Definition at line 14 of file Py8PtAndLxyGun.cc.

14 {}

Member Function Documentation

◆ classname()

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

Implements gen::Py8InterfaceBase.

Definition at line 199 of file Py8PtAndLxyGun.cc.

199 { return "Py8PtAndLxyGun"; }

◆ generatePartonsAndHadronize()

bool gen::Py8PtAndLxyGun::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 59 of file Py8PtAndLxyGun.cc.

References funct::abs(), mps_setup::append, funct::cos(), PVValHelper::dxy, PVValHelper::dz, PVValHelper::eta, gen::BaseHadronizer::event(), gen::Py8GunBase::evtGenDecay(), JetChargeProducer_cfi::exp, fAddAntiParticle, fConeH, fConeRadius, fDistanceToAPEX, fDxyMax, fDzMax, gen::P8RndmEngine::flat(), fLxyBackFraction, fLxyMax, fLxyMin, fLzMax, fLzOppositeFraction, gen::Py8InterfaceBase::fMasterGen, fMaxEta, gen::Py8GunBase::fMaxPhi, fMaxPt, fMinEta, gen::Py8GunBase::fMinPhi, fMinPt, gen::Py8GunBase::fPartIDs, mps_fire::i, dqmiolumiharvest::j, dqm-mbProfile::log, M_PI, EgHLTOffHistBins_cfi::mass, EgammaObjectsElectrons_cfi::particleID, DiDispStaMuonMonitor_cfi::pt, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, gen::Py8InterfaceBase::randomEngine(), Validation_hcalonly_cfi::sign, funct::sin(), findQualityFiles::size, mathSSE::sqrt(), funct::tan(), metsig::tau, theta(), protons_cff::time, and gen::Py8InterfaceBase::toHepMC.

59  {
60  fMasterGen->event.reset();
61 
62  for (size_t i = 0; i < fPartIDs.size(); i++) {
63  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
64 
65  double phi = 0;
66  double dxy = 0;
67  double pt = 0;
68  double eta = 0;
69  double px = 0;
70  double py = 0;
71  double pz = 0;
72  double mass = 0;
73  double ee = 0;
74  double vx = 0;
75  double vy = 0;
76  double vz = 0;
77  double lxy = 0;
78 
79  bool passLoop = false;
80  while (!passLoop) {
81  bool passDxy = false;
82  bool passLz = false;
83  bool passDz = false;
84 
86  pt = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
87  px = pt * cos(phi);
88  py = pt * sin(phi);
89 
90  lxy = (fLxyMax - fLxyMin) * randomEngine().flat() + fLxyMin;
91 
92  int sign = 1;
93  for (int i = 0; i < 10000; i++) {
94  double vphi = 2 * M_PI * randomEngine().flat();
95  vx = lxy * cos(vphi);
96  vy = lxy * sin(vphi);
97 
98  dxy = -vx * sin(phi) + vy * cos(phi);
99 
100  sign = 1;
101  if (fLxyBackFraction > 0 && randomEngine().flat() <= fLxyBackFraction) {
102  sign = -1;
103  }
104  if ((std::abs(dxy) < fDxyMax || fDxyMax < 0) && sign * (vx * px + vy * py) > 0) {
105  passDxy = true;
106  break;
107  }
108  }
109 
110  eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
111  double theta = 2. * atan(exp(-eta));
112 
113  mass = (fMasterGen->particleData).m0(particleID);
114 
115  double pp = pt / sin(theta); // sqrt( ee*ee - mass*mass );
116  ee = sqrt(pp * pp + mass * mass);
117 
118  pz = pp * cos(theta);
119 
120  float coneTheta = fConeRadius / fConeH;
121  for (int j = 0; j < 100; j++) {
122  vz = fLzMax * randomEngine().flat(); // this is abs(vz)
123  float v0 = vz - fDistanceToAPEX;
124  if (v0 <= 0 || lxy * lxy / (coneTheta * coneTheta) > v0 * v0) {
125  passLz = true;
126  break;
127  }
128  }
129 
131  sign *= -1;
132  if (sign * pz < 0)
133  vz = -vz;
134 
135  double dz = vz - (vx * cos(phi) + vy * sin(phi)) / tan(theta);
136  if (std::abs(dz) < fDzMax || fDzMax < 0) {
137  passDz = true;
138  }
139 
140  passLoop = (passDxy && passLz && passDz);
141  if (passLoop)
142  break;
143  }
144 
145  float time = sqrt(vx * vx + vy * vy + vz * vz);
146 
147  if (!((fMasterGen->particleData).isParticle(particleID))) {
149  }
150  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
151  (fMasterGen->event).append(particleID, 23, 101, 0, px, py, pz, ee, mass);
152  else if (std::abs(particleID) == 21) // gluons
153  (fMasterGen->event).append(21, 23, 101, 102, px, py, pz, ee, mass);
154  // other
155  else {
156  (fMasterGen->event).append(particleID, 1, 0, 0, px, py, pz, ee, mass);
157  int eventSize = (fMasterGen->event).size() - 1;
158  // -log(flat) = exponential distribution
159  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
160  (fMasterGen->event)[eventSize].tau(tauTmp);
161  }
162  (fMasterGen->event).back().vProd(vx, vy, vz, time);
163 
164  // Here also need to add anti-particle (if any)
165  // otherwise just add a 2nd particle of the same type
166  // (for example, gamma).
167  // Added anti-particle has momentum opposite to corresponding
168  // particle, (px,py,pz)=>(-px,-py,-pz), and production vertex
169  // symmetric wrt (0,0,0), (vx, vy, vz)=>(-vx, -vy, -vz).
170  //
171  if (fAddAntiParticle) {
172  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
173  (fMasterGen->event).append(-particleID, 23, 0, 101, -px, -py, -pz, ee, mass);
174  } else if (std::abs(particleID) == 21) { // gluons
175  (fMasterGen->event).append(21, 23, 102, 101, -px, -py, -pz, ee, mass);
176  } else {
177  if ((fMasterGen->particleData).isParticle(-particleID)) {
178  (fMasterGen->event).append(-particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
179  } else {
180  (fMasterGen->event).append(particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
181  }
182  int eventSize = (fMasterGen->event).size() - 1;
183  // -log(flat) = exponential distribution
184  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
185  (fMasterGen->event)[eventSize].tau(tauTmp);
186  }
187  (fMasterGen->event).back().vProd(-vx, -vy, -vz, time);
188  }
189  }
190 
191  if (!fMasterGen->next())
192  return false;
193  evtGenDecay();
194 
195  event() = std::make_unique<HepMC::GenEvent>();
196  return toHepMC.fill_next_event(fMasterGen->event, event().get());
197  }
size
Write out results.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double flat() override
Definition: P8RndmEngine.cc:7
void evtGenDecay()
Definition: Py8GunBase.cc:144
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.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
#define M_PI
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Geom::Theta< T > theta() const
Definition: event.py:1

Member Data Documentation

◆ fAddAntiParticle

bool gen::Py8PtAndLxyGun::fAddAntiParticle
private

Definition at line 25 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fConeH

double gen::Py8PtAndLxyGun::fConeH
private

Definition at line 32 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fConeRadius

double gen::Py8PtAndLxyGun::fConeRadius
private

Definition at line 31 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fDistanceToAPEX

double gen::Py8PtAndLxyGun::fDistanceToAPEX
private

Definition at line 33 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fDxyMax

double gen::Py8PtAndLxyGun::fDxyMax
private

Definition at line 26 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fDzMax

double gen::Py8PtAndLxyGun::fDzMax
private

Definition at line 27 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fLxyBackFraction

double gen::Py8PtAndLxyGun::fLxyBackFraction
private

Definition at line 34 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fLxyMax

double gen::Py8PtAndLxyGun::fLxyMax
private

Definition at line 29 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fLxyMin

double gen::Py8PtAndLxyGun::fLxyMin
private

Definition at line 28 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fLzMax

double gen::Py8PtAndLxyGun::fLzMax
private

Definition at line 30 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fLzOppositeFraction

double gen::Py8PtAndLxyGun::fLzOppositeFraction
private

Definition at line 35 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fMaxEta

double gen::Py8PtAndLxyGun::fMaxEta
private

Definition at line 22 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fMaxPt

double gen::Py8PtAndLxyGun::fMaxPt
private

Definition at line 24 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fMinEta

double gen::Py8PtAndLxyGun::fMinEta
private

Definition at line 21 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().

◆ fMinPt

double gen::Py8PtAndLxyGun::fMinPt
private

Definition at line 23 of file Py8PtAndLxyGun.cc.

Referenced by generatePartonsAndHadronize(), and Py8PtAndLxyGun().