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:
VertexAssociatorBase

Public Member Functions

reco::VertexRecoToSimCollection associateRecoToSim (edm::Handle< edm::View< reco::Vertex > > &, edm::Handle< TrackingVertexCollection > &, const edm::Event &, reco::RecoToSimCollection &) const
 
reco::VertexSimToRecoCollection associateSimToReco (edm::Handle< edm::View< reco::Vertex > > &, edm::Handle< TrackingVertexCollection > &, const edm::Event &, reco::SimToRecoCollection &) const
 
 VertexAssociatorByTracks (const edm::ParameterSet &)
 
 ~VertexAssociatorByTracks ()
 
- Public Member Functions inherited from VertexAssociatorBase
 VertexAssociatorBase ()
 
virtual ~VertexAssociatorBase ()
 

Private Attributes

const edm::ParameterSetconfig_
 
double R2SMatchedRecoRatio_
 
double R2SMatchedSimRatio_
 
double S2RMatchedRecoRatio_
 
double S2RMatchedSimRatio_
 
TrackingParticleSelector selector_
 
reco::TrackBase::TrackQuality trackQuality_
 

Detailed Description

Definition at line 11 of file VertexAssociatorByTracks.h.

Constructor & Destructor Documentation

VertexAssociatorByTracks::VertexAssociatorByTracks ( const edm::ParameterSet config)
explicit

Definition at line 24 of file VertexAssociatorByTracks.cc.

References edm::ParameterSet::getParameter(), reco::TrackBase::qualityByName(), R2SMatchedRecoRatio_, R2SMatchedSimRatio_, S2RMatchedRecoRatio_, S2RMatchedSimRatio_, selector_, AlCaHLTBitMon_QueryRunRegistry::string, and trackQuality_.

24  : config_(config)
25 {
26  R2SMatchedSimRatio_ = config.getParameter<double>("R2SMatchedSimRatio");
27  R2SMatchedRecoRatio_ = config.getParameter<double>("R2SMatchedRecoRatio");
28  S2RMatchedSimRatio_ = config.getParameter<double>("S2RMatchedSimRatio");
29  S2RMatchedRecoRatio_ = config.getParameter<double>("S2RMatchedRecoRatio");
30 
31  std::string trackQualityType = config.getParameter<std::string>("trackQuality");
33 
34  edm::ParameterSet param = config.getParameter<edm::ParameterSet>("trackingParticleSelector");
35 
37  param.getParameter<double>("ptMinTP"),
38  param.getParameter<double>("minRapidityTP"),
39  param.getParameter<double>("maxRapidityTP"),
40  param.getParameter<double>("tipTP"),
41  param.getParameter<double>("lipTP"),
42  param.getParameter<int>("minHitTP"),
43  param.getParameter<bool>("signalOnlyTP"),
44  param.getParameter<bool>("chargedOnlyTP"),
45  param.getParameter<bool>("stableOnlyTP"),
46  param.getParameter<std::vector<int> >("pdgIdTP")
47  );
48 }
T getParameter(std::string const &) const
reco::TrackBase::TrackQuality trackQuality_
TrackingParticleSelector selector_
SingleObjectSelector< TrackingParticleCollection,::TrackingParticleSelector > TrackingParticleSelector
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:114
const edm::ParameterSet & config_
VertexAssociatorByTracks::~VertexAssociatorByTracks ( )

Definition at line 52 of file VertexAssociatorByTracks.cc.

53 {
54  //do cleanup here
55 }

Member Function Documentation

reco::VertexRecoToSimCollection VertexAssociatorByTracks::associateRecoToSim ( edm::Handle< edm::View< reco::Vertex > > &  recoVertexes,
edm::Handle< TrackingVertexCollection > &  trackingVertexes,
const edm::Event event,
reco::RecoToSimCollection associator 
) const
virtual

Implements VertexAssociatorBase.

Definition at line 58 of file VertexAssociatorByTracks.cc.

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

64 {
65  reco::VertexRecoToSimCollection outputCollection;
66 
67  std::map<TrackingVertexRef,std::pair<double, std::size_t> > matches;
68 
69  std::cout << "reco::VertexCollection size = " << recoVertexes->size();
70  std::cout << " ; TrackingVertexCollection size = " << trackingVertexes->size() << std::endl << std::endl;
71 
72  // Loop over RecoVertex
73  for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex)
74  {
75  reco::VertexBaseRef recoVertex(recoVertexes, recoIndex);
76 
77  matches.clear();
78 
79  double recoDaughterWeight = 0.;
80 
81  // Loop over daughter tracks of RecoVertex
82  for (
83  reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
84  recoDaughter != recoVertex->tracks_end();
85  ++recoDaughter
86  )
87  {
88  // Check the quality of the RecoDoughters
89  if ( !(*recoDaughter)->quality(trackQuality_) ) continue;
90 
91  // Check for association for the given RecoDaughter
92  if ( associator.numberOfAssociations(*recoDaughter) > 0 )
93  {
94  std::vector<std::pair<TrackingParticleRef,double> > associations = associator[*recoDaughter];
95 
96  // Loop over TrackingParticles associated with RecoDaughter
97  for (
98  std::vector<std::pair<TrackingParticleRef,double> >::const_iterator association = associations.begin();
99  association != associations.end();
100  ++association
101  )
102  {
103  // Get a reference to parent vertex of TrackingParticle associated to the RecoDaughter
104  TrackingVertexRef trackingVertex = association->first->parentVertex();
105  // Store matched RecoDaughter to the trackingVertex
106  matches[trackingVertex].first += recoVertex->trackWeight(*recoDaughter);
107  matches[trackingVertex].second++;
108  }
109  }
110 
111  recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
112  } //loop over daughter tracks of RecoVertex
113 
114  std::size_t assoIndex = 0;
115 
116  // Loop over map between TrackingVertexes and matched RecoDaugther
117  for (
118  std::map<TrackingVertexRef,std::pair<double, std::size_t> >::const_iterator match = matches.begin();
119  match != matches.end();
120  ++match
121  )
122  {
123  // Getting the TrackingVertex information
124  TrackingVertexRef trackingVertex = match->first;
125  double matchedDaughterWeight = match->second.first;
126  std::size_t matchedDaughterCounter = match->second.second;
127 
128  // Count for only reconstructible SimDaughters
129  std::size_t simDaughterCounter = 0;
130 
131  for (
132  TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
133  simDaughter != trackingVertex->daughterTracks_end();
134  ++simDaughter
135  )
136  if ( selector_(**simDaughter) ) simDaughterCounter++;
137 
138  // Sanity condition in case that reconstructable condition is too tight
139  if ( simDaughterCounter < matchedDaughterCounter )
140  simDaughterCounter = matchedDaughterCounter;
141 
142  // Condition over S2RMatchedSimRatio
143  if ( (double)matchedDaughterCounter/simDaughterCounter < R2SMatchedSimRatio_ ) continue;
144 
145  double quality = (double)matchedDaughterWeight/recoDaughterWeight;
146 
147  // Condition over R2SMatchedRecoRatio
148  if (quality < R2SMatchedRecoRatio_) continue;
149 
150  outputCollection.insert(recoVertex, std::make_pair(trackingVertex, quality));
151 
152  std::cout << "R2S: INDEX " << assoIndex;
153  std::cout << " ; SimDaughterCounter = " << simDaughterCounter;
154  std::cout << " ; RecoDaughterWeight = " << recoDaughterWeight;
155  std::cout << " ; MatchedDaughterCounter = " << matchedDaughterCounter;
156  std::cout << " ; MatchedDaughterWeight = " << matchedDaughterWeight;
157  std::cout << " ; quality = " << quality << std::endl;
158 
159  assoIndex++;
160  }
161  } // Loop on RecoVertex
162 
163  std::cout << std::endl;
164  std::cout << "RecoToSim OUTPUT COLLECTION: outputCollection.size() = " << outputCollection.size() << std::endl << std::endl;
165 
166  return outputCollection;
167 }
reco::TrackBase::TrackQuality trackQuality_
TrackingParticleSelector selector_
size_type numberOfAssociations(const key_type &k) const
number of associations to a key
size_type size() const
map size
void insert(const key_type &k, const data_type &v)
insert an association
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
tuple cout
Definition: gather_cfg.py:121
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 ( edm::Handle< edm::View< reco::Vertex > > &  recoVertexes,
edm::Handle< TrackingVertexCollection > &  trackingVertexes,
const edm::Event event,
reco::SimToRecoCollection associator 
) const
virtual

Implements VertexAssociatorBase.

Definition at line 171 of file VertexAssociatorByTracks.cc.

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

177 {
178  reco::VertexSimToRecoCollection outputCollection;
179 
180  // Loop over TrackingVertexes
181  std::map<std::size_t,std::pair<double, std::size_t> > matches;
182 
183  // Loop over TrackingVertexes
184  for (std::size_t simIndex = 0; simIndex < trackingVertexes->size(); ++simIndex)
185  {
186  TrackingVertexRef trackingVertex(trackingVertexes, simIndex);
187 
188  matches.clear();
189 
190  std::size_t simDaughterCounter = 0;
191 
192  // Loop over daughter tracks of TrackingVertex
193  for (
194  TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
195  simDaughter != trackingVertex->daughterTracks_end();
196  ++simDaughter
197  )
198  {
199  // Select only reconstructible SimDaughters
200  if ( !selector_(**simDaughter) ) continue;
201 
202  // Check for association for the given RecoDaughter
203  if ( associator.numberOfAssociations(*simDaughter) > 0 )
204  {
205  std::vector<std::pair<reco::TrackBaseRef, double> > associations = associator[*simDaughter];
206 
207  // Loop over RecoTracks associated with TrackingParticle
208  for (
209  std::vector<std::pair<reco::TrackBaseRef,double> >::const_iterator association = associations.begin();
210  association != associations.end();
211  ++association
212  )
213  {
214  reco::TrackBaseRef recoTrack = association->first;
215 
216  for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex)
217  {
218  reco::VertexBaseRef recoVertex(recoVertexes, recoIndex);
219 
220  for (
221  reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
222  recoDaughter != recoVertex->tracks_end();
223  ++recoDaughter
224  )
225  {
226  // Store matched RecoDaughter to the RecoVertex
227  if (
228  recoDaughter->id() == recoTrack.id() &&
229  recoDaughter->key() == recoTrack.key()
230  )
231  {
232  matches[recoIndex].first += recoVertex->trackWeight(*recoDaughter);
233  matches[recoIndex].second++;
234  }
235  }
236  }
237  } // loop a recotracks
238  }
239 
240  simDaughterCounter++;
241  } // loop a simDaughter
242 
243  std::size_t assoIndex = 0;
244 
245  // Loop over map, set score, add to outputCollection
246  for (
247  std::map<std::size_t,std::pair<double,std::size_t> >::const_iterator match = matches.begin();
248  match != matches.end();
249  ++match
250  )
251  {
252  // Getting the TrackingVertex information
253  reco::VertexBaseRef recoVertex(recoVertexes, match->first);
254  double matchedDaughterWeight = match->second.first;
255  std::size_t matchedDaughterCounter = match->second.second;
256 
257  double recoDaughterWeight = 0.;
258 
259  // Weighted count of those tracks with a given quality of the RecoDoughters
260  for (
261  reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
262  recoDaughter != recoVertex->tracks_end();
263  ++recoDaughter
264  )
265  if ( (*recoDaughter)->quality(trackQuality_) )
266  recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
267 
268  // Sanity condition in case that track quality condition is too tight
269  if ( recoDaughterWeight < matchedDaughterWeight )
270  recoDaughterWeight = matchedDaughterWeight;
271 
272  // Condition over S2RMatchedRecoRatio
273  if ( matchedDaughterWeight/recoDaughterWeight < S2RMatchedRecoRatio_ ) continue;
274 
275  double quality = (double)matchedDaughterCounter/simDaughterCounter;
276 
277  // Condition over S2RMatchedSimRatio
278  if (quality < S2RMatchedSimRatio_) continue;
279 
280  outputCollection.insert(trackingVertex, std::make_pair(recoVertex, quality));
281 
282  std::cout << "R2S: INDEX " << assoIndex;
283  std::cout << " ; SimDaughterCounter = " << simDaughterCounter;
284  std::cout << " ; RecoDaughterWeight = " << recoDaughterWeight;
285  std::cout << " ; MatchedDaughterCounter = " << matchedDaughterCounter;
286  std::cout << " ; MatchedDaughterWeight = " << matchedDaughterWeight;
287  std::cout << " ; quality = " << quality << std::endl;
288 
289  assoIndex++;
290  }
291  } // loop over TrackingVertexes
292 
293  std::cout << "SimToReco OUTPUT COLLECTION: outputCollection.size() = " << outputCollection.size() << std::endl << std::endl;
294 
295  return outputCollection;
296 }
reco::TrackBase::TrackQuality trackQuality_
ProductID id() const
Definition: RefToBase.h:220
TrackingParticleSelector selector_
size_type numberOfAssociations(const key_type &k) const
number of associations to a key
size_t key() const
Definition: RefToBase.h:228
size_type size() const
map size
void insert(const key_type &k, const data_type &v)
insert an association
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
tuple cout
Definition: gather_cfg.py:121
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::ParameterSet& VertexAssociatorByTracks::config_
private

Definition at line 40 of file VertexAssociatorByTracks.h.

double VertexAssociatorByTracks::R2SMatchedRecoRatio_
private

Definition at line 43 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and VertexAssociatorByTracks().

double VertexAssociatorByTracks::R2SMatchedSimRatio_
private

Definition at line 42 of file VertexAssociatorByTracks.h.

Referenced by associateRecoToSim(), and VertexAssociatorByTracks().

double VertexAssociatorByTracks::S2RMatchedRecoRatio_
private

Definition at line 45 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco(), and VertexAssociatorByTracks().

double VertexAssociatorByTracks::S2RMatchedSimRatio_
private

Definition at line 44 of file VertexAssociatorByTracks.h.

Referenced by associateSimToReco(), and VertexAssociatorByTracks().

TrackingParticleSelector VertexAssociatorByTracks::selector_
private
reco::TrackBase::TrackQuality VertexAssociatorByTracks::trackQuality_
private