CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Pythia6JetGun.cc
Go to the documentation of this file.
1 
2 #include <iostream>
3 
4 #include "Pythia6JetGun.h"
5 
7 
10 
12 
14 
15 using namespace edm;
16 using namespace gen;
17 
18 Pythia6JetGun::Pythia6JetGun( const ParameterSet& pset ) :
19  Pythia6ParticleGun(pset),
20  fMinEta(0.), fMaxEta(0.),
21  fMinE(0.), fMaxE(0.),
22  fMinP(0.), fMaxP(0.)
23 {
24 
25  ParameterSet pgun_params =
26  pset.getParameter<ParameterSet>("PGunParameters");
27  fMinEta = pgun_params.getParameter<double>("MinEta");
28  fMaxEta = pgun_params.getParameter<double>("MaxEta");
29  fMinE = pgun_params.getParameter<double>("MinE");
30  fMaxE = pgun_params.getParameter<double>("MaxE");
31  fMinP = pgun_params.getParameter<double>("MinP");
32  fMaxP = pgun_params.getParameter<double>("MaxP");
33 
34 }
35 
37 {
38 }
39 
40 void Pythia6JetGun::generateEvent(CLHEP::HepRandomEngine*)
41 {
42  Pythia6Service::InstanceWrapper guard(fPy6Service); // grab Py6 instance
43 
44  // now actualy, start cooking up the event gun
45  //
46 
47  // 1st, primary vertex
48  //
49  HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
50 
51  // here re-create fEvt (memory)
52  //
53  fEvt = new HepMC::GenEvent() ;
54 
55  int ip=1;
56  double totPx = 0.;
57  double totPy = 0.;
58  double totPz = 0.;
59  double totE = 0.;
60  double totM = 0.;
61  double phi, eta, the, ee, pp;
62  int dum = 0;
63  for ( size_t i=0; i<fPartIDs.size(); i++ )
64  {
65  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py6 !!!
66  int py6PID = HepPID::translatePDTtoPythia( particleID );
67 
68  // internal numbers
69  //
70  phi = 2. * M_PI * pyr_(&dum);
71  the = std::acos( -1. + 2.*pyr_(&dum) );
72 
73  // from input
74  //
75  ee = (fMaxE-fMinE)*pyr_(&dum)+fMinE;
76 
77  // fill p(ip,5) (in PYJETS) with mass value right now,
78  // because the (hardcoded) mstu(10)=1 will make py1ent
79  // pick the mass from there
80  double mass = pymass_(py6PID);
81  pyjets.p[4][ip-1]=mass;
82 
83  // add entry to py6
84  //
85  py1ent_(ip, py6PID, ee, the, phi);
86 
87  // values for computing total mass
88  //
89  totPx += pyjets.p[0][ip-1];
90  totPy += pyjets.p[1][ip-1];
91  totPz += pyjets.p[2][ip-1];
92  totE += pyjets.p[3][ip-1];
93 
94  ip++;
95 
96  } // end forming up py6 record of the jet
97 
98 
99  // compute total mass
100  //
101  totM = std::sqrt( totE*totE - (totPx*totPx+totPy*totPy+totPz*totPz) );
102 
103  //now the boost (from input params)
104  //
105  pp = (fMaxP-fMinP)*pyr_(&dum)+fMinP;
106  ee = std::sqrt( totM*totM + pp*pp );
107 
108  //the boost direction (from input params)
109  //
110  phi = (fMaxPhi-fMinPhi)*pyr_(&dum)+fMinPhi;
111  eta = (fMaxEta-fMinEta)*pyr_(&dum)+fMinEta;
112  the = 2.*atan(exp(-eta));
113 
114  double betaX = pp/ee * std::sin(the) * std::cos(phi);
115  double betaY = pp/ee * std::sin(the) * std::sin(phi);
116  double betaZ = pp/ee * std::cos(the);
117 
118  // boost all particles
119  // the first 2 params (-1) tell to boost all (fisrt-to-last),
120  // and the next 2 params (0.) tell no rotation
121  //
122  int first=-1, last=-1;
123  double rothe=0, rophi=0.;
124 
125  pyrobo_( first, last, rothe, rophi, betaX, betaY, betaZ );
126 
127  // event should be formed from boosted record !!!
128  // that's why additional loop
129  //
130  for ( int i=0; i<pyjets.n; i++ )
131  {
132  HepMC::FourVector p(pyjets.p[0][i],pyjets.p[1][i],pyjets.p[2][i],pyjets.p[3][i]) ;
133  HepMC::GenParticle* Part =
134  new HepMC::GenParticle(p,
135  HepPID::translatePythiatoPDT( pyjets.k[1][i] ),
136  1);
137  Part->suggest_barcode( i+1 ) ;
138  Vtx->add_particle_out(Part);
139  }
140  fEvt->add_vertex(Vtx);
141 
142  // run pythia
143  pyexec_();
144 
145  return;
146 }
147 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
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
HepMC::GenEvent * fEvt
Definition: Pythia6Gun.h:71
def gen
run2 Cosmic #### Run 256259 @ 0T 2015C### Run 272133 @ 3.8T 2016B###
void generateEvent(CLHEP::HepRandomEngine *)
double p[5][pyjets_maxn]
double fMinPhi
Definition: Pythia6Gun.h:66
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double pymass_(int &)
std::vector< int > fPartIDs
#define M_PI
double fMaxPhi
Definition: Pythia6Gun.h:67
Geom::Phi< T > phi() const
virtual ~Pythia6JetGun()
void py1ent_(int &ip, int &kf, double &pe, double &the, double &phi)
Pythia6Service * fPy6Service
Definition: Pythia6Gun.h:61
double pyr_(int *idummy)
void pyrobo_(int &, int &, double &, double &, double &, double &, double &)
void pyexec_()