Main Page
Namespaces
Namespace List
Namespace Members
All
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
o
p
q
r
s
t
u
v
w
z
Enumerator
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
:
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
Enumerator
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Properties
_
a
d
e
f
l
m
o
p
s
t
u
v
Related Functions
:
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
Package Documentation
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
GeneratorInterface
Pythia8Interface
plugins
Py8PtGun.cc
Go to the documentation of this file.
1
#include "
GeneratorInterface/Core/interface/GeneratorFilter.h
"
2
#include "
GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h
"
3
4
#include "
GeneratorInterface/Pythia8Interface/interface/Py8GunBase.h
"
5
6
namespace
gen
{
7
8
class
Py8PtGun
:
public
Py8GunBase
{
9
public
:
10
Py8PtGun
(
edm::ParameterSet
const
&);
11
~Py8PtGun
()
override
{}
12
13
bool
generatePartonsAndHadronize
()
override
;
14
const
char
*
classname
()
const override
;
15
16
private
:
17
// PtGun particle(s) characteristics
18
double
fMinEta
;
19
double
fMaxEta
;
20
double
fMinPt
;
21
double
fMaxPt
;
22
bool
fAddAntiParticle
;
23
};
24
25
// implementation
26
//
27
Py8PtGun::Py8PtGun
(
edm::ParameterSet
const
& ps) :
Py8GunBase
(ps) {
28
// ParameterSet defpset ;
29
edm::ParameterSet
pgun_params = ps.
getParameter
<
edm::ParameterSet
>(
"PGunParameters"
);
// , defpset ) ;
30
fMinEta
= pgun_params.
getParameter
<
double
>(
"MinEta"
);
// ,-2.2);
31
fMaxEta
= pgun_params.
getParameter
<
double
>(
"MaxEta"
);
// , 2.2);
32
fMinPt
= pgun_params.
getParameter
<
double
>(
"MinPt"
);
// , 0.);
33
fMaxPt
= pgun_params.
getParameter
<
double
>(
"MaxPt"
);
// , 0.);
34
fAddAntiParticle
= pgun_params.
getParameter
<
bool
>(
"AddAntiParticle"
);
//, false) ;
35
}
36
37
bool
Py8PtGun::generatePartonsAndHadronize
() {
38
fMasterGen
->event.reset();
39
40
for
(
size_t
i
= 0;
i
<
fPartIDs
.size();
i
++) {
41
int
particleID
=
fPartIDs
[
i
];
// this is PDG - need to convert to Py8 ???
42
43
double
phi = (
fMaxPhi
-
fMinPhi
) *
randomEngine
().
flat
() +
fMinPhi
;
44
double
eta
= (
fMaxEta
-
fMinEta
) *
randomEngine
().
flat
() +
fMinEta
;
45
double
the = 2. * atan(
exp
(-
eta
));
46
47
double
pt
= (
fMaxPt
-
fMinPt
) *
randomEngine
().
flat
() +
fMinPt
;
48
49
double
mass
= (
fMasterGen
->particleData).m0(
particleID
);
50
51
double
pp
=
pt
/
sin
(the);
// sqrt( ee*ee - mass*mass );
52
double
ee =
sqrt
(
pp
*
pp
+
mass
*
mass
);
53
54
double
px
=
pt
*
cos
(phi);
55
double
py
=
pt
*
sin
(phi);
56
double
pz =
pp
*
cos
(the);
57
58
if
(!((
fMasterGen
->particleData).isParticle(
particleID
))) {
59
particleID
=
std::abs
(
particleID
);
60
}
61
if
(1 <=
std::abs
(
particleID
) &&
std::abs
(
particleID
) <= 6)
// quarks
62
(
fMasterGen
->event).
append
(
particleID
, 23, 101, 0,
px
,
py
, pz, ee,
mass
);
63
else
if
(
std::abs
(
particleID
) == 21)
// gluons
64
(
fMasterGen
->event).append(21, 23, 101, 102,
px
,
py
, pz, ee,
mass
);
65
// other
66
else
{
67
(
fMasterGen
->event).
append
(
particleID
, 1, 0, 0,
px
,
py
, pz, ee,
mass
);
68
int
eventSize = (
fMasterGen
->event).
size
() - 1;
69
// -log(flat) = exponential distribution
70
double
tauTmp = -(
fMasterGen
->event)[eventSize].tau0() *
log
(
randomEngine
().flat());
71
(
fMasterGen
->event)[eventSize].
tau
(tauTmp);
72
}
73
74
// Here also need to add anti-particle (if any)
75
// otherwise just add a 2nd particle of the same type
76
// (for example, gamma)
77
//
78
if
(
fAddAntiParticle
) {
79
if
(1 <=
std::abs
(
particleID
) &&
std::abs
(
particleID
) <= 6) {
// quarks
80
(
fMasterGen
->event).
append
(-
particleID
, 23, 0, 101, -
px
, -
py
, -pz, ee,
mass
);
81
}
else
if
(
std::abs
(
particleID
) == 21) {
// gluons
82
(
fMasterGen
->event).
append
(21, 23, 102, 101, -
px
, -
py
, -pz, ee,
mass
);
83
}
else
{
84
if
((
fMasterGen
->particleData).isParticle(-
particleID
)) {
85
(
fMasterGen
->event).
append
(-
particleID
, 1, 0, 0, -
px
, -
py
, -pz, ee,
mass
);
86
}
else
{
87
(
fMasterGen
->event).
append
(
particleID
, 1, 0, 0, -
px
, -
py
, -pz, ee,
mass
);
88
}
89
int
eventSize = (
fMasterGen
->event).
size
() - 1;
90
// -log(flat) = exponential distribution
91
double
tauTmp = -(
fMasterGen
->event)[eventSize].tau0() *
log
(
randomEngine
().flat());
92
(
fMasterGen
->event)[eventSize].
tau
(tauTmp);
93
}
94
}
95
}
96
97
if
(!
fMasterGen
->next())
98
return
false
;
99
evtGenDecay
();
100
101
event
().reset(
new
HepMC::GenEvent
);
102
return
toHepMC
.fill_next_event(
fMasterGen
->event,
event
().
get
());
103
}
104
105
const
char
*
Py8PtGun::classname
()
const
{
return
"Py8PtGun"
; }
106
107
typedef
edm::GeneratorFilter<gen::Py8PtGun, gen::ExternalDecayDriver>
Pythia8PtGun
;
108
109
}
// namespace gen
110
111
using
gen::Pythia8PtGun
;
112
DEFINE_FWK_MODULE
(
Pythia8PtGun
);
gen::Py8PtGun::fMinEta
double fMinEta
Definition:
Py8PtGun.cc:18
gen::Py8InterfaceBase::fMasterGen
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition:
Py8InterfaceBase.h:44
mps_fire.i
i
Definition:
mps_fire.py:355
gen::Py8GunBase::fMinPhi
double fMinPhi
Definition:
Py8GunBase.h:58
gen::Py8PtGun::~Py8PtGun
~Py8PtGun() override
Definition:
Py8PtGun.cc:11
metsig::tau
Definition:
SignAlgoResolutions.h:49
DiDispStaMuonMonitor_cfi.pt
pt
Definition:
DiDispStaMuonMonitor_cfi.py:39
gen::P8RndmEngine::flat
double flat() override
Definition:
P8RndmEngine.cc:7
multPhiCorr_741_25nsDY_cfi.py
py
Definition:
multPhiCorr_741_25nsDY_cfi.py:12
ExternalDecayDriver.h
gen::Py8PtGun::fMaxEta
double fMaxEta
Definition:
Py8PtGun.cc:19
gen::Py8PtGun::fAddAntiParticle
bool fAddAntiParticle
Definition:
Py8PtGun.cc:22
gen::Py8PtGun::fMaxPt
double fMaxPt
Definition:
Py8PtGun.cc:21
HepMC::GenEvent
Definition:
hepmc_rootio.cc:9
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
EgammaObjectsElectrons_cfi.particleID
particleID
Definition:
EgammaObjectsElectrons_cfi.py:4
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
edm::GeneratorFilter
Definition:
GeneratorFilter.h:40
gen::Py8GunBase::evtGenDecay
void evtGenDecay()
Definition:
Py8GunBase.cc:144
gen::Py8PtGun::classname
const char * classname() const override
Definition:
Py8PtGun.cc:105
PVValHelper::eta
Definition:
PVValidationHelpers.h:69
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
gen::Py8PtGun::fMinPt
double fMinPt
Definition:
Py8PtGun.cc:20
gen
Definition:
PythiaDecays.h:13
gen::Py8InterfaceBase::randomEngine
P8RndmEngine & randomEngine()
Definition:
Py8InterfaceBase.h:41
edm::ParameterSet
Definition:
ParameterSet.h:36
GeneratorFilter.h
mps_setup.append
append
Definition:
mps_setup.py:85
gen::Py8GunBase::fPartIDs
std::vector< int > fPartIDs
Definition:
Py8GunBase.h:57
get
#define get
gen::Py8PtGun::generatePartonsAndHadronize
bool generatePartonsAndHadronize() override
Definition:
Py8PtGun.cc:37
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
multPhiCorr_741_25nsDY_cfi.px
px
Definition:
multPhiCorr_741_25nsDY_cfi.py:10
gen::BaseHadronizer::event
std::unique_ptr< HepMC::GenEvent > & event()
Definition:
BaseHadronizer.h:86
gen::Py8PtGun
Definition:
Py8PtGun.cc:8
gen::Py8PtGun::Py8PtGun
Py8PtGun(edm::ParameterSet const &)
Definition:
Py8PtGun.cc:27
gen::Py8GunBase
Definition:
Py8GunBase.h:40
gen::Py8InterfaceBase::toHepMC
HepMC::Pythia8ToHepMC toHepMC
Definition:
Py8InterfaceBase.h:46
EgHLTOffHistBins_cfi.mass
mass
Definition:
EgHLTOffHistBins_cfi.py:34
gen::Pythia8PtGun
edm::GeneratorFilter< gen::Py8PtGun, gen::ExternalDecayDriver > Pythia8PtGun
Definition:
Py8PtGun.cc:107
dqm-mbProfile.log
log
Definition:
dqm-mbProfile.py:17
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
createTree.pp
pp
Definition:
createTree.py:17
JetChargeProducer_cfi.exp
exp
Definition:
JetChargeProducer_cfi.py:6
gen::Py8GunBase::fMaxPhi
double fMaxPhi
Definition:
Py8GunBase.h:59
Py8GunBase.h
event
Definition:
event.py:1
findQualityFiles.size
size
Write out results.
Definition:
findQualityFiles.py:443
Generated for CMSSW Reference Manual by
1.8.16