CMS 3D CMS Logo

Py8PtAndDxyGun.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
7 
8 namespace gen {
9 
10  class Py8PtAndDxyGun : public Py8GunBase {
11  public:
13  ~Py8PtAndDxyGun() override {}
14 
15  bool generatePartonsAndHadronize() override;
16  const char* classname() const override;
17 
18  private:
19  // PtAndDxyGun particle(s) characteristics
20  double fMinEta;
21  double fMaxEta;
22  double fMinPt;
23  double fMaxPt;
25  double fDxyMin;
26  double fDxyMax;
27  double fLxyMax;
28  double fLzMax;
29  double fConeRadius;
30  double fConeH;
32  };
33 
34  // implementation
35  //
37  // ParameterSet defpset ;
38  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
39  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
40  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
41  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
42  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 0.);
43  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
44  fDxyMin = pgun_params.getParameter<double>("dxyMin");
45  fDxyMax = pgun_params.getParameter<double>("dxyMax");
46  fLxyMax = pgun_params.getParameter<double>("LxyMax");
47  fLzMax = pgun_params.getParameter<double>("LzMax");
48  fConeRadius = pgun_params.getParameter<double>("ConeRadius");
49  fConeH = pgun_params.getParameter<double>("ConeH");
50  fDistanceToAPEX = pgun_params.getParameter<double>("DistanceToAPEX");
51  }
52 
54  fMasterGen->event.reset();
55 
56  for (size_t i = 0; i < fPartIDs.size(); i++) {
57  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
58 
59  double phi = 0;
60  double dxy = 0;
61  double pt = 0;
62  double eta = 0;
63  double px = 0;
64  double py = 0;
65  double pz = 0;
66  double mass = 0;
67  double ee = 0;
68  double vx = 0;
69  double vy = 0;
70  double vz = 0;
71  double lxy = 0;
72 
73  bool passLoop = false;
74  while (!passLoop) {
75  bool passLxy = false;
76  bool passLz = false;
77 
78  phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
79  dxy = (fDxyMax - fDxyMin) * randomEngine().flat() + fDxyMin;
80  float dxysign = randomEngine().flat() - 0.5;
81  if (dxysign < 0)
82  dxy = -dxy;
83 
84  pt = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
85  px = pt * cos(phi);
86  py = pt * sin(phi);
87 
88  for (int i = 0; i < 10000; i++) {
89  vx = 2 * fLxyMax * randomEngine().flat() - fLxyMax;
90  vy = (pt * dxy + vx * py) / px;
91  lxy = sqrt(vx * vx + vy * vy);
92  if (lxy < std::abs(fLxyMax) && (vx * px + vy * py) > 0) {
93  passLxy = true;
94  break;
95  }
96  }
97 
98  eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
99  double theta = 2. * atan(exp(-eta));
100 
101  mass = (fMasterGen->particleData).m0(particleID);
102 
103  double pp = pt / sin(theta); // sqrt( ee*ee - mass*mass );
104  ee = sqrt(pp * pp + mass * mass);
105 
106  pz = pp * cos(theta);
107 
108  float coneTheta = fConeRadius / fConeH;
109  for (int j = 0; j < 100; j++) {
110  vz = fLzMax * randomEngine().flat(); // this is abs(vz)
111  float v0 = vz - fDistanceToAPEX;
112  if (v0 <= 0 || lxy * lxy / (coneTheta * coneTheta) > v0 * v0) {
113  passLz = true;
114  break;
115  }
116  }
117  if (pz < 0)
118  vz = -vz;
119  passLoop = (passLxy && passLz);
120 
121  if (passLoop)
122  break;
123  }
124 
125  float time = sqrt(vx * vx + vy * vy + vz * vz);
126 
127  if (!((fMasterGen->particleData).isParticle(particleID))) {
128  particleID = std::abs(particleID);
129  }
130  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
131  (fMasterGen->event).append(particleID, 23, 101, 0, px, py, pz, ee, mass);
132  else if (std::abs(particleID) == 21) // gluons
133  (fMasterGen->event).append(21, 23, 101, 102, px, py, pz, ee, mass);
134  // other
135  else {
136  (fMasterGen->event).append(particleID, 1, 0, 0, px, py, pz, ee, mass);
137  int eventSize = (fMasterGen->event).size() - 1;
138  // -log(flat) = exponential distribution
139  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
140  (fMasterGen->event)[eventSize].tau(tauTmp);
141  }
142  (fMasterGen->event).back().vProd(vx, vy, vz, time);
143 
144  // Here also need to add anti-particle (if any)
145  // otherwise just add a 2nd particle of the same type
146  // (for example, gamma)
147  //
148  if (fAddAntiParticle) {
149  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
150  (fMasterGen->event).append(-particleID, 23, 0, 101, -px, -py, -pz, ee, mass);
151  } else if (std::abs(particleID) == 21) { // gluons
152  (fMasterGen->event).append(21, 23, 102, 101, -px, -py, -pz, ee, mass);
153  } else {
154  if ((fMasterGen->particleData).isParticle(-particleID)) {
155  (fMasterGen->event).append(-particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
156  } else {
157  (fMasterGen->event).append(particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
158  }
159  int eventSize = (fMasterGen->event).size() - 1;
160  // -log(flat) = exponential distribution
161  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
162  (fMasterGen->event)[eventSize].tau(tauTmp);
163  }
164  (fMasterGen->event).back().vProd(-vx, -vy, -vz, time);
165  }
166  }
167 
168  if (!fMasterGen->next())
169  return false;
170  evtGenDecay();
171 
172  event() = std::make_unique<HepMC::GenEvent>();
173  return toHepMC.fill_next_event(fMasterGen->event, event().get());
174  }
175 
176  const char* Py8PtAndDxyGun::classname() const { return "Py8PtAndDxyGun"; }
177 
179 
180 } // namespace gen
181 
size
Write out results.
T getParameter(std::string const &) const
edm::GeneratorFilter< gen::Py8PtAndDxyGun, gen::ExternalDecayDriver > Pythia8PtAndDxyGun
const char * classname() const override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:60
Geom::Theta< T > theta() const
double flat() override
Definition: P8RndmEngine.cc:7
bool generatePartonsAndHadronize() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void evtGenDecay()
Definition: Py8GunBase.cc:159
Py8PtAndDxyGun(edm::ParameterSet const &)
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
~Py8PtAndDxyGun() override
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:61
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: event.py:1