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 
7 
8 namespace gen {
9 
10 class Py8EGun : public Py8GunBase {
11 
12  public:
13 
14  Py8EGun( edm::ParameterSet const& );
15  ~Py8EGun() {}
16 
18  const char* classname() const;
19 
20  private:
21 
22  // EGun particle(s) characteristics
23  double fMinEta;
24  double fMaxEta;
25  double fMinE ;
26  double fMaxE ;
28 
29 };
30 
31 // implementation
32 //
34  : Py8GunBase(ps)
35 {
36 
37  // ParameterSet defpset ;
38  edm::ParameterSet pgun_params =
39  ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
40  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
41  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
42  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
43  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 0.);
44  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
45 
46 }
47 
49 {
50 
51  fMasterGen->event.reset();
52 
53  for ( size_t i=0; i<fPartIDs.size(); i++ )
54  {
55 
56  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
57 
58  // FIXME !!!
59  // Ouch, it's using bare randomEngine pointer - that's NOT safe.
60  // Need to hold a pointer somewhere properly !!!
61  //
62  double phi = (fMaxPhi-fMinPhi) * randomEngine->flat() + fMinPhi;
63  double ee = (fMaxE-fMinE) * randomEngine->flat() + fMinE;
64  double eta = (fMaxEta-fMinEta) * randomEngine->flat() + fMinEta;
65  double the = 2.*atan(exp(-eta));
66 
67  double mass = (fMasterGen->particleData).m0( particleID );
68 
69  double pp = sqrt( ee*ee - mass*mass );
70  double px = pp * sin(the) * cos(phi);
71  double py = pp * sin(the) * sin(phi);
72  double pz = pp * cos(the);
73 
74  if ( !((fMasterGen->particleData).isParticle( particleID )) )
75  {
76  particleID = std::fabs(particleID) ;
77  }
78  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
79 
80 // Here also need to add anti-particle (if any)
81 // otherwise just add a 2nd particle of the same type
82 // (for example, gamma)
83 //
84  if ( fAddAntiParticle )
85  {
86  if ( (fMasterGen->particleData).isParticle( -particleID ) )
87  {
88  (fMasterGen->event).append( -particleID, 1, 0, 0, px, py, pz, ee, mass );
89  }
90  else
91  {
92  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
93  }
94  }
95 
96  }
97 
98  if ( !fMasterGen->next() ) return false;
99 
100  event().reset(new HepMC::GenEvent);
101  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
102 
103 }
104 
105 const char* Py8EGun::classname() const
106 {
107  return "Py8EGun";
108 }
109 
111 
112 } // end namespace
113 
114 using gen::Pythia8EGun;
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool fAddAntiParticle
Definition: Py8EGun.cc:27
std::auto_ptr< Pythia8::Pythia > fMasterGen
tuple pp
Definition: createTree.py:15
double fMinE
Definition: Py8EGun.cc:25
CLHEP::HepRandomEngine * randomEngine
Definition: PYR.cc:4
#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:66
const char * classname() const
Definition: Py8EGun.cc:105
double fMinEta
Definition: Py8EGun.cc:23
T eta() const
T sqrt(T t)
Definition: SSEVec.h:46
double fMaxEta
Definition: Py8EGun.cc:24
bool generatePartonsAndHadronize()
Definition: Py8EGun.cc:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::GeneratorFilter< gen::Py8EGun, gen::ExternalDecayDriver > Pythia8EGun
Definition: Py8EGun.cc:110
std::vector< int > fPartIDs
Definition: Py8GunBase.h:65
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:60
double fMaxE
Definition: Py8EGun.cc:26
double fMaxPhi
Definition: Py8GunBase.h:67
tuple mass
Definition: scaleCards.py:27
HepMC::Pythia8ToHepMC toHepMC
Py8EGun(edm::ParameterSet const &)
Definition: Py8EGun.cc:33
Definition: DDAxes.h:10