CMS 3D CMS Logo

calculateVertexSharedTracks.cc
Go to the documentation of this file.
2 
3 unsigned int calculateVertexSharedTracks(const reco::Vertex &recoV,
4  const TrackingVertex &simV,
5  const reco::RecoToSimCollection &trackRecoToSimAssociation) {
6  unsigned int sharedTracks = 0;
7  for (auto iTrack = recoV.tracks_begin(); iTrack != recoV.tracks_end(); ++iTrack) {
8  auto found = trackRecoToSimAssociation.find(*iTrack);
9 
10  if (found == trackRecoToSimAssociation.end())
11  continue;
12 
13  // matched TP equal to any TP of sim vertex => increase counter
14  for (const auto &tp : found->val) {
15  if (std::find_if(simV.daughterTracks_begin(), simV.daughterTracks_end(), [&](const TrackingParticleRef &vtp) {
16  return tp.first == vtp;
17  }) != simV.daughterTracks_end()) {
18  sharedTracks += 1;
19  break;
20  }
21  }
22  }
23 
24  return sharedTracks;
25 }
26 
28  const reco::Vertex &recoV,
29  const reco::SimToRecoCollection &trackSimToRecoAssociation) {
30  unsigned int sharedTracks = 0;
31  for (auto iTP = simV.daughterTracks_begin(); iTP != simV.daughterTracks_end(); ++iTP) {
32  auto found = trackSimToRecoAssociation.find(*iTP);
33 
34  if (found == trackSimToRecoAssociation.end())
35  continue;
36 
37  // matched track equal to any track of reco vertex => increase counter
38  for (const auto &tk : found->val) {
39  if (std::find_if(recoV.tracks_begin(), recoV.tracks_end(), [&](const reco::TrackBaseRef &vtk) {
40  return tk.first == vtk;
41  }) != recoV.tracks_end()) {
42  sharedTracks += 1;
43  break;
44  }
45  }
46  }
47 
48  return sharedTracks;
49 }
tp_iterator daughterTracks_begin() const
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
const_iterator end() const
last iterator over the map (read only)
unsigned int calculateVertexSharedTracks(const reco::Vertex &recoV, const TrackingVertex &simV, const reco::RecoToSimCollection &trackRecoToSimAssociation)
const_iterator find(const key_type &k) const
find element with specified reference key
tp_iterator daughterTracks_end() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76