CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  for (size_t i = 0; i < fPartIDs.size(); i++) {
43  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
44 
45  double phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
46  double eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
47  double the = 2. * atan(exp(-eta));
48 
49  double pt = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
50 
51  double mass = (fMasterGen->particleData).m0(particleID);
52 
53  double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
54  double ee = sqrt(pp * pp + mass * mass);
55 
56  double px = pt * cos(phi);
57  double py = pt * sin(phi);
58  double pz = pp * cos(the);
59 
60  if (!((fMasterGen->particleData).isParticle(particleID))) {
61  particleID = std::abs(particleID);
62  }
63  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
64  (fMasterGen->event).append(particleID, 23, 101, 0, px, py, pz, ee, mass);
65  else if (std::abs(particleID) == 21) // gluons
66  (fMasterGen->event).append(21, 23, 101, 102, px, py, pz, ee, mass);
67  // other
68  else {
69  (fMasterGen->event).append(particleID, 1, 0, 0, px, py, pz, ee, mass);
70  int eventSize = (fMasterGen->event).size() - 1;
71  // -log(flat) = exponential distribution
72  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
73  (fMasterGen->event)[eventSize].tau(tauTmp);
74  }
75 
76  // Here also need to add anti-particle (if any)
77  // otherwise just add a 2nd particle of the same type
78  // (for example, gamma)
79  //
80  if (fAddAntiParticle) {
81  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
82  (fMasterGen->event).append(-particleID, 23, 0, 101, -px, -py, -pz, ee, mass);
83  } else if (std::abs(particleID) == 21) { // gluons
84  (fMasterGen->event).append(21, 23, 102, 101, -px, -py, -pz, ee, mass);
85  } else {
86  if ((fMasterGen->particleData).isParticle(-particleID)) {
87  (fMasterGen->event).append(-particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
88  } else {
89  (fMasterGen->event).append(particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
90  }
91  int eventSize = (fMasterGen->event).size() - 1;
92  // -log(flat) = exponential distribution
93  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
94  (fMasterGen->event)[eventSize].tau(tauTmp);
95  }
96  }
97  }
98 
99  if (!fMasterGen->next())
100  return false;
101  evtGenDecay();
102 
103  event() = std::make_unique<HepMC::GenEvent>();
104  return toHepMC.fill_next_event(fMasterGen->event, event().get());
105  }
106 
107  const char* Py8PtGun::classname() const { return "Py8PtGun"; }
108 
110 
111 } // namespace gen
112 
113 using gen::Pythia8PtGun;
edm::GeneratorFilter< gen::Py8PtGun, gen::ExternalDecayDriver > Pythia8PtGun
Definition: Py8PtGun.cc:109
static std::vector< std::string > checklist log
tuple pp
Definition: createTree.py:17
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
boost::dynamic_bitset append(const boost::dynamic_bitset<> &bs1, const boost::dynamic_bitset<> &bs2)
this method takes two bitsets bs1 and bs2 and returns result of bs2 appended to the end of bs1 ...
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:58
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
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:107
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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
tuple size
Write out results.