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
CommonTools
ParticleFlow
src
PFPileUpAlgo.cc
Go to the documentation of this file.
1
#include "
CommonTools/ParticleFlow/interface/PFPileUpAlgo.h
"
2
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidate.h
"
3
#include "
DataFormats/VertexReco/interface/Vertex.h
"
4
5
void
PFPileUpAlgo::process
(
const
PFCollection
&
pfCandidates
,
const
reco::VertexCollection
&
vertices
) {
6
pfCandidatesFromVtx_
.clear();
7
pfCandidatesFromPU_
.clear();
8
9
for
(
unsigned
i
= 0;
i
<
pfCandidates
.size();
i
++) {
10
const
reco::PFCandidate
&
cand
= *(
pfCandidates
[
i
]);
11
12
int
ivertex;
13
14
switch
(
cand
.particleId()) {
15
case
reco::PFCandidate::h
:
16
ivertex =
chargedHadronVertex
(
vertices
,
cand
);
17
break
;
18
default
:
19
continue
;
20
}
21
22
// no associated vertex, or primary vertex
23
// not pile-up
24
if
(ivertex == -1 || ivertex == 0) {
25
if
(
verbose_
)
26
std::cout
<<
"VTX "
<<
i
<<
" "
<< *(
pfCandidates
[
i
]) << std::endl;
27
pfCandidatesFromVtx_
.push_back(
pfCandidates
[
i
]);
28
}
else
{
29
if
(
verbose_
)
30
std::cout
<<
"PU "
<<
i
<<
" "
<< *(
pfCandidates
[
i
]) << std::endl;
31
// associated to a vertex
32
pfCandidatesFromPU_
.push_back(
pfCandidates
[
i
]);
33
}
34
}
35
}
36
37
int
PFPileUpAlgo::chargedHadronVertex
(
const
reco::VertexCollection
&
vertices
,
const
reco::PFCandidate
&
pfcand
)
const
{
38
auto
const
&
track
=
pfcand
.trackRef();
39
size_t
iVertex = 0;
40
unsigned
int
index
= 0;
41
unsigned
int
nFoundVertex = 0;
42
float
bestweight = 0;
43
for
(
auto
const
&
vtx
:
vertices
) {
44
float
w
=
vtx
.trackWeight(
track
);
45
//select the vertex for which the track has the highest weight
46
if
(
w
> bestweight) {
47
bestweight =
w
;
48
iVertex =
index
;
49
nFoundVertex++;
50
}
51
++
index
;
52
}
53
54
if
(nFoundVertex > 0) {
55
if
(nFoundVertex != 1)
56
edm::LogWarning
(
"TrackOnTwoVertex"
) <<
"a track is shared by at least two verteces. Used to be an assert"
;
57
return
iVertex;
58
}
59
// no vertex found with this track.
60
61
// optional: as a secondary solution, associate the closest vertex in z
62
if
(
checkClosestZVertex_
) {
63
double
dzmin = 10000;
64
double
ztrack =
pfcand
.vertex().z();
65
bool
foundVertex =
false
;
66
index
= 0;
67
for
(
auto
iv =
vertices
.begin(); iv !=
vertices
.end(); ++iv, ++
index
) {
68
double
dz
= fabs(ztrack - iv->z());
69
if
(
dz
< dzmin) {
70
dzmin =
dz
;
71
iVertex =
index
;
72
foundVertex =
true
;
73
}
74
}
75
76
if
(foundVertex)
77
return
iVertex;
78
}
79
80
return
-1;
81
}
zmumugammaAnalyzer_cfi.pfCandidates
pfCandidates
Definition:
zmumugammaAnalyzer_cfi.py:11
mps_fire.i
i
Definition:
mps_fire.py:355
PFCandidate.h
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition:
VertexFwd.h:9
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
reco::PFCandidate::h
Definition:
PFCandidate.h:45
PFPileUpAlgo::process
void process(const PFCollection &pfCandidates, const reco::VertexCollection &vertices)
Definition:
PFPileUpAlgo.cc:5
pfDeepBoostedJetPreprocessParams_cfi.pfcand
pfcand
Definition:
pfDeepBoostedJetPreprocessParams_cfi.py:8
w
const double w
Definition:
UKUtility.cc:23
badGlobalMuonTaggersAOD_cff.vtx
vtx
Definition:
badGlobalMuonTaggersAOD_cff.py:5
Vertex.h
edm::LogWarning
Definition:
MessageLogger.h:141
cand
Definition:
decayParser.h:34
PFPileUpAlgo.h
PVValHelper::dz
Definition:
PVValidationHelpers.h:50
PFPileUpAlgo::pfCandidatesFromPU_
PFCollection pfCandidatesFromPU_
Definition:
PFPileUpAlgo.h:46
PFPileUpAlgo::pfCandidatesFromVtx_
PFCollection pfCandidatesFromVtx_
Definition:
PFPileUpAlgo.h:45
AlignmentPI::index
index
Definition:
AlignmentPayloadInspectorHelper.h:46
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition:
PFCandidate.h:40
HLT_2018_cff.track
track
Definition:
HLT_2018_cff.py:10352
PFPileUpAlgo::chargedHadronVertex
int chargedHadronVertex(const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
Definition:
PFPileUpAlgo.cc:37
PFPileUpAlgo::verbose_
bool verbose_
verbose ?
Definition:
PFPileUpAlgo.h:43
PFPileUpAlgo::checkClosestZVertex_
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition:
PFPileUpAlgo.h:40
PFPileUpAlgo::PFCollection
std::vector< edm::FwdPtr< reco::PFCandidate > > PFCollection
Definition:
PFPileUpAlgo.h:14
pwdgSkimBPark_cfi.vertices
vertices
Definition:
pwdgSkimBPark_cfi.py:7
Generated for CMSSW Reference Manual by
1.8.16