CMS 3D CMS Logo

Py8EGun.cc
Go to the documentation of this file.
1 
2 #include <memory>
3 
6 
8 
9 namespace gen {
10 
11  class Py8EGun : public Py8GunBase {
12  public:
13  Py8EGun(edm::ParameterSet const&);
14  ~Py8EGun() override {}
15 
16  bool generatePartonsAndHadronize() override;
17  const char* classname() const override;
18 
19  private:
20  // EGun particle(s) characteristics
21  double fMinEta;
22  double fMaxEta;
23  double fMinE;
24  double fMaxE;
26  };
27 
28  // implementation
29  //
31  // ParameterSet defpset ;
32  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
33  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
34  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
35  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
36  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 0.);
37  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
38  }
39 
41  fMasterGen->event.reset();
42 
43  int NTotParticles = fPartIDs.size();
44  if (fAddAntiParticle)
45  NTotParticles *= 2;
46 
47  // energy below is dummy, it is not used
48  (fMasterGen->event).append(990, -11, 0, 0, 2, 1 + NTotParticles, 0, 0, 0., 0., 0., 15000., 15000.);
49 
50  int colorindex = 101;
51 
52  for (size_t i = 0; i < fPartIDs.size(); i++) {
53  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
54  if ((std::abs(particleID) <= 6 || particleID == 21) && !(fAddAntiParticle)) {
55  throw cms::Exception("PythiaError") << "Attempting to generate quarks or gluons without setting "
56  "AddAntiParticle to true. This will not handle color properly."
57  << std::endl;
58  }
59 
60  double phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
61  double ee = (fMaxE - fMinE) * randomEngine().flat() + fMinE;
62  double eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
63  double the = 2. * atan(exp(-eta));
64 
65  double mass = (fMasterGen->particleData).m0(particleID);
66 
67  double pp = sqrt(ee * ee - mass * mass);
68  double px = pp * sin(the) * cos(phi);
69  double py = pp * sin(the) * sin(phi);
70  double pz = pp * cos(the);
71 
72  if (!((fMasterGen->particleData).isParticle(particleID))) {
73  particleID = std::fabs(particleID);
74  }
75  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
76  (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, colorindex, 0, px, py, pz, ee, mass);
77  if (!fAddAntiParticle)
78  colorindex += 1;
79  } else if (std::abs(particleID) == 21) { // gluons
80  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, colorindex, colorindex + 1, px, py, pz, ee, mass);
81  if (!fAddAntiParticle) {
82  colorindex += 2;
83  }
84  }
85  // other
86  else {
87  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
88  int eventSize = (fMasterGen->event).size() - 1;
89  // -log(flat) = exponential distribution
90  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
91  (fMasterGen->event)[eventSize].tau(tauTmp);
92  }
93 
94  // Here also need to add anti-particle (if any)
95  // otherwise just add a 2nd particle of the same type
96  // (for example, gamma)
97  //
98  if (fAddAntiParticle) {
99  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
100  (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, colorindex, -px, -py, -pz, ee, mass);
101  colorindex += 1;
102  } else if (std::abs(particleID) == 21) { // gluons
103  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, colorindex + 1, colorindex, -px, -py, -pz, ee, mass);
104  colorindex += 2;
105  } else {
106  if ((fMasterGen->particleData).isParticle(-particleID)) {
107  (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
108  } else {
109  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
110  }
111  int eventSize = (fMasterGen->event).size() - 1;
112  // -log(flat) = exponential distribution
113  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
114  (fMasterGen->event)[eventSize].tau(tauTmp);
115  }
116  }
117  }
118 
119  if (!fMasterGen->next())
120  return false;
121  evtGenDecay();
122 
123  event() = std::make_unique<HepMC::GenEvent>();
124  return toHepMC.fill_next_event(fMasterGen->event, event().get());
125  }
126 
127  const char* Py8EGun::classname() const { return "Py8EGun"; }
128 
130 
131 } // namespace gen
132 
133 using gen::Pythia8EGun;
size
Write out results.
bool fAddAntiParticle
Definition: Py8EGun.cc:25
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
double fMinE
Definition: Py8EGun.cc:23
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double fMinEta
Definition: Py8EGun.cc:21
double flat() override
Definition: P8RndmEngine.cc:7
void evtGenDecay()
Definition: Py8GunBase.cc:144
T sqrt(T t)
Definition: SSEVec.h:19
double fMaxEta
Definition: Py8EGun.cc:22
bool generatePartonsAndHadronize() override
Definition: Py8EGun.cc:40
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
edm::GeneratorFilter< gen::Py8EGun, gen::ExternalDecayDriver > Pythia8EGun
Definition: Py8EGun.cc:129
std::unique_ptr< HepMC::GenEvent > & event()
std::vector< int > fPartIDs
Definition: Py8GunBase.h:57
const char * classname() const override
Definition: Py8EGun.cc:127
P8RndmEngine & randomEngine()
double fMaxE
Definition: Py8EGun.cc:24
double fMaxPhi
Definition: Py8GunBase.h:59
~Py8EGun() override
Definition: Py8EGun.cc:14
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
Py8EGun(edm::ParameterSet const &)
Definition: Py8EGun.cc:30
Definition: event.py:1