CMS 3D CMS Logo

Py8PtAndLxyGun.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <algorithm>
3 
6 
8 
9 namespace gen {
10 
11  class Py8PtAndLxyGun : public Py8GunBase {
12  public:
14  ~Py8PtAndLxyGun() override {}
15 
16  bool generatePartonsAndHadronize() override;
17  const char* classname() const override;
18 
19  private:
20  // PtAndLxyGun particle(s) characteristics
21  double fMinEta;
22  double fMaxEta;
23  double fMinPt;
24  double fMaxPt;
26  double fDxyMax;
27  double fDzMax;
28  double fLxyMin;
29  double fLxyMax;
30  double fLzMax;
31  double fConeRadius;
32  double fConeH;
36  };
37 
38  // implementation
39  //
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  }
58 
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 
85  phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
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  }
198 
199  const char* Py8PtAndLxyGun::classname() const { return "Py8PtAndLxyGun"; }
200 
202 
203 } // namespace gen
204 
size
Write out results.
edm::GeneratorFilter< gen::Py8PtAndLxyGun, gen::ExternalDecayDriver > Pythia8PtAndLxyGun
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
~Py8PtAndLxyGun() override
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
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::unique_ptr< HepMC::GenEvent > & event()
std::vector< int > fPartIDs
Definition: Py8GunBase.h:57
#define M_PI
const char * classname() const override
P8RndmEngine & randomEngine()
bool generatePartonsAndHadronize() override
Py8PtAndLxyGun(edm::ParameterSet const &)
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Geom::Theta< T > theta() const
Definition: event.py:1