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  int NTotParticles = fPartIDs.size();
57  if (fAddAntiParticle)
58  NTotParticles *= 2;
59 
60  // energy below is dummy, it is not used
61  (fMasterGen->event).append(990, -11, 0, 0, 2, 1 + NTotParticles, 0, 0, 0., 0., 0., 15000., 15000.);
62 
63  for (size_t i = 0; i < fPartIDs.size(); i++) {
64  int particleID = fPartIDs[i]; // this is PDG
65 
66  double phi = 0;
67  double dxy = 0;
68  double pt = 0;
69  double eta = 0;
70  double px = 0;
71  double py = 0;
72  double pz = 0;
73  double mass = 0;
74  double ee = 0;
75  double vx = 0;
76  double vy = 0;
77  double vz = 0;
78  double lxy = 0;
79 
80  bool passLoop = false;
81  while (!passLoop) {
82  bool passLxy = false;
83  bool passLz = false;
84 
85  phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
87  float dxysign = randomEngine().flat() - 0.5;
88  if (dxysign < 0)
89  dxy = -dxy;
90 
91  pt = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
92  px = pt * cos(phi);
93  py = pt * sin(phi);
94 
95  for (int i = 0; i < 10000; i++) {
96  vx = 2 * fLxyMax * randomEngine().flat() - fLxyMax;
97  vy = (pt * dxy + vx * py) / px;
98  lxy = sqrt(vx * vx + vy * vy);
99  if (lxy < std::abs(fLxyMax) && (vx * px + vy * py) > 0) {
100  passLxy = true;
101  break;
102  }
103  }
104 
105  eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
106  double theta = 2. * atan(exp(-eta));
107 
108  mass = (fMasterGen->particleData).m0(particleID);
109 
110  double pp = pt / sin(theta); // sqrt( ee*ee - mass*mass );
111  ee = sqrt(pp * pp + mass * mass);
112 
113  pz = pp * cos(theta);
114 
115  float coneTheta = fConeRadius / fConeH;
116  for (int j = 0; j < 100; j++) {
117  vz = fLzMax * randomEngine().flat(); // this is abs(vz)
118  float v0 = vz - fDistanceToAPEX;
119  if (v0 <= 0 || lxy * lxy / (coneTheta * coneTheta) > v0 * v0) {
120  passLz = true;
121  break;
122  }
123  }
124  if (pz < 0)
125  vz = -vz;
126  passLoop = (passLxy && passLz);
127 
128  if (passLoop)
129  break;
130  }
131 
132  float time = sqrt(vx * vx + vy * vy + vz * vz);
133 
134  if (!((fMasterGen->particleData).isParticle(particleID))) {
136  }
137  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
138  (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, 101, 0, px, py, pz, ee, mass);
139  else if (std::abs(particleID) == 21) // gluons
140  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 101, 102, px, py, pz, ee, mass);
141  // other
142  else {
143  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
144  int eventSize = (fMasterGen->event).size() - 1;
145  // -log(flat) = exponential distribution
146  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
147  (fMasterGen->event)[eventSize].tau(tauTmp);
148  }
149  (fMasterGen->event).back().vProd(vx, vy, vz, time);
150 
151  // Here also need to add anti-particle (if any)
152  // otherwise just add a 2nd particle of the same type
153  // (for example, gamma)
154  //
155  if (fAddAntiParticle) {
156  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
157  (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, 101, -px, -py, -pz, ee, mass);
158  } else if (std::abs(particleID) == 21) { // gluons
159  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 102, 101, -px, -py, -pz, ee, mass);
160  } else {
161  if ((fMasterGen->particleData).isParticle(-particleID)) {
162  (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
163  } else {
164  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
165  }
166  int eventSize = (fMasterGen->event).size() - 1;
167  // -log(flat) = exponential distribution
168  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
169  (fMasterGen->event)[eventSize].tau(tauTmp);
170  }
171  (fMasterGen->event).back().vProd(-vx, -vy, -vz, time);
172  }
173  }
174 
175  if (!fMasterGen->next())
176  return false;
177  evtGenDecay();
178 
179  event() = std::make_unique<HepMC::GenEvent>();
180  return toHepMC.fill_next_event(fMasterGen->event, event().get());
181  }
182 
183  const char* Py8PtAndDxyGun::classname() const { return "Py8PtAndDxyGun"; }
184 
186 
187 } // namespace gen
188 
size
Write out results.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::GeneratorFilter< gen::Py8PtAndDxyGun, gen::ExternalDecayDriver > Pythia8PtAndDxyGun
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double flat() override
Definition: P8RndmEngine.cc:7
bool generatePartonsAndHadronize() override
void evtGenDecay()
Definition: Py8GunBase.cc:144
Py8PtAndDxyGun(edm::ParameterSet const &)
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.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
~Py8PtAndDxyGun() override
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
const char * classname() const override
Geom::Theta< T > theta() const
Definition: event.py:1