CMS 3D CMS Logo

Py8PtExpGun.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
7 
8 namespace gen {
9 
10  class Py8PtExpGun : public Py8GunBase {
11  public:
13  ~Py8PtExpGun() override {}
14 
15  bool generatePartonsAndHadronize() override;
16  const char* classname() const override;
17 
18  private:
19  // PtExpGun 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  //-log(flat) = exponential distribution
57  // need the /10.0 and the min with 1.0 to make sure pt doesn't go too high
58  // 10.0 chosen to give last pt bin (overflow bin) a reasonable (not unnaturally high) content
59  double pt = (std::min(-1 / 10.0 * log(randomEngine().flat()), 1.0)) * (fMaxPt - fMinPt) + fMinPt;
60 
61  double mass = (fMasterGen->particleData).m0(particleID);
62 
63  double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
64  double ee = sqrt(pp * pp + mass * mass);
65 
66  double px = pt * cos(phi);
67  double py = pt * sin(phi);
68  double pz = pp * cos(the);
69 
70  if (!((fMasterGen->particleData).isParticle(particleID))) {
72  }
73  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
74  (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, 101, 0, px, py, pz, ee, mass);
75  else if (std::abs(particleID) == 21) // gluons
76  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 101, 102, px, py, pz, ee, mass);
77  // other
78  else {
79  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
80  int eventSize = (fMasterGen->event).size() - 1;
81  // -log(flat) = exponential distribution
82  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
83  (fMasterGen->event)[eventSize].tau(tauTmp);
84  }
85 
86  // Here also need to add anti-particle (if any)
87  // otherwise just add a 2nd particle of the same type
88  // (for example, gamma)
89  //
90  if (fAddAntiParticle) {
91  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
92  (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, 101, -px, -py, -pz, ee, mass);
93  } else if (std::abs(particleID) == 21) { // gluons
94  (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 102, 101, -px, -py, -pz, ee, mass);
95  } else {
96  if ((fMasterGen->particleData).isParticle(-particleID)) {
97  (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
98  } else {
99  (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
100  }
101  int eventSize = (fMasterGen->event).size() - 1;
102  // -log(flat) = exponential distribution
103  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
104  (fMasterGen->event)[eventSize].tau(tauTmp);
105  }
106  }
107  }
108 
109  if (!fMasterGen->next())
110  return false;
111  evtGenDecay();
112 
113  event() = std::make_unique<HepMC::GenEvent>();
114  return toHepMC.fill_next_event(fMasterGen->event, event().get());
115  }
116 
117  const char* Py8PtExpGun::classname() const { return "Py8PtExpGun"; }
118 
120 
121 } // namespace gen
122 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const char * classname() const override
Definition: Py8PtExpGun.cc:117
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
double flat() override
Definition: P8RndmEngine.cc:7
Py8PtExpGun(edm::ParameterSet const &)
Definition: Py8PtExpGun.cc:29
void evtGenDecay()
Definition: Py8GunBase.cc:144
edm::GeneratorFilter< gen::Py8PtExpGun, gen::ExternalDecayDriver > Pythia8PtExpGun
Definition: Py8PtExpGun.cc:119
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
~Py8PtExpGun() override
Definition: Py8PtExpGun.cc:13
P8RndmEngine & randomEngine()
double fMaxPhi
Definition: Py8GunBase.h:59
HepMC::Pythia8ToHepMC toHepMC
std::unique_ptr< Pythia8::Pythia > fMasterGen
bool generatePartonsAndHadronize() override
Definition: Py8PtExpGun.cc:39
Definition: event.py:1