CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Py8EGun.cc
Go to the documentation of this file.
1 
4 
6 
7 namespace gen {
8 
9 class Py8EGun : public Py8GunBase {
10 
11  public:
12 
13  Py8EGun( edm::ParameterSet const& );
14  ~Py8EGun() {}
15 
16  bool generatePartonsAndHadronize() override;
17  const char* classname() const override;
18 
19  private:
20 
21  // EGun particle(s) characteristics
22  double fMinEta;
23  double fMaxEta;
24  double fMinE ;
25  double fMaxE ;
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  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
42  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 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 ee = (fMaxE-fMinE) * randomEngine().flat() + fMinE;
59  double eta = (fMaxEta-fMinEta) * randomEngine().flat() + fMinEta;
60  double the = 2.*atan(exp(-eta));
61 
62  double mass = (fMasterGen->particleData).m0( particleID );
63 
64  double pp = sqrt( ee*ee - mass*mass );
65  double px = pp * sin(the) * cos(phi);
66  double py = pp * sin(the) * sin(phi);
67  double pz = pp * cos(the);
68 
69  if ( !((fMasterGen->particleData).isParticle( particleID )) )
70  {
71  particleID = std::fabs(particleID) ;
72  }
73  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
74 
75 // Here also need to add anti-particle (if any)
76 // otherwise just add a 2nd particle of the same type
77 // (for example, gamma)
78 //
79  if ( fAddAntiParticle )
80  {
81  if ( (fMasterGen->particleData).isParticle( -particleID ) )
82  {
83  (fMasterGen->event).append( -particleID, 1, 0, 0, px, py, pz, ee, mass );
84  }
85  else
86  {
87  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
88  }
89  }
90 
91  }
92 
93  if ( !fMasterGen->next() ) return false;
94 
95  event().reset(new HepMC::GenEvent);
96  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
97 
98 }
99 
100 const char* Py8EGun::classname() const
101 {
102  return "Py8EGun";
103 }
104 
106 
107 } // end namespace
108 
109 using gen::Pythia8EGun;
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool fAddAntiParticle
Definition: Py8EGun.cc:26
std::auto_ptr< Pythia8::Pythia > fMasterGen
tuple pp
Definition: createTree.py:15
double fMinE
Definition: Py8EGun.cc:24
#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
double fMinEta
Definition: Py8EGun.cc:22
T eta() const
virtual double flat() override
Definition: P8RndmEngine.cc:7
T sqrt(T t)
Definition: SSEVec.h:48
double fMaxEta
Definition: Py8EGun.cc:23
bool generatePartonsAndHadronize() override
Definition: Py8EGun.cc:47
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::GeneratorFilter< gen::Py8EGun, gen::ExternalDecayDriver > Pythia8EGun
Definition: Py8EGun.cc:105
std::vector< int > fPartIDs
Definition: Py8GunBase.h:74
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:69
const char * classname() const override
Definition: Py8EGun.cc:100
P8RndmEngine & randomEngine()
double fMaxE
Definition: Py8EGun.cc:25
double fMaxPhi
Definition: Py8GunBase.h:76
HepMC::Pythia8ToHepMC toHepMC
Py8EGun(edm::ParameterSet const &)
Definition: Py8EGun.cc:32
Definition: DDAxes.h:10