CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFPileUpAlgo.cc
Go to the documentation of this file.
4 
6  const reco::VertexCollection & vertices) {
7 
8  pfCandidatesFromVtx_.clear();
9  pfCandidatesFromPU_.clear();
10 
11  for( unsigned i=0; i<pfCandidates.size(); i++ ) {
12 
13  const reco::PFCandidate& cand = * ( pfCandidates[i] );
14 
15  int ivertex;
16 
17  switch( cand.particleId() ) {
19  ivertex = chargedHadronVertex( vertices, cand );
20  break;
21  default:
22  continue;
23  }
24 
25  // no associated vertex, or primary vertex
26  // not pile-up
27  if( ivertex == -1 ||
28  ivertex == 0 ) {
29  if(verbose_)
30  std::cout<<"VTX "<<i<<" "<< *(pfCandidates[i])<<std::endl;
31  pfCandidatesFromVtx_.push_back( pfCandidates[i] );
32  } else {
33  if(verbose_)
34  std::cout<<"PU "<<i<<" "<< *(pfCandidates[i])<<std::endl;
35  // associated to a vertex
36  pfCandidatesFromPU_.push_back( pfCandidates[i] );
37  }
38  }
39 }
40 
41 
42 int
44 
45 
46  reco::TrackBaseRef trackBaseRef( pfcand.trackRef() );
47 
48  size_t iVertex = 0;
49  unsigned index=0;
50  unsigned nFoundVertex = 0;
51  typedef reco::VertexCollection::const_iterator IV;
53  float bestweight=0;
54  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
55 
56  const reco::Vertex& vtx = *iv;
57 
58  // loop on tracks in vertices
59  for(IT iTrack=vtx.tracks_begin();
60  iTrack!=vtx.tracks_end(); ++iTrack) {
61 
62  const reco::TrackBaseRef& baseRef = *iTrack;
63 
64  // one of the tracks in the vertex is the same as
65  // the track considered in the function
66  if(baseRef == trackBaseRef ) {
67  float w = vtx.trackWeight(baseRef);
68  //select the vertex for which the track has the highest weight
69  if (w > bestweight){
70  bestweight=w;
71  iVertex=index;
72  nFoundVertex++;
73  }
74  }
75  }
76  }
77 
78  if (nFoundVertex>0){
79  if (nFoundVertex!=1)
80  edm::LogWarning("TrackOnTwoVertex")<<"a track is shared by at least two verteces. Used to be an assert";
81  return iVertex;
82  }
83  // no vertex found with this track.
84 
85  // optional: as a secondary solution, associate the closest vertex in z
86  if ( checkClosestZVertex_ ) {
87 
88  double dzmin = 10000;
89  double ztrack = pfcand.vertex().z();
90  bool foundVertex = false;
91  index = 0;
92  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
93 
94  double dz = fabs(ztrack - iv->z());
95  if(dz<dzmin) {
96  dzmin = dz;
97  iVertex = index;
98  foundVertex = true;
99  }
100  }
101 
102  if( foundVertex )
103  return iVertex;
104 
105  }
106 
107 
108  return -1 ;
109 }
110 
int i
Definition: DBlmapReader.cc:9
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:349
PFCollection pfCandidatesFromVtx_
Definition: PFPileUpAlgo.h:50
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
bool verbose_
verbose ?
Definition: PFPileUpAlgo.h:48
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:563
std::vector< LinkConnSpec >::const_iterator IT
PFCollection pfCandidatesFromPU_
Definition: PFPileUpAlgo.h:51
int chargedHadronVertex(const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
Definition: PFPileUpAlgo.cc:43
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:44
std::vector< edm::FwdPtr< reco::PFCandidate > > PFCollection
Definition: PFPileUpAlgo.h:16
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:35
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
void process(const PFCollection &pfCandidates, const reco::VertexCollection &vertices)
Definition: PFPileUpAlgo.cc:5
tuple cout
Definition: gather_cfg.py:121
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
virtual ParticleType particleId() const
Definition: PFCandidate.h:347