CMS 3D CMS Logo

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

#include <VertexAssociatorByTracks.h>

Inheritance diagram for VertexAssociatorByTracks:
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
 
 VertexAssociatorByTracks (const edm::EDProductGetter *productGetter, double R2SMatchedSimRatio, double R2SMatchedRecoRatio, double S2RMatchedSimRatio, double S2RMatchedRecoRatio, const TrackingParticleSelector *selector, reco::TrackBase::TrackQuality trackQuality, const reco::RecoToSimCollection *trackRecoToSimAssociation, const reco::SimToRecoCollection *trackSimToRecoAssociation)
 
 ~VertexAssociatorByTracks () override
 
- Public Member Functions inherited from reco::VertexToTrackingVertexAssociatorBaseImpl
 VertexToTrackingVertexAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~VertexToTrackingVertexAssociatorBaseImpl ()
 Destructor. More...
 

Private Attributes

const edm::EDProductGetterproductGetter_
 
const double R2SMatchedRecoRatio_
 
const double R2SMatchedSimRatio_
 
const double S2RMatchedRecoRatio_
 
const double S2RMatchedSimRatio_
 
const TrackingParticleSelectorselector_
 
const reco::TrackBase::TrackQuality trackQuality_
 
const reco::RecoToSimCollectiontrackRecoToSimAssociation_
 
const reco::SimToRecoCollectiontrackSimToRecoAssociation_
 

Detailed Description

Definition at line 9 of file VertexAssociatorByTracks.h.

Constructor & Destructor Documentation

VertexAssociatorByTracks::VertexAssociatorByTracks ( const edm::EDProductGetter productGetter,
double  R2SMatchedSimRatio,
double  R2SMatchedRecoRatio,
double  S2RMatchedSimRatio,
double  S2RMatchedRecoRatio,
const TrackingParticleSelector selector,
reco::TrackBase::TrackQuality  trackQuality,
const reco::RecoToSimCollection trackRecoToSimAssociation,
const reco::SimToRecoCollection trackSimToRecoAssociation 
)

Definition at line 21 of file VertexAssociatorByTracks.cc.

30  : productGetter_(productGetter),
31  R2SMatchedSimRatio_(R2SMatchedSimRatio),
32  R2SMatchedRecoRatio_(R2SMatchedRecoRatio),
33  S2RMatchedSimRatio_(S2RMatchedSimRatio),
34  S2RMatchedRecoRatio_(S2RMatchedRecoRatio),
35  selector_(selector),
37  trackRecoToSimAssociation_(trackRecoToSimAssociation),
38  trackSimToRecoAssociation_(trackSimToRecoAssociation) {}
const TrackingParticleSelector * selector_
const reco::RecoToSimCollection * trackRecoToSimAssociation_
const reco::TrackBase::TrackQuality trackQuality_
const reco::SimToRecoCollection * trackSimToRecoAssociation_
const edm::EDProductGetter * productGetter_
VertexAssociatorByTracks::~VertexAssociatorByTracks ( )
override

Definition at line 41 of file VertexAssociatorByTracks.cc.

41 {}

Member Function Documentation

reco::VertexRecoToSimCollection VertexAssociatorByTracks::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 43 of file VertexAssociatorByTracks.cc.

References TrackValidation_cff::association, genericTrackCandidates_cff::associations, edm::AssociationMap< Tag >::insert(), LogDebug, LogTrace, genParticles_cff::map, match(), patCandidatesForDimuonsSequences_cff::matches, edm::AssociationMap< Tag >::numberOfAssociations(), productGetter_, qcdUeDQM_cfi::quality, R2SMatchedRecoRatio_, R2SMatchedSimRatio_, ecalDetailedTimeRecHit_cfi::recoVertex, selector_, edm::AssociationMap< Tag >::size(), trackQuality_, trackRecoToSimAssociation_, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), and reco::Vertex::trackWeight().

45  {
47 
48  std::map<TrackingVertexRef, std::pair<double, std::size_t>> matches;
49 
50  LogDebug("VertexAssociation") << "reco::VertexCollection size = " << recoVertexes->size()
51  << " ; TrackingVertexCollection size = " << trackingVertexes->size() << std::endl;
52 
53  // Loop over RecoVertex
54  for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex) {
55  reco::VertexBaseRef recoVertex(recoVertexes, recoIndex);
56 
57  matches.clear();
58 
59  double recoDaughterWeight = 0.;
60 
61  // Loop over daughter tracks of RecoVertex
62  for (reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
63  recoDaughter != recoVertex->tracks_end();
64  ++recoDaughter) {
65  // Check the quality of the RecoDoughters
66  if (!(*recoDaughter)->quality(trackQuality_))
67  continue;
68 
69  // Check for association for the given RecoDaughter
70  if (trackRecoToSimAssociation_->numberOfAssociations(*recoDaughter) > 0) {
71  std::vector<std::pair<TrackingParticleRef, double>> associations = (*trackRecoToSimAssociation_)[*recoDaughter];
72 
73  // Loop over TrackingParticles associated with RecoDaughter
74  for (std::vector<std::pair<TrackingParticleRef, double>>::const_iterator association = associations.begin();
75  association != associations.end();
76  ++association) {
77  // Get a reference to parent vertex of TrackingParticle associated to
78  // the RecoDaughter
79  TrackingVertexRef trackingVertex = association->first->parentVertex();
80  // Store matched RecoDaughter to the trackingVertex
81  matches[trackingVertex].first += recoVertex->trackWeight(*recoDaughter);
82  matches[trackingVertex].second++;
83  }
84  }
85 
86  recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
87  } // loop over daughter tracks of RecoVertex
88 
89  std::size_t assoIndex = 0;
90 
91  // Loop over map between TrackingVertexes and matched RecoDaugther
92  for (std::map<TrackingVertexRef, std::pair<double, std::size_t>>::const_iterator match = matches.begin();
93  match != matches.end();
94  ++match) {
95  // Getting the TrackingVertex information
96  TrackingVertexRef trackingVertex = match->first;
97  double matchedDaughterWeight = match->second.first;
98  std::size_t matchedDaughterCounter = match->second.second;
99 
100  // Count for only reconstructible SimDaughters
101  std::size_t simDaughterCounter = 0;
102 
103  for (TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
104  simDaughter != trackingVertex->daughterTracks_end();
105  ++simDaughter)
106  if ((*selector_)(**simDaughter))
107  simDaughterCounter++;
108 
109  // Sanity condition in case that reconstructable condition is too tight
110  if (simDaughterCounter < matchedDaughterCounter)
111  simDaughterCounter = matchedDaughterCounter;
112 
113  // Condition over S2RMatchedSimRatio
114  if ((double)matchedDaughterCounter / simDaughterCounter < R2SMatchedSimRatio_)
115  continue;
116 
117  double quality = (double)matchedDaughterWeight / recoDaughterWeight;
118 
119  // Condition over R2SMatchedRecoRatio
120  if (quality < R2SMatchedRecoRatio_)
121  continue;
122 
123  outputCollection.insert(recoVertex, std::make_pair(trackingVertex, quality));
124 
125  LogTrace("VertexAssociation") << "R2S: INDEX " << assoIndex << " ; SimDaughterCounter = " << simDaughterCounter
126  << " ; RecoDaughterWeight = " << recoDaughterWeight
127  << " ; MatchedDaughterCounter = " << matchedDaughterCounter
128  << " ; MatchedDaughterWeight = " << matchedDaughterWeight
129  << " ; quality = " << quality;
130 
131  assoIndex++;
132  }
133  } // Loop on RecoVertex
134 
135  LogTrace("VertexAssociation") << "\nRecoToSim OUTPUT COLLECTION: outputCollection.size() = "
136  << outputCollection.size() << std::endl;
137 
138  return outputCollection;
139 }
#define LogDebug(id)
const TrackingParticleSelector * selector_
size_type numberOfAssociations(const key_type &k) const
number of associations to a key
const reco::RecoToSimCollection * trackRecoToSimAssociation_
#define LogTrace(id)
const reco::TrackBase::TrackQuality trackQuality_
const edm::EDProductGetter * productGetter_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
reco::VertexSimToRecoCollection VertexAssociatorByTracks::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 141 of file VertexAssociatorByTracks.cc.

References TrackValidation_cff::association, genericTrackCandidates_cff::associations, edm::RefToBase< T >::id(), edm::AssociationMap< Tag >::insert(), edm::RefToBase< T >::key(), LogTrace, genParticles_cff::map, match(), patCandidatesForDimuonsSequences_cff::matches, edm::AssociationMap< Tag >::numberOfAssociations(), productGetter_, qcdUeDQM_cfi::quality, ecalDetailedTimeRecHit_cfi::recoVertex, S2RMatchedRecoRatio_, S2RMatchedSimRatio_, selector_, edm::AssociationMap< Tag >::size(), trackQuality_, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), trackSimToRecoAssociation_, and reco::Vertex::trackWeight().

143  {
145 
146  // Loop over TrackingVertexes
147  std::map<std::size_t, std::pair<double, std::size_t>> matches;
148 
149  // Loop over TrackingVertexes
150  for (std::size_t simIndex = 0; simIndex < trackingVertexes->size(); ++simIndex) {
151  TrackingVertexRef trackingVertex(trackingVertexes, simIndex);
152 
153  matches.clear();
154 
155  std::size_t simDaughterCounter = 0;
156 
157  // Loop over daughter tracks of TrackingVertex
158  for (TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
159  simDaughter != trackingVertex->daughterTracks_end();
160  ++simDaughter) {
161  // Select only reconstructible SimDaughters
162  if (!(*selector_)(**simDaughter))
163  continue;
164 
165  // Check for association for the given RecoDaughter
166  if (trackSimToRecoAssociation_->numberOfAssociations(*simDaughter) > 0) {
167  std::vector<std::pair<reco::TrackBaseRef, double>> associations = (*trackSimToRecoAssociation_)[*simDaughter];
168 
169  // Loop over RecoTracks associated with TrackingParticle
170  for (std::vector<std::pair<reco::TrackBaseRef, double>>::const_iterator association = associations.begin();
171  association != associations.end();
172  ++association) {
173  reco::TrackBaseRef recoTrack = association->first;
174 
175  for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex) {
176  reco::VertexBaseRef recoVertex(recoVertexes, recoIndex);
177 
178  for (reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
179  recoDaughter != recoVertex->tracks_end();
180  ++recoDaughter) {
181  // Store matched RecoDaughter to the RecoVertex
182  if (recoDaughter->id() == recoTrack.id() && recoDaughter->key() == recoTrack.key()) {
183  matches[recoIndex].first += recoVertex->trackWeight(*recoDaughter);
184  matches[recoIndex].second++;
185  }
186  }
187  }
188  } // loop a recotracks
189  }
190 
191  simDaughterCounter++;
192  } // loop a simDaughter
193 
194  std::size_t assoIndex = 0;
195 
196  // Loop over map, set score, add to outputCollection
197  for (std::map<std::size_t, std::pair<double, std::size_t>>::const_iterator match = matches.begin();
198  match != matches.end();
199  ++match) {
200  // Getting the TrackingVertex information
201  reco::VertexBaseRef recoVertex(recoVertexes, match->first);
202  double matchedDaughterWeight = match->second.first;
203  std::size_t matchedDaughterCounter = match->second.second;
204 
205  double recoDaughterWeight = 0.;
206 
207  // Weighted count of those tracks with a given quality of the
208  // RecoDoughters
209  for (reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
210  recoDaughter != recoVertex->tracks_end();
211  ++recoDaughter)
212  if ((*recoDaughter)->quality(trackQuality_))
213  recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
214 
215  // Sanity condition in case that track quality condition is too tight
216  if (recoDaughterWeight < matchedDaughterWeight)
217  recoDaughterWeight = matchedDaughterWeight;
218 
219  // Condition over S2RMatchedRecoRatio
220  if (matchedDaughterWeight / recoDaughterWeight < S2RMatchedRecoRatio_)
221  continue;
222 
223  double quality = (double)matchedDaughterCounter / simDaughterCounter;
224 
225  // Condition over S2RMatchedSimRatio
226  if (quality < S2RMatchedSimRatio_)
227  continue;
228 
229  outputCollection.insert(trackingVertex, std::make_pair(recoVertex, quality));
230 
231  LogTrace("VertexAssociation") << "R2S: INDEX " << assoIndex << " ; SimDaughterCounter = " << simDaughterCounter
232  << " ; RecoDaughterWeight = " << recoDaughterWeight
233  << " ; MatchedDaughterCounter = " << matchedDaughterCounter
234  << " ; MatchedDaughterWeight = " << matchedDaughterWeight
235  << " ; quality = " << quality;
236 
237  assoIndex++;
238  }
239  } // loop over TrackingVertexes
240 
241  LogTrace("VertexAssociation") << "SimToReco OUTPUT COLLECTION: outputCollection.size() = " << outputCollection.size()
242  << std::endl;
243 
244  return outputCollection;
245 }
const TrackingParticleSelector * selector_
ProductID id() const
Definition: RefToBase.h:214
size_type numberOfAssociations(const key_type &k) const
number of associations to a key
size_t key() const
Definition: RefToBase.h:219
#define LogTrace(id)
const reco::TrackBase::TrackQuality trackQuality_
const reco::SimToRecoCollection * trackSimToRecoAssociation_
const edm::EDProductGetter * productGetter_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37

Member Data Documentation

const edm::EDProductGetter* VertexAssociatorByTracks::productGetter_
private

Definition at line 32 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByTracks::R2SMatchedRecoRatio_
private

Definition at line 35 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const double VertexAssociatorByTracks::R2SMatchedSimRatio_
private

Definition at line 34 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const double VertexAssociatorByTracks::S2RMatchedRecoRatio_
private

Definition at line 37 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().

const double VertexAssociatorByTracks::S2RMatchedSimRatio_
private

Definition at line 36 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().

const TrackingParticleSelector* VertexAssociatorByTracks::selector_
private

Definition at line 39 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const reco::TrackBase::TrackQuality VertexAssociatorByTracks::trackQuality_
private

Definition at line 40 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const reco::RecoToSimCollection* VertexAssociatorByTracks::trackRecoToSimAssociation_
private

Definition at line 42 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const reco::SimToRecoCollection* VertexAssociatorByTracks::trackSimToRecoAssociation_
private

Definition at line 43 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().