GeneratorInterface
Pythia6Interface
plugins
Pythia6EGun.cc
Go to the documentation of this file.
1
2
#include <iostream>
3
4
#include "
Pythia6EGun.h
"
5
6
#include "
FWCore/Utilities/interface/Exception.h
"
7
8
#include "
FWCore/Framework/interface/EDProducer.h
"
9
#include "
FWCore/Framework/interface/EventSetup.h
"
10
11
#include "
SimDataFormats/GeneratorProducts/interface/HepMCProduct.h
"
12
13
#include "
FWCore/Framework/interface/MakerMacros.h
"
14
15
using namespace
edm
;
16
using namespace
gen
;
17
18
Pythia6EGun::Pythia6EGun(
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
fMinE
= pgun_params.
getParameter
<
double
>(
"MinE"
);
// , 0.);
24
fMaxE
= pgun_params.
getParameter
<
double
>(
"MaxE"
);
// , 0.);
25
fAddAntiParticle
= pgun_params.
getParameter
<
bool
>(
"AddAntiParticle"
);
//, false) ;
26
}
27
28
Pythia6EGun::~Pythia6EGun
() {}
29
30
void
Pythia6EGun::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
ee = 0, the = 0,
eta
= 0;
50
double
mass
=
pymass_
(
particleID
);
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
ee = (
fMaxE
-
fMinE
) *
pyr_
(&dum) +
fMinE
;
59
eta
= (
fMaxEta
-
fMinEta
) *
pyr_
(&dum) +
fMinEta
;
60
the = 2. * atan(
exp
(-
eta
));
61
62
py1ent_
(ip, py6PID, ee, the, phi);
63
64
double
px
= pyjets.p[0][ip - 1];
// pt*cos(phi) ;
65
double
py
= pyjets.p[1][ip - 1];
// pt*sin(phi) ;
66
double
pz = pyjets.p[2][ip - 1];
// mom*cos(the) ;
67
68
HepMC::FourVector
p
(
px
,
py
, pz, ee);
69
HepMC::GenParticle
* Part =
new
HepMC::GenParticle
(
p
,
particleID
, 1);
70
Part->suggest_barcode(ip);
71
Vtx->add_particle_out(Part);
72
73
if
(
fAddAntiParticle
) {
74
ip = ip + 1;
75
HepMC::GenParticle
* APart =
addAntiParticle
(ip,
particleID
, ee,
eta
, phi);
76
if
(APart)
77
Vtx->add_particle_out(APart);
78
}
79
ip++;
80
}
81
82
fEvt
->add_vertex(Vtx);
83
84
// run pythia
85
pyexec_
();
86
87
return
;
88
}
89
90
DEFINE_FWK_MODULE
(
Pythia6EGun
);
gen::Pythia6EGun::fMaxEta
double fMaxEta
Definition:
Pythia6EGun.h:23
mps_fire.i
i
Definition:
mps_fire.py:428
EDProducer.h
gen::pymass_
double pymass_(int &)
gen::Pythia6EGun::generateEvent
void generateEvent(CLHEP::HepRandomEngine *) override
Definition:
Pythia6EGun.cc:30
gen::Pythia6ParticleGun
Definition:
Pythia6ParticleGun.h:28
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:59
gen::Pythia6Gun::fEvt
HepMC::GenEvent * fEvt
Definition:
Pythia6Gun.h:66
HepMC::GenEvent
Definition:
hepmc_rootio.cc:9
gen::pyexec_
void pyexec_()
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
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
gen::Pythia6EGun::fAddAntiParticle
bool fAddAntiParticle
Definition:
Pythia6EGun.h:26
gen::FortranInstance::InstanceWrapper
Definition:
FortranInstance.h:54
gen::Pythia6EGun::fMinE
double fMinE
Definition:
Pythia6EGun.h:24
gen
Definition:
PythiaDecays.h:13
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::Pythia6EGun::~Pythia6EGun
~Pythia6EGun() override
Definition:
Pythia6EGun.cc:28
gen::Pythia6Gun::fMinPhi
double fMinPhi
Definition:
Pythia6Gun.h:61
gen::Pythia6Gun::fPy6Service
Pythia6Service * fPy6Service
Definition:
Pythia6Gun.h:56
gen::Pythia6EGun::fMinEta
double fMinEta
Definition:
Pythia6EGun.h:22
multPhiCorr_741_25nsDY_cfi.px
px
Definition:
multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition:
GenParticle.py:18
Pythia6EGun.h
EgHLTOffHistBins_cfi.mass
mass
Definition:
EgHLTOffHistBins_cfi.py:34
gen::Pythia6EGun::fMaxE
double fMaxE
Definition:
Pythia6EGun.h:25
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
Exception.h
HepMCProduct.h
JetChargeProducer_cfi.exp
exp
Definition:
JetChargeProducer_cfi.py:6
gen::Pythia6EGun
Definition:
Pythia6EGun.h:12
muonDTDigis_cfi.pset
pset
Definition:
muonDTDigis_cfi.py:27
Generated for CMSSW Reference Manual by
1.8.16