CMS 3D CMS Logo

VertexAssociatorByPositionAndTracks.cc
Go to the documentation of this file.
3 #include "CLHEP/Units/GlobalSystemOfUnits.h"
4 
6 
8  double absZ,
9  double sigmaZ,
10  double maxRecoZ,
11  double absT,
12  double sigmaT,
13  double maxRecoT,
14  double sharedTrackFraction,
15  const reco::RecoToSimCollection *trackRecoToSimAssociation,
16  const reco::SimToRecoCollection *trackSimToRecoAssociation):
17  productGetter_(productGetter),
18  absZ_(absZ),
19  sigmaZ_(sigmaZ),
20  maxRecoZ_(maxRecoZ),
21  absT_(absT),
22  sigmaT_(sigmaT),
23  maxRecoT_(maxRecoT),
24  sharedTrackFraction_(sharedTrackFraction),
25  trackRecoToSimAssociation_(trackRecoToSimAssociation),
26  trackSimToRecoAssociation_(trackSimToRecoAssociation)
27 {}
28 
30  double absZ,
31  double sigmaZ,
32  double maxRecoZ,
33  double sharedTrackFraction,
34  const reco::RecoToSimCollection *trackRecoToSimAssociation,
35  const reco::SimToRecoCollection *trackSimToRecoAssociation):
36  productGetter_(productGetter),
37  absZ_(absZ),
38  sigmaZ_(sigmaZ),
39  maxRecoZ_(maxRecoZ),
40  absT_(std::numeric_limits<double>::max()),
41  sigmaT_(std::numeric_limits<double>::max()),
42  maxRecoT_(std::numeric_limits<double>::max()),
43  sharedTrackFraction_(sharedTrackFraction),
44  trackRecoToSimAssociation_(trackRecoToSimAssociation),
45  trackSimToRecoAssociation_(trackSimToRecoAssociation)
46 {}
47 
49 
51  const edm::Handle<TrackingVertexCollection>& tVCH) const {
53 
54  const edm::View<reco::Vertex>& recoVertices = *vCH;
55  const TrackingVertexCollection& simVertices = *tVCH;
56 
57  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateRecoToSim(): associating "
58  << recoVertices.size() << " reco::Vertices to" << simVertices.size() << " TrackingVertices";
59 
60  for(size_t iReco=0; iReco != recoVertices.size(); ++iReco) {
61  const reco::Vertex& recoVertex = recoVertices[iReco];
62 
63  // skip fake vertices
64  if(std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() || recoVertex.ndof() < 0.)
65  continue;
66 
67  LogTrace("VertexAssociation") << " reco::Vertex at Z " << recoVertex.z();
68 
69  int current_event = -1;
70  for(size_t iSim=0; iSim != simVertices.size(); ++iSim) {
71  const TrackingVertex& simVertex = simVertices[iSim];
72 
73  // Associate only to primary vertices of the in-time pileup
74  // events (BX=0, first vertex in each of the events)
75  if(simVertex.eventId().bunchCrossing() != 0) continue;
76  if(simVertex.eventId().event() != current_event) {
77  current_event = simVertex.eventId().event();
78  }
79  else {
80  continue;
81  }
82 
83  LogTrace("VertexAssociation") << " Considering TrackingVertex at Z " << simVertex.position().z();
84 
85  // recoVertex.t() == 0. is a special value
86  // need to change this to std::numeric_limits<double>::max() or something more clear
87  const bool useTiming = ( absT_ != std::numeric_limits<double>::max() && recoVertex.t() != 0. );
88  if( useTiming ) {
89  LogTrace("VertexAssociation") << " and T " << recoVertex.t()*CLHEP::second << std::endl;
90  }
91 
92  const double tdiff = std::abs(recoVertex.t() - simVertex.position().t()*CLHEP::second);
93  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
94  if( zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_ &&
95  ( !useTiming || ( tdiff < absT_ && tdiff / recoVertex.tError() < sigmaT_ ) ) ) {
96  auto sharedTracks = calculateVertexSharedTracks(recoVertex, simVertex, *trackRecoToSimAssociation_);
97  auto fraction = double(sharedTracks)/recoVertex.tracksSize();
98  if(sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
99  LogTrace("VertexAssociation") << " Matched with significance " << zdiff/recoVertex.zError() << " " << tdiff/recoVertex.tError()
100  << " shared tracks " << sharedTracks << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles " << simVertex.nDaughterTracks();
101 
102  ret.insert(reco::VertexBaseRef(vCH, iReco), std::make_pair(TrackingVertexRef(tVCH, iSim), sharedTracks));
103  }
104  }
105  }
106  }
107 
108  ret.post_insert();
109 
110  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateRecoToSim(): finished";
111 
112  return ret;
113 }
114 
116  const edm::Handle<TrackingVertexCollection>& tVCH) const {
118 
119  const edm::View<reco::Vertex>& recoVertices = *vCH;
120  const TrackingVertexCollection& simVertices = *tVCH;
121 
122  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateSimToReco(): associating "
123  << simVertices.size() << " TrackingVertices to " << recoVertices.size() << " reco::Vertices";
124 
125  int current_event = -1;
126  for(size_t iSim=0; iSim != simVertices.size(); ++iSim) {
127  const TrackingVertex& simVertex = simVertices[iSim];
128 
129  // Associate only primary vertices of the in-time pileup
130  // events (BX=0, first vertex in each of the events)
131  if(simVertex.eventId().bunchCrossing() != 0) continue;
132  if(simVertex.eventId().event() != current_event) {
133  current_event = simVertex.eventId().event();
134  }
135  else {
136  continue;
137  }
138 
139  LogTrace("VertexAssociation") << " TrackingVertex at Z " << simVertex.position().z();
140 
141  for(size_t iReco=0; iReco != recoVertices.size(); ++iReco) {
142  const reco::Vertex& recoVertex = recoVertices[iReco];
143 
144  // skip fake vertices
145  if(std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() || recoVertex.ndof() < 0.)
146  continue;
147 
148  LogTrace("VertexAssociation") << " Considering reco::Vertex at Z " << recoVertex.z();
149  const bool useTiming = ( absT_ != std::numeric_limits<double>::max() && recoVertex.t() != 0. );
150  if( useTiming ) {
151  LogTrace("VertexAssociation") << " and T " << recoVertex.t()*CLHEP::second << std::endl;
152  }
153 
154  const double tdiff = std::abs(recoVertex.t() - simVertex.position().t()*CLHEP::second);
155  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
156  if( zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_ &&
157  ( !useTiming || ( tdiff < absT_ && tdiff / recoVertex.tError() < sigmaT_ ) ) ) {
158  auto sharedTracks = calculateVertexSharedTracks(simVertex, recoVertex, *trackSimToRecoAssociation_);
159  auto fraction = double(sharedTracks)/recoVertex.tracksSize();
160  if(sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
161  LogTrace("VertexAssociation") << " Matched with significance " << zdiff/recoVertex.zError() << " " << tdiff/recoVertex.tError()
162  << " shared tracks " << sharedTracks << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles " << simVertex.nDaughterTracks();
163 
164  ret.insert(TrackingVertexRef(tVCH, iSim), std::make_pair(reco::VertexBaseRef(vCH, iReco), sharedTracks));
165  }
166  }
167  }
168  }
169 
170  ret.post_insert();
171 
172  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateSimToReco(): finished";
173 
174  return ret;
175 }
#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:123
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:68
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
U second(std::pair< T, U > const &p)
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
z coordinate
Definition: Vertex.h:115
VertexAssociatorByPositionAndTracks(const edm::EDProductGetter *productGetter, double absZ, double sigmaZ, double maxRecoZ, double absT, double sigmaT, double maxRecoT, double sharedTrackFraction, const reco::RecoToSimCollection *trackRecoToSimAssociation, const reco::SimToRecoCollection *trackSimToRecoAssociation)
#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:105
edm::Ref< TrackingVertexCollection > TrackingVertexRef
bool isFake() const
Definition: Vertex.h:72
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_
double tError() const
error on t
Definition: Vertex.h:125
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
double t() const
t coordinate
Definition: Vertex.h:117
const LorentzVector & position() const