CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
VertexAssociatorByPositionAndTracks Class Reference

#include <VertexAssociatorByPositionAndTracks.h>

Inheritance diagram for VertexAssociatorByPositionAndTracks:
reco::VertexToTrackingVertexAssociatorBaseImpl

Public Member Functions

reco::VertexRecoToSimCollection associateRecoToSim (const edm::Handle< edm::View< reco::Vertex >> &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const override
 
reco::VertexSimToRecoCollection associateSimToReco (const edm::Handle< edm::View< reco::Vertex >> &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const override
 
 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)
 
 VertexAssociatorByPositionAndTracks (const edm::EDProductGetter *productGetter, double absZ, double sigmaZ, double maxRecoZ, double sharedTrackFraction, const reco::RecoToSimCollection *trackRecoToSimAssociation, const reco::SimToRecoCollection *trackSimToRecoAssociation)
 
 ~VertexAssociatorByPositionAndTracks () override
 
- Public Member Functions inherited from reco::VertexToTrackingVertexAssociatorBaseImpl
 VertexToTrackingVertexAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~VertexToTrackingVertexAssociatorBaseImpl ()
 Destructor. More...
 

Private Attributes

const double absT_
 
const double absZ_
 
const double maxRecoT_
 
const double maxRecoZ_
 
const edm::EDProductGetterproductGetter_
 
const double sharedTrackFraction_
 
const double sigmaT_
 
const double sigmaZ_
 
const reco::RecoToSimCollectiontrackRecoToSimAssociation_
 
const reco::SimToRecoCollectiontrackSimToRecoAssociation_
 

Detailed Description

This class associates reco::Vertices and TrackingVertices by their position (maximum distance in Z should be smaller than absZ and sigmaZ*zError of reco::Vertex), and (optionally) by the fraction of tracks shared by reco::Vertex and TrackingVertex divided by the number of tracks in reco::Vertex. This fraction is always used as the quality in the association, i.e. multiple associations are sorted by it in descending order.

Definition at line 16 of file VertexAssociatorByPositionAndTracks.h.

Constructor & Destructor Documentation

VertexAssociatorByPositionAndTracks::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 
)
VertexAssociatorByPositionAndTracks::VertexAssociatorByPositionAndTracks ( const edm::EDProductGetter productGetter,
double  absZ,
double  sigmaZ,
double  maxRecoZ,
double  sharedTrackFraction,
const reco::RecoToSimCollection trackRecoToSimAssociation,
const reco::SimToRecoCollection trackSimToRecoAssociation 
)
VertexAssociatorByPositionAndTracks::~VertexAssociatorByPositionAndTracks ( )
override

Definition at line 48 of file VertexAssociatorByPositionAndTracks.cc.

48 {}

Member Function Documentation

reco::VertexRecoToSimCollection VertexAssociatorByPositionAndTracks::associateRecoToSim ( const edm::Handle< edm::View< reco::Vertex >> &  vCH,
const edm::Handle< TrackingVertexCollection > &  tVCH 
) const
overridevirtual

compare reco to sim the handle of reco::Vertex and TrackingVertex collections

Implements reco::VertexToTrackingVertexAssociatorBaseImpl.

Definition at line 50 of file VertexAssociatorByPositionAndTracks.cc.

References funct::abs(), absT_, absZ_, EncodedEventId::bunchCrossing(), calculateVertexSharedTracks(), EncodedEventId::event(), TrackingVertex::eventId(), HLT_2018_cff::fraction, edm::AssociationMap< Tag >::insert(), reco::Vertex::isFake(), reco::Vertex::isValid(), LogDebug, LogTrace, SiStripPI::max, maxRecoZ_, TrackingVertex::nDaughterTracks(), reco::Vertex::ndof(), TrackingVertex::position(), edm::AssociationMap< Tag >::post_insert(), productGetter_, ecalDetailedTimeRecHit_cfi::recoVertex, runTheMatrix::ret, edm::second(), sharedTrackFraction_, sigmaT_, sigmaZ_, ecalDetailedTimeRecHit_cfi::simVertex, HGCalValidator_cfi::simVertices, edm::View< T >::size(), reco::Vertex::t(), reco::Vertex::tError(), trackRecoToSimAssociation_, reco::Vertex::tracksSize(), RecoVertex_phase2_timing_cff::useTiming, reco::Vertex::z(), and reco::Vertex::zError().

51  {
53 
54  const edm::View<reco::Vertex> &recoVertices = *vCH;
56 
57  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::"
58  "associateRecoToSim(): associating "
59  << recoVertices.size() << " reco::Vertices to" << simVertices.size()
60  << " TrackingVertices";
61 
62  // filter sim PVs
63  std::vector<size_t> simPVindices;
64  simPVindices.reserve(recoVertices.size());
65  {
66  int current_event = -1;
67  for (size_t iSim = 0; iSim != simVertices.size(); ++iSim) {
68  const TrackingVertex &simVertex = simVertices[iSim];
69 
70  // Associate only to primary vertices of the in-time pileup
71  // events (BX=0, first vertex in each of the events)
72  if (simVertex.eventId().bunchCrossing() != 0)
73  continue;
74  if (simVertex.eventId().event() != current_event) {
75  current_event = simVertex.eventId().event();
76  simPVindices.push_back(iSim);
77  }
78  }
79  }
80 
81  for (size_t iReco = 0; iReco != recoVertices.size(); ++iReco) {
82  const reco::Vertex &recoVertex = recoVertices[iReco];
83 
84  // skip fake vertices
85  if (std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() || recoVertex.ndof() < 0.)
86  continue;
87 
88  LogTrace("VertexAssociation") << " reco::Vertex at Z " << recoVertex.z();
89 
90  for (const size_t iSim : simPVindices) {
91  const TrackingVertex &simVertex = simVertices[iSim];
92  LogTrace("VertexAssociation") << " Considering TrackingVertex at Z " << simVertex.position().z();
93 
94  // recoVertex.t() == 0. is a special value
95  // need to change this to std::numeric_limits<double>::max() or something
96  // more clear
97  const bool useTiming = (absT_ != std::numeric_limits<double>::max() && recoVertex.t() != 0.);
98  if (useTiming) {
99  LogTrace("VertexAssociation") << " and T " << recoVertex.t() * CLHEP::second << std::endl;
100  }
101 
102  const double tdiff = std::abs(recoVertex.t() - simVertex.position().t() * CLHEP::second);
103  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
104  if (zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_ &&
105  (!useTiming || (tdiff < absT_ && tdiff / recoVertex.tError() < sigmaT_))) {
106  auto sharedTracks = calculateVertexSharedTracks(recoVertex, simVertex, *trackRecoToSimAssociation_);
107  auto fraction = double(sharedTracks) / recoVertex.tracksSize();
108  if (sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
109  LogTrace("VertexAssociation") << " Matched with significance " << zdiff / recoVertex.zError() << " "
110  << tdiff / recoVertex.tError() << " shared tracks " << sharedTracks
111  << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles "
112  << simVertex.nDaughterTracks();
113 
114  ret.insert(reco::VertexBaseRef(vCH, iReco), std::make_pair(TrackingVertexRef(tVCH, iSim), sharedTracks));
115  }
116  }
117  }
118  }
119 
120  ret.post_insert();
121 
122  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateRecoToSim(): finished";
123 
124  return ret;
125 }
#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:127
ret
prodAgent to be discontinued
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:71
size_type size() const
const reco::RecoToSimCollection * trackRecoToSimAssociation_
U second(std::pair< T, U > const &p)
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:119
#define LogTrace(id)
double ndof() const
Definition: Vertex.h:109
edm::Ref< TrackingVertexCollection > TrackingVertexRef
bool isFake() const
Definition: Vertex.h:75
unsigned int nDaughterTracks() const
std::vector< TrackingVertex > TrackingVertexCollection
const EncodedEventId & eventId() const
unsigned int calculateVertexSharedTracks(const reco::Vertex &recoV, const TrackingVertex &simV, const reco::RecoToSimCollection &trackRecoToSimAssociation)
double tError() const
error on t
Definition: Vertex.h:129
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
double t() const
t coordinate
Definition: Vertex.h:121
const LorentzVector & position() const
reco::VertexSimToRecoCollection VertexAssociatorByPositionAndTracks::associateSimToReco ( const edm::Handle< edm::View< reco::Vertex >> &  vCH,
const edm::Handle< TrackingVertexCollection > &  tVCH 
) const
overridevirtual

compare reco to sim the handle of reco::Vertex and TrackingVertex collections

Implements reco::VertexToTrackingVertexAssociatorBaseImpl.

Definition at line 127 of file VertexAssociatorByPositionAndTracks.cc.

References funct::abs(), absT_, absZ_, EncodedEventId::bunchCrossing(), calculateVertexSharedTracks(), EncodedEventId::event(), TrackingVertex::eventId(), HLT_2018_cff::fraction, edm::AssociationMap< Tag >::insert(), reco::Vertex::isFake(), reco::Vertex::isValid(), LogDebug, LogTrace, SiStripPI::max, maxRecoZ_, TrackingVertex::nDaughterTracks(), reco::Vertex::ndof(), TrackingVertex::position(), edm::AssociationMap< Tag >::post_insert(), productGetter_, ecalDetailedTimeRecHit_cfi::recoVertex, runTheMatrix::ret, edm::second(), sharedTrackFraction_, sigmaT_, sigmaZ_, ecalDetailedTimeRecHit_cfi::simVertex, HGCalValidator_cfi::simVertices, edm::View< T >::size(), reco::Vertex::t(), reco::Vertex::tError(), trackSimToRecoAssociation_, reco::Vertex::tracksSize(), RecoVertex_phase2_timing_cff::useTiming, reco::Vertex::z(), and reco::Vertex::zError().

128  {
130 
131  const edm::View<reco::Vertex> &recoVertices = *vCH;
132  const TrackingVertexCollection &simVertices = *tVCH;
133 
134  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::"
135  "associateSimToReco(): associating "
136  << simVertices.size() << " TrackingVertices to " << recoVertices.size()
137  << " reco::Vertices";
138 
139  int current_event = -1;
140  for (size_t iSim = 0; iSim != simVertices.size(); ++iSim) {
141  const TrackingVertex &simVertex = simVertices[iSim];
142 
143  // Associate only primary vertices of the in-time pileup
144  // events (BX=0, first vertex in each of the events)
145  if (simVertex.eventId().bunchCrossing() != 0)
146  continue;
147  if (simVertex.eventId().event() != current_event) {
148  current_event = simVertex.eventId().event();
149  } else {
150  continue;
151  }
152 
153  LogTrace("VertexAssociation") << " TrackingVertex at Z " << simVertex.position().z();
154 
155  for (size_t iReco = 0; iReco != recoVertices.size(); ++iReco) {
156  const reco::Vertex &recoVertex = recoVertices[iReco];
157 
158  // skip fake vertices
159  if (std::abs(recoVertex.z()) > maxRecoZ_ || recoVertex.isFake() || !recoVertex.isValid() ||
160  recoVertex.ndof() < 0.)
161  continue;
162 
163  LogTrace("VertexAssociation") << " Considering reco::Vertex at Z " << recoVertex.z();
164  const bool useTiming = (absT_ != std::numeric_limits<double>::max() && recoVertex.t() != 0.);
165  if (useTiming) {
166  LogTrace("VertexAssociation") << " and T " << recoVertex.t() * CLHEP::second << std::endl;
167  }
168 
169  const double tdiff = std::abs(recoVertex.t() - simVertex.position().t() * CLHEP::second);
170  const double zdiff = std::abs(recoVertex.z() - simVertex.position().z());
171  if (zdiff < absZ_ && zdiff / recoVertex.zError() < sigmaZ_ &&
172  (!useTiming || (tdiff < absT_ && tdiff / recoVertex.tError() < sigmaT_))) {
173  auto sharedTracks = calculateVertexSharedTracks(simVertex, recoVertex, *trackSimToRecoAssociation_);
174  auto fraction = double(sharedTracks) / recoVertex.tracksSize();
175  if (sharedTrackFraction_ < 0 || fraction > sharedTrackFraction_) {
176  LogTrace("VertexAssociation") << " Matched with significance " << zdiff / recoVertex.zError() << " "
177  << tdiff / recoVertex.tError() << " shared tracks " << sharedTracks
178  << " reco Tracks " << recoVertex.tracksSize() << " TrackingParticles "
179  << simVertex.nDaughterTracks();
180 
181  ret.insert(TrackingVertexRef(tVCH, iSim), std::make_pair(reco::VertexBaseRef(vCH, iReco), sharedTracks));
182  }
183  }
184  }
185  }
186 
187  ret.post_insert();
188 
189  LogDebug("VertexAssociation") << "VertexAssociatorByPositionAndTracks::associateSimToReco(): finished";
190 
191  return ret;
192 }
#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:127
ret
prodAgent to be discontinued
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:71
size_type size() const
U second(std::pair< T, U > const &p)
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:119
#define LogTrace(id)
double ndof() const
Definition: Vertex.h:109
edm::Ref< TrackingVertexCollection > TrackingVertexRef
bool isFake() const
Definition: Vertex.h:75
unsigned int nDaughterTracks() const
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:129
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
double t() const
t coordinate
Definition: Vertex.h:121
const LorentzVector & position() const

Member Data Documentation

const double VertexAssociatorByPositionAndTracks::absT_
private

Definition at line 53 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByPositionAndTracks::absZ_
private

Definition at line 50 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByPositionAndTracks::maxRecoT_
private

Definition at line 55 of file VertexAssociatorByPositionAndTracks.h.

const double VertexAssociatorByPositionAndTracks::maxRecoZ_
private

Definition at line 52 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const edm::EDProductGetter* VertexAssociatorByPositionAndTracks::productGetter_
private

Definition at line 48 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByPositionAndTracks::sharedTrackFraction_
private

Definition at line 56 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByPositionAndTracks::sigmaT_
private

Definition at line 54 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByPositionAndTracks::sigmaZ_
private

Definition at line 51 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const reco::RecoToSimCollection* VertexAssociatorByPositionAndTracks::trackRecoToSimAssociation_
private

Definition at line 58 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateRecoToSim().

const reco::SimToRecoCollection* VertexAssociatorByPositionAndTracks::trackSimToRecoAssociation_
private

Definition at line 59 of file VertexAssociatorByPositionAndTracks.h.

Referenced by associateSimToReco().