CMS 3D CMS Logo

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