CMS 3D CMS Logo

Py8PtotGun.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
7 
8 namespace gen {
9 
10  class Py8PtotGun : public Py8GunBase {
11  public:
13  ~Py8PtotGun() override {}
14 
15  bool generatePartonsAndHadronize() override;
16  const char* classname() const override;
17 
18  private:
19  // Ptot Gun particle(s) characteristics
20  double fMinEta;
21  double fMaxEta;
22  double fMinPtot;
23  double fMaxPtot;
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  fMinPtot = pgun_params.getParameter<double>("MinPtot"); // , 0.);
35  fMaxPtot = pgun_params.getParameter<double>("MaxPtot"); // , 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
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 pp = (fMaxPtot - fMinPtot) * randomEngine().flat() + fMinPtot;
57 
58  double mass = (fMasterGen->particleData).m0(particleID);
59 
60  double pt = pp * sin(the);
61  double ee = sqrt(pp * pp + mass * mass);
62  double px = pt * cos(phi);
63  double py = pt * sin(phi);
64  double pz = pp * cos(the);
65 
66  if (!((fMasterGen->particleData).isParticle(particleID))) {
68  }
69  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
70  (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, 101, 0, px, py, pz, ee, mass);
71  else if (std::abs(particleID) == 21) // gluons
72  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 101, 102, px, py, pz, ee, mass);
73  // other
74  else {
75  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
76  int eventSize = (fMasterGen->event).size() - 1;
77  // -log(flat) = exponential distribution
78  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
79  (fMasterGen->event)[eventSize].tau(tauTmp);
80  }
81 
82  // Here also need to add anti-particle (if any)
83  // otherwise just add a 2nd particle of the same type
84  // (for example, gamma)
85  //
86  if (fAddAntiParticle) {
87  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
88  (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, 101, -px, -py, -pz, ee, mass);
89  } else if (std::abs(particleID) == 21) { // gluons
90  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 102, 101, -px, -py, -pz, ee, mass);
91  } else {
92  if ((fMasterGen->particleData).isParticle(-particleID)) {
93  (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
94  } else {
95  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
96  }
97  int eventSize = (fMasterGen->event).size() - 1;
98  // -log(flat) = exponential distribution
99  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
100  (fMasterGen->event)[eventSize].tau(tauTmp);
101  }
102  }
103  }
104 
105  if (!fMasterGen->next())
106  return false;
107  evtGenDecay();
108 
109  event() = std::make_unique<HepMC::GenEvent>();
110  return toHepMC.fill_next_event(fMasterGen->event, event().get());
111  }
112 
113  const char* Py8PtotGun::classname() const { return "Py8PtotGun"; }
114 
116 
117 } // namespace gen
118 
119 using gen::Pythia8PtotGun;
size
Write out results.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
~Py8PtotGun() override
Definition: Py8PtotGun.cc:13
edm::GeneratorFilter< gen::Py8PtotGun, gen::ExternalDecayDriver > Pythia8PtotGun
Definition: Py8PtotGun.cc:115
double flat() override
Definition: P8RndmEngine.cc:7
bool generatePartonsAndHadronize() override
Definition: Py8PtotGun.cc:39
bool fAddAntiParticle
Definition: Py8PtotGun.cc:24
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
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
Py8PtotGun(edm::ParameterSet const &)
Definition: Py8PtotGun.cc:29
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
const char * classname() const override
Definition: Py8PtotGun.cc:113
Definition: event.py:1