CMS 3D CMS Logo

Pythia6PtGun.cc
Go to the documentation of this file.
1 
2 #include <iostream>
3 
4 #include "Pythia6PtGun.h"
5 
7 
10 
12 
14 
15 using namespace edm;
16 using namespace gen;
17 
18 Pythia6PtGun::Pythia6PtGun(const ParameterSet& pset) : Pythia6ParticleGun(pset) {
19  // ParameterSet defpset ;
20  ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters"); //, defpset ) ;
21  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
22  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
23  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 20.);
24  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 420.);
25  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
26 }
27 
29 
30 void Pythia6PtGun::generateEvent(CLHEP::HepRandomEngine*) {
31  Pythia6Service::InstanceWrapper guard(fPy6Service); // grab Py6 instance
32 
33  // now actualy, start cooking up the event gun
34  //
35 
36  // 1st, primary vertex
37  //
38  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
39 
40  // here re-create fEvt (memory)
41  //
42  fEvt = new HepMC::GenEvent();
43 
44  int ip = 1;
45  for (size_t i = 0; i < fPartIDs.size(); i++) {
46  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py6 !!!
47  int py6PID = HepPID::translatePDTtoPythia(particleID);
48  int dum = 0;
49  double pt = 0, mom = 0, ee = 0, the = 0, eta = 0;
50  double mass = pymass_(py6PID);
51 
52  // fill p(ip,5) (in PYJETS) with mass value right now,
53  // because the (hardcoded) mstu(10)=1 will make py1ent
54  // pick the mass from there
55  pyjets.p[4][ip - 1] = mass;
56 
57  double phi = (fMaxPhi - fMinPhi) * pyr_(&dum) + fMinPhi;
58 
59  eta = (fMaxEta - fMinEta) * pyr_(&dum) + fMinEta;
60 
61  the = 2. * atan(exp(-eta));
62 
63  pt = (fMaxPt - fMinPt) * pyr_(&dum) + fMinPt;
64 
65  mom = pt / sin(the);
66  ee = sqrt(mom * mom + mass * mass);
67 
68  py1ent_(ip, py6PID, ee, the, phi);
69 
70  double px = pyjets.p[0][ip - 1]; // pt*cos(phi) ;
71  double py = pyjets.p[1][ip - 1]; // pt*sin(phi) ;
72  double pz = pyjets.p[2][ip - 1]; // mom*cos(the) ;
73 
74  HepMC::FourVector p(px, py, pz, ee);
76  Part->suggest_barcode(ip);
77  Vtx->add_particle_out(Part);
78 
79  if (fAddAntiParticle) {
80  ip = ip + 1;
81  HepMC::GenParticle* APart = addAntiParticle(ip, particleID, ee, eta, phi);
82  if (APart)
83  Vtx->add_particle_out(APart);
84  }
85  ip++;
86  }
87 
88  fEvt->add_vertex(Vtx);
89 
90  // run pythia
91  pyexec_();
92 
93  return;
94 }
95 
mps_fire.i
i
Definition: mps_fire.py:428
gen::Pythia6PtGun::fMaxPt
double fMaxPt
Definition: Pythia6PtGun.h:25
EDProducer.h
gen::Pythia6PtGun::fAddAntiParticle
bool fAddAntiParticle
Definition: Pythia6PtGun.h:26
gen::pymass_
double pymass_(int &)
gen::Pythia6ParticleGun
Definition: Pythia6ParticleGun.h:28
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
gen::pyr_
double pyr_(int *idummy)
Definition: Pythia6Service.cc:60
gen::Pythia6Gun::fEvt
HepMC::GenEvent * fEvt
Definition: Pythia6Gun.h:66
gen::Pythia6PtGun::fMinPt
double fMinPt
Definition: Pythia6PtGun.h:24
gen::Pythia6PtGun::generateEvent
void generateEvent(CLHEP::HepRandomEngine *) override
Definition: Pythia6PtGun.cc:30
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
gen::pyexec_
void pyexec_()
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
gen::py1ent_
void py1ent_(int &ip, int &kf, double &pe, double &the, double &phi)
MakerMacros.h
EgammaObjectsElectrons_cfi.particleID
particleID
Definition: EgammaObjectsElectrons_cfi.py:4
gen::Pythia6Gun::addAntiParticle
HepMC::GenParticle * addAntiParticle(int &, int &, double &, double &, double &)
Definition: Pythia6Gun.cc:198
Pythia6PtGun.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PVValHelper::eta
Definition: PVValidationHelpers.h:70
gen::p
double p[5][pyjets_maxn]
Definition: Cascade2Hadronizer.cc:76
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
gen::FortranInstance::InstanceWrapper
Definition: FortranInstance.h:54
gen
Definition: PythiaDecays.h:13
gen::Pythia6PtGun::~Pythia6PtGun
~Pythia6PtGun() override
Definition: Pythia6PtGun.cc:28
gen::Pythia6PtGun
Definition: Pythia6PtGun.h:12
gen::Pythia6Gun::fMaxPhi
double fMaxPhi
Definition: Pythia6Gun.h:62
gen::Pythia6ParticleGun::fPartIDs
std::vector< int > fPartIDs
Definition: Pythia6ParticleGun.h:36
edm::ParameterSet
Definition: ParameterSet.h:47
gen::Pythia6Gun::fMinPhi
double fMinPhi
Definition: Pythia6Gun.h:61
gen::Pythia6Gun::fPy6Service
Pythia6Service * fPy6Service
Definition: Pythia6Gun.h:56
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
gen::Pythia6PtGun::fMaxEta
double fMaxEta
Definition: Pythia6PtGun.h:23
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
gen::Pythia6PtGun::fMinEta
double fMinEta
Definition: Pythia6PtGun.h:22
HepMCProduct.h
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27