CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Py8JetGun.cc
Go to the documentation of this file.
1 
4 
6 
7 namespace gen {
8 
9 class Py8JetGun : public Py8GunBase {
10 
11  public:
12 
13  Py8JetGun( edm::ParameterSet const& );
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 fMinP ;
25  double fMaxP ;
26  double fMinE ;
27  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  fMinP = pgun_params.getParameter<double>("MinP"); // , 0.);
43  fMaxP = pgun_params.getParameter<double>("MaxP"); // , 0.);
44  fMinE = pgun_params.getParameter<double>("MinE"); // , 0.);
45  fMaxE = pgun_params.getParameter<double>("MaxE"); // , 0.);
46 
47 }
48 
50 {
51 
52  fMasterGen->event.reset();
53 
54  double totPx = 0.;
55  double totPy = 0.;
56  double totPz = 0.;
57  double totE = 0.;
58  double totM = 0.;
59  double phi, eta, the, ee, pp;
60 
61  for ( size_t i=0; i<fPartIDs.size(); i++ )
62  {
63 
64  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
65 
66  phi = 2. * M_PI * randomEngine().flat() ;
67  the = acos( -1. + 2.*randomEngine().flat() );
68 
69  // from input
70  //
71  ee = (fMaxE-fMinE)*randomEngine().flat() + fMinE;
72 
73  double mass = (fMasterGen->particleData).m0( particleID );
74 
75  pp = sqrt( ee*ee - mass*mass );
76 
77  double px = pp * sin(the) * cos(phi);
78  double py = pp * sin(the) * sin(phi);
79  double pz = pp * cos(the);
80 
81  if ( !((fMasterGen->particleData).isParticle( particleID )) )
82  {
83  particleID = std::fabs(particleID) ;
84  }
85  (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
86 
87  // values for computing total mass
88  //
89  totPx += px;
90  totPy += py;
91  totPz += pz;
92  totE += ee;
93 
94  }
95 
96  totM = sqrt( totE*totE - (totPx*totPx+totPy*totPy+totPz*totPz) );
97 
98  //now the boost (from input params)
99  //
100  pp = (fMaxP-fMinP)*randomEngine().flat() + fMinP;
101  ee = sqrt( totM*totM + pp*pp );
102 
103  //the boost direction (from input params)
104  //
105  phi = (fMaxPhi-fMinPhi)*randomEngine().flat() + fMinPhi;
106  eta = (fMaxEta-fMinEta)*randomEngine().flat() + fMinEta;
107  the = 2.*atan(exp(-eta));
108 
109  double betaX = pp/ee * std::sin(the) * std::cos(phi);
110  double betaY = pp/ee * std::sin(the) * std::sin(phi);
111  double betaZ = pp/ee * std::cos(the);
112 
113  // boost all particles
114  //
115  (fMasterGen->event).bst( betaX, betaY, betaZ );
116 
117  if ( !fMasterGen->next() ) return false;
118 
119  event().reset(new HepMC::GenEvent);
120  return toHepMC.fill_next_event( fMasterGen->event, event().get() );
121 
122 }
123 
124 const char* Py8JetGun::classname() const
125 {
126  return "Py8JetGun";
127 }
128 
130 
131 } // end namespace
132 
133 using gen::Pythia8JetGun;
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
double fMinEta
Definition: Py8JetGun.cc:22
double fMinE
Definition: Py8JetGun.cc:26
tuple pp
Definition: createTree.py:15
const char * classname() const override
Definition: Py8JetGun.cc:124
#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 fMaxP
Definition: Py8JetGun.cc:25
T eta() const
virtual double flat() override
Definition: P8RndmEngine.cc:7
double fMaxE
Definition: Py8JetGun.cc:27
Py8JetGun(edm::ParameterSet const &)
Definition: Py8JetGun.cc:33
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
#define M_PI
P8RndmEngine & randomEngine()
edm::GeneratorFilter< gen::Py8JetGun, gen::ExternalDecayDriver > Pythia8JetGun
Definition: Py8JetGun.cc:129
double fMaxPhi
Definition: Py8GunBase.h:76
HepMC::Pythia8ToHepMC toHepMC
double fMinP
Definition: Py8JetGun.cc:24
bool generatePartonsAndHadronize() override
Definition: Py8JetGun.cc:49
double fMaxEta
Definition: Py8JetGun.cc:23
Definition: DDAxes.h:10