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 
7 
8 namespace gen {
9 
10 class Py8PtGun : public Py8GunBase {
11 
12  public:
13 
14  Py8PtGun( edm::ParameterSet const& );
15  ~Py8PtGun() {}
16 
18  const char* classname() const;
19 
20  private:
21 
22  // PtGun particle(s) characteristics
23  double fMinEta;
24  double fMaxEta;
25  double fMinPt ;
26  double fMaxPt ;
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  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
43  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 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 eta = (fMaxEta-fMinEta) * randomEngine->flat() + fMinEta;
64  double the = 2.*atan(exp(-eta));
65 
66  double pt = (fMaxPt-fMinPt) * randomEngine->flat() + fMinPt;
67 
68  double mass = (fMasterGen->particleData).m0( particleID );
69 
70  double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
71  double ee = sqrt( pp*pp + mass*mass );
72 
73  double px = pt * cos(phi);
74  double py = pt * sin(phi);
75  double pz = pp * cos(the);
76 
77  if ( !((fMasterGen->particleData).isParticle( particleID )) )
78  {
79  particleID = std::fabs(particleID) ;
80  }
81  if( 1<= fabs(particleID) && fabs(particleID) <= 6) // quarks
82  (fMasterGen->event).append( particleID, 23, 101, 0, px, py, pz, ee, mass );
83  else if (fabs(particleID) == 21) // gluons
84  (fMasterGen->event).append( 21, 23, 101, 102, px, py, pz, ee, mass );
85  else // other
86  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
87 
88 // Here also need to add anti-particle (if any)
89 // otherwise just add a 2nd particle of the same type
90 // (for example, gamma)
91 //
92  if ( fAddAntiParticle )
93  {
94  if( 1 <= fabs(particleID) && fabs(particleID) <= 6){ // quarks
95  (fMasterGen->event).append( -particleID, 23, 0, 101, -px, -py, -pz, ee, mass );
96  }
97  else if (fabs(particleID) == 21){ // gluons
98  (fMasterGen->event).append( 21, 23, 102, 101, -px, -py, -pz, ee, mass );
99  }
100  else if ( (fMasterGen->particleData).isParticle( -particleID ) )
101  {
102  (fMasterGen->event).append( -particleID, 1, 0, 0, -px, -py, -pz, ee, mass );
103  }
104  else
105  {
106  (fMasterGen->event).append( particleID, 1, 0, 0, -px, -py, -pz, ee, mass );
107  }
108  }
109 
110  }
111 
112 
113  if ( !fMasterGen->next() ) return false;
114 
115  event().reset(new HepMC::GenEvent);
116  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
117 
118 }
119 
120 const char* Py8PtGun::classname() const
121 {
122  return "Py8PtGun";
123 }
124 
126 
127 } // end namespace
128 
129 using gen::Pythia8PtGun;
edm::GeneratorFilter< gen::Py8PtGun, gen::ExternalDecayDriver > Pythia8PtGun
Definition: Py8PtGun.cc:125
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
tuple pp
Definition: createTree.py:15
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
T eta() const
bool fAddAntiParticle
Definition: Py8PtGun.cc:27
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< int > fPartIDs
Definition: Py8GunBase.h:65
std::auto_ptr< HepMC::GenEvent > & event()
Definition: Py8GunBase.h:60
double fMinPt
Definition: Py8PtGun.cc:25
const char * classname() const
Definition: Py8PtGun.cc:120
double fMinEta
Definition: Py8PtGun.cc:23
double fMaxPhi
Definition: Py8GunBase.h:67
tuple mass
Definition: scaleCards.py:27
Py8PtGun(edm::ParameterSet const &)
Definition: Py8PtGun.cc:33
HepMC::Pythia8ToHepMC toHepMC
double fMaxPt
Definition: Py8PtGun.cc:26
double fMaxEta
Definition: Py8PtGun.cc:24
bool generatePartonsAndHadronize()
Definition: Py8PtGun.cc:48
Definition: DDAxes.h:10