CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Py8PtGun.cc
Go to the documentation of this file.
1 
4 
6 
7 namespace gen {
8 
9 class Py8PtGun : public Py8GunBase {
10 
11  public:
12 
13  Py8PtGun( edm::ParameterSet const& );
14  ~Py8PtGun() {}
15 
16  bool generatePartonsAndHadronize() override;
17  const char* classname() const override;
18 
19  private:
20 
21  // PtGun particle(s) characteristics
22  double fMinEta;
23  double fMaxEta;
24  double fMinPt ;
25  double fMaxPt ;
27 
28 };
29 
30 // implementation
31 //
33  : Py8GunBase(ps)
34 {
35 
36  // ParameterSet defpset ;
37  edm::ParameterSet pgun_params =
38  ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
39  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
40  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
41  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
42  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 0.);
43  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
44 
45 }
46 
48 {
49 
50  fMasterGen->event.reset();
51 
52  for ( size_t i=0; i<fPartIDs.size(); i++ )
53  {
54 
55  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
56 
57  double phi = (fMaxPhi-fMinPhi) * randomEngine().flat() + fMinPhi;
58  double eta = (fMaxEta-fMinEta) * randomEngine().flat() + fMinEta;
59  double the = 2.*atan(exp(-eta));
60 
61  double pt = (fMaxPt-fMinPt) * randomEngine().flat() + fMinPt;
62 
63  double mass = (fMasterGen->particleData).m0( particleID );
64 
65  double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
66  double ee = sqrt( pp*pp + mass*mass );
67 
68  double px = pt * cos(phi);
69  double py = pt * sin(phi);
70  double pz = pp * cos(the);
71 
72  if ( !((fMasterGen->particleData).isParticle( particleID )) )
73  {
74  particleID = std::fabs(particleID) ;
75  }
76  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
77 
78 // Here also need to add anti-particle (if any)
79 // otherwise just add a 2nd particle of the same type
80 // (for example, gamma)
81 //
82  if ( fAddAntiParticle )
83  {
84  if ( (fMasterGen->particleData).isParticle( -particleID ) )
85  {
86  (fMasterGen->event).append( -particleID, 1, 0, 0, px, py, pz, ee, mass );
87  }
88  else
89  {
90  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
91  }
92  }
93 
94  }
95 
96  if ( !fMasterGen->next() ) return false;
97 
98  event().reset(new HepMC::GenEvent);
99  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
100 
101 }
102 
103 const char* Py8PtGun::classname() const
104 {
105  return "Py8PtGun";
106 }
107 
109 
110 } // end namespace
111 
112 using gen::Pythia8PtGun;
edm::GeneratorFilter< gen::Py8PtGun, gen::ExternalDecayDriver > Pythia8PtGun
Definition: Py8PtGun.cc:108
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
tuple pp
Definition: createTree.py:15
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double fMinPhi
Definition: Py8GunBase.h:75
T eta() const
virtual double flat() override
Definition: P8RndmEngine.cc:7
bool fAddAntiParticle
Definition: Py8PtGun.cc:26
const char * classname() const override
Definition: Py8PtGun.cc:103
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< int > fPartIDs
Definition: Py8GunBase.h:74
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:69
P8RndmEngine & randomEngine()
double fMinPt
Definition: Py8PtGun.cc:24
double fMinEta
Definition: Py8PtGun.cc:22
bool generatePartonsAndHadronize() override
Definition: Py8PtGun.cc:47
double fMaxPhi
Definition: Py8GunBase.h:76
Py8PtGun(edm::ParameterSet const &)
Definition: Py8PtGun.cc:32
HepMC::Pythia8ToHepMC toHepMC
double fMaxPt
Definition: Py8PtGun.cc:25
double fMaxEta
Definition: Py8PtGun.cc:23
Definition: DDAxes.h:10