CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
VertexAssociatorByTracks Class Reference

#include <VertexAssociatorByTracks.h>

Inheritance diagram for VertexAssociatorByTracks:
reco::VertexToTrackingVertexAssociatorBaseImpl

Public Member Functions

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 More...
 
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 More...
 
 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 ()
 
- 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.

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

Definition at line 42 of file VertexAssociatorByTracks.cc.

42 {}

Member Function Documentation

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

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

Implements reco::VertexToTrackingVertexAssociatorBaseImpl.

Definition at line 44 of file VertexAssociatorByTracks.cc.

References edm::AssociationMap< Tag >::insert(), LogDebug, LogTrace, python.multivaluedict::map(), match(), matches, edm::AssociationMap< Tag >::numberOfAssociations(), productGetter_, HLT_25ns14e33_v1_cff::quality, R2SMatchedRecoRatio_, R2SMatchedSimRatio_, selector_, edm::AssociationMap< Tag >::size(), trackQuality_, trackRecoToSimAssociation_, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), and reco::Vertex::trackWeight().

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

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

Implements reco::VertexToTrackingVertexAssociatorBaseImpl.

Definition at line 154 of file VertexAssociatorByTracks.cc.

References edm::RefToBase< T >::id(), edm::AssociationMap< Tag >::insert(), edm::RefToBase< T >::key(), LogTrace, python.multivaluedict::map(), match(), matches, edm::AssociationMap< Tag >::numberOfAssociations(), productGetter_, HLT_25ns14e33_v1_cff::quality, S2RMatchedRecoRatio_, S2RMatchedSimRatio_, selector_, edm::AssociationMap< Tag >::size(), trackQuality_, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), trackSimToRecoAssociation_, and reco::Vertex::trackWeight().

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

Member Data Documentation

const edm::EDProductGetter* VertexAssociatorByTracks::productGetter_
private

Definition at line 34 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double VertexAssociatorByTracks::R2SMatchedRecoRatio_
private

Definition at line 37 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const double VertexAssociatorByTracks::R2SMatchedSimRatio_
private

Definition at line 36 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const double VertexAssociatorByTracks::S2RMatchedRecoRatio_
private

Definition at line 39 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().

const double VertexAssociatorByTracks::S2RMatchedSimRatio_
private

Definition at line 38 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().

const TrackingParticleSelector* VertexAssociatorByTracks::selector_
private

Definition at line 41 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

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

Definition at line 42 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const reco::RecoToSimCollection* VertexAssociatorByTracks::trackRecoToSimAssociation_
private

Definition at line 44 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim().

const reco::SimToRecoCollection* VertexAssociatorByTracks::trackSimToRecoAssociation_
private

Definition at line 45 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco().