CMS 3D CMS Logo

PFPileUpAlgo.cc
Go to the documentation of this file.
5 
7  pfCandidatesFromVtx_.clear();
8  pfCandidatesFromPU_.clear();
9 
10  for (unsigned i = 0; i < pfCandidates.size(); i++) {
11  const reco::PFCandidate& cand = *(pfCandidates[i]);
12 
13  int ivertex;
14 
15  switch (cand.particleId()) {
17  ivertex = chargedHadronVertex(vertices, cand);
18  break;
19  default:
20  continue;
21  }
22 
23  // no associated vertex, or primary vertex
24  // not pile-up
25  if ((ivertex == -1 || ivertex == 0 ||
26  (fNumOfPUVtxsForCharged_ > 0 && !vertices.empty() && ivertex <= int(fNumOfPUVtxsForCharged_) &&
27  std::abs(cand.vertex().z() - vertices[0].z()) < fDzCutForChargedFromPUVtxs_))) {
28  if (verbose_)
29  std::cout << "VTX " << i << " " << *(pfCandidates[i]) << std::endl;
31  } else {
32  if (verbose_)
33  std::cout << "PU " << i << " " << *(pfCandidates[i]) << std::endl;
34  // associated to a vertex
36  }
37  }
38 }
39 
41  auto const& track = pfcand.trackRef();
42  size_t iVertex = 0;
43  unsigned int index = 0;
44  unsigned int nFoundVertex = 0;
45  float bestweight = 0;
46  for (auto const& vtx : vertices) {
47  float w = vtx.trackWeight(track);
48  //select the vertex for which the track has the highest weight
49  if (w > bestweight) {
50  bestweight = w;
51  iVertex = index;
52  nFoundVertex++;
53  }
54  ++index;
55  }
56 
57  if (nFoundVertex > 0) {
58  if (nFoundVertex != 1)
59  edm::LogWarning("TrackOnTwoVertex") << "a track is shared by at least two verteces. Used to be an assert";
60  return iVertex;
61  }
62  // no vertex found with this track.
63 
64  // optional: as a secondary solution, associate the closest vertex in z
66  double dzmin = 10000;
67  double ztrack = pfcand.vertex().z();
68  bool foundVertex = false;
69  index = 0;
70  for (auto iv = vertices.begin(); iv != vertices.end(); ++iv, ++index) {
71  double dz = fabs(ztrack - iv->z());
72  if (dz < dzmin) {
73  dzmin = dz;
74  iVertex = index;
75  foundVertex = true;
76  }
77  }
78 
79  if (foundVertex)
80  return iVertex;
81  }
82 
83  return -1;
84 }
double fDzCutForChargedFromPUVtxs_
Definition: PFPileUpAlgo.h:44
int32_t *__restrict__ iv
std::vector< edm::FwdPtr< reco::PFCandidate > > PFCollection
Definition: PFPileUpAlgo.h:14
T w() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int chargedHadronVertex(const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
Definition: PFPileUpAlgo.cc:40
unsigned int fNumOfPUVtxsForCharged_
Definition: PFPileUpAlgo.h:43
PFCollection pfCandidatesFromVtx_
Definition: PFPileUpAlgo.h:49
bool verbose_
verbose ?
Definition: PFPileUpAlgo.h:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PFCollection pfCandidatesFromPU_
Definition: PFPileUpAlgo.h:50
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:42
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
void process(const PFCollection &pfCandidates, const reco::VertexCollection &vertices)
Definition: PFPileUpAlgo.cc:6
Log< level::Warning, false > LogWarning