CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VertexAssociatorByPositionAndTracks.cc
Go to the documentation of this file.
3 
5 
7  double absZ,
8  double sigmaZ,
9  double maxRecoZ,
10  double sharedTrackFraction,
11  const reco::RecoToSimCollection *trackRecoToSimAssociation,
12  const reco::SimToRecoCollection *trackSimToRecoAssociation):
13  productGetter_(productGetter),
14  absZ_(absZ),
15  sigmaZ_(sigmaZ),
16  maxRecoZ_(maxRecoZ),
17  sharedTrackFraction_(sharedTrackFraction),
18  trackRecoToSimAssociation_(trackRecoToSimAssociation),
19  trackSimToRecoAssociation_(trackSimToRecoAssociation)
20 {}
21 
23 
25  const edm::Handle<TrackingVertexCollection>& tVCH) const {
27 
28  const edm::View<reco::Vertex>& recoVertices = *vCH;
29  const TrackingVertexCollection& simVertices = *tVCH;
30 
31  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateRecoToSim(): associating "
32  << recoVertices.size() << " reco::Vertices to" << simVertices.size() << " TrackingVertices";
33 
34  for(size_t iReco=0; iReco != recoVertices.size(); ++iReco) {
35  const reco::Vertex& recoVertex = recoVertices[iReco];
36 
37  // skip fake vertices
38  if(std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() || recoVertex.ndof() < 0.)
39  continue;
40 
41  LogTrace("VertexAssociation") << " reco::Vertex at Z " << recoVertex.z();
42 
43  int current_event = -1;
44  for(size_t iSim=0; iSim != simVertices.size(); ++iSim) {
45  const TrackingVertex& simVertex = simVertices[iSim];
46 
47  // Associate only to primary vertices of the in-time pileup
48  // events (BX=0, first vertex in each of the events)
49  if(simVertex.eventId().bunchCrossing() != 0) continue;
50  if(simVertex.eventId().event() != current_event) {
51  current_event = simVertex.eventId().event();
52  }
53  else {
54  continue;
55  }
56 
57  LogTrace("VertexAssociation") << " Considering TrackingVertex at Z " << simVertex.position().z();
58 
59  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
60  if(zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_) {
61  auto sharedTracks = calculateVertexSharedTracks(recoVertex, simVertex, *trackRecoToSimAssociation_);
62  auto fraction = double(sharedTracks)/recoVertex.tracksSize();
63  if(sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
64  LogTrace("VertexAssociation") << " Matched with significance " << zdiff/recoVertex.zError()
65  << " shared tracks " << sharedTracks << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles " << simVertex.nDaughterTracks();
66 
67  ret.insert(reco::VertexBaseRef(vCH, iReco), std::make_pair(TrackingVertexRef(tVCH, iSim), sharedTracks));
68  }
69  }
70  }
71  }
72 
73  ret.post_insert();
74 
75  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateRecoToSim(): finished";
76 
77  return ret;
78 }
79 
81  const edm::Handle<TrackingVertexCollection>& tVCH) const {
83 
84  const edm::View<reco::Vertex>& recoVertices = *vCH;
85  const TrackingVertexCollection& simVertices = *tVCH;
86 
87  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateSimToReco(): associating "
88  << simVertices.size() << " TrackingVertices to " << recoVertices.size() << " reco::Vertices";
89 
90  int current_event = -1;
91  for(size_t iSim=0; iSim != simVertices.size(); ++iSim) {
92  const TrackingVertex& simVertex = simVertices[iSim];
93 
94  // Associate only primary vertices of the in-time pileup
95  // events (BX=0, first vertex in each of the events)
96  if(simVertex.eventId().bunchCrossing() != 0) continue;
97  if(simVertex.eventId().event() != current_event) {
98  current_event = simVertex.eventId().event();
99  }
100  else {
101  continue;
102  }
103 
104  LogTrace("VertexAssociation") << " TrackingVertex at Z " << simVertex.position().z();
105 
106  for(size_t iReco=0; iReco != recoVertices.size(); ++iReco) {
107  const reco::Vertex& recoVertex = recoVertices[iReco];
108 
109  // skip fake vertices
110  if(std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() || recoVertex.ndof() < 0.)
111  continue;
112 
113  LogTrace("VertexAssociation") << " Considering reco::Vertex at Z " << recoVertex.z();
114 
115  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
116  if(zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_) {
117  auto sharedTracks = calculateVertexSharedTracks(simVertex, recoVertex, *trackSimToRecoAssociation_);
118  auto fraction = double(sharedTracks)/recoVertex.tracksSize();
119  if(sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
120  LogTrace("VertexAssociation") << " Matched with significance " << zdiff/recoVertex.zError()
121  << " shared tracks " << sharedTracks << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles " << simVertex.nDaughterTracks();
122 
123  ret.insert(TrackingVertexRef(tVCH, iSim), std::make_pair(reco::VertexBaseRef(vCH, iReco), sharedTracks));
124  }
125  }
126  }
127  }
128 
129  ret.post_insert();
130 
131  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateSimToReco(): finished";
132 
133  return ret;
134 }
#define LogDebug(id)
int event() const
get the contents of the subdetector field (should be protected?)
double zError() const
error on z
Definition: Vertex.h:118
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:60
size_type size() const
const reco::RecoToSimCollection * trackRecoToSimAssociation_
EDProductGetter const * productGetter(std::atomic< void const * > const &iCache)
virtual reco::VertexSimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Vertex > > &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const
compare reco to sim the handle of reco::Vertex and TrackingVertex collections
VertexAssociatorByPositionAndTracks(const edm::EDProductGetter *productGetter, double absZ, double sigmaZ, double maxRecoZ, double sharedTrackFraction, const reco::RecoToSimCollection *trackRecoToSimAssociation, const reco::SimToRecoCollection *trackSimToRecoAssociation)
void post_insert()
post insert action
int bunchCrossing() const
get the detector field from this detid
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
y coordinate
Definition: Vertex.h:112
#define LogTrace(id)
virtual reco::VertexRecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Vertex > > &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const
compare reco to sim the handle of reco::Vertex and TrackingVertex collections
double ndof() const
Definition: Vertex.h:102
edm::Ref< TrackingVertexCollection > TrackingVertexRef
bool isFake() const
Definition: Vertex.h:64
unsigned int nDaughterTracks() const
void insert(const key_type &k, const data_type &v)
insert an association
std::vector< TrackingVertex > TrackingVertexCollection
const EncodedEventId & eventId() const
unsigned int calculateVertexSharedTracks(const reco::Vertex &recoV, const TrackingVertex &simV, const reco::RecoToSimCollection &trackRecoToSimAssociation)
const reco::SimToRecoCollection * trackSimToRecoAssociation_
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:34
const LorentzVector & position() const