CMS 3D CMS Logo

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

This class analyses the reconstruction quality for a given track. More...

#include <TrackQuality.h>

Classes

struct  Layer
 

Public Types

typedef std::vector
< TrackingParticleRef
SimParticleTrail
 

Public Member Functions

void evaluate (SimParticleTrail const &, reco::TrackBaseRef const &)
 Compute information about the track reconstruction quality. More...
 
const Layerlayer (unsigned int index) const
 Return information about the given layer by index. More...
 
void newEvent (const edm::Event &, const edm::EventSetup &)
 Pre-process event information (for accessing reconstruction information) More...
 
unsigned int numberOfLayers () const
 Return the number of layers with simulated and/or reconstructed hits. More...
 
 TrackQuality (const edm::ParameterSet &)
 Constructor by pset. More...
 

Private Attributes

std::auto_ptr
< TrackerHitAssociator
associator_
 
const edm::ParameterSet associatorPSet_
 
std::vector< Layerlayers_
 

Detailed Description

This class analyses the reconstruction quality for a given track.

Definition at line 26 of file TrackQuality.h.

Member Typedef Documentation

Definition at line 29 of file TrackQuality.h.

Constructor & Destructor Documentation

TrackQuality::TrackQuality ( const edm::ParameterSet config)

Constructor by pset.

Definition at line 182 of file TrackQuality.cc.

182  :
183  associatorPSet_(config.getParameter<edm::ParameterSet>("hitAssociator"))
184 {
185 }
T getParameter(std::string const &) const
const edm::ParameterSet associatorPSet_
Definition: TrackQuality.h:91

Member Function Documentation

void TrackQuality::evaluate ( SimParticleTrail const &  spt,
reco::TrackBaseRef const &  tr 
)

Compute information about the track reconstruction quality.

Definition at line 192 of file TrackQuality.cc.

References associator_, TrackingRecHit::bad, TrackQuality::Layer::Bad, gather_cfg::cout, TrackQuality::Layer::Dead, DetId::det(), getDetLayer(), TrackQuality::Layer::Good, TrackQuality::Layer::hits, i, TrackingRecHit::inactive, TrackQuality::Layer::layer, layer(), layers_, TrackQuality::Layer::Misassoc, TrackQuality::Layer::Missed, DetId::Muon, TrackQuality::Layer::Noise, pos, TrackQuality::Layer::Hit::recHitId, reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), TrackQuality::Layer::Shared, findQualityFiles::size, TrackQuality::Layer::Hit::state, TrackQuality::Layer::subDet, and TrackQuality::Layer::Unknown.

Referenced by TrackClassifier::qualityInformation().

194 {
195  std::vector<MatchedHit> matchedHits;
196 
197  // iterate over reconstructed hits
198  for (trackingRecHit_iterator hit = tr->recHitsBegin();
199  hit != tr->recHitsEnd(); ++hit)
200  {
201  // on which module the hit lies
202  DetId detId = (*hit)->geographicalId();
203 
204  // FIXME: check for double-sided modules?
205 
206  // didn't find a hit on that module
207  if (!(*hit)->isValid())
208  {
209  MatchedHit matchedHit;
210  matchedHit.detId = detId;
211  matchedHit.simTrackId = NonMatchedTrackId;
212  // check why hit wasn't valid and propagate information
213  switch ((*hit)->getType())
214  {
216  matchedHit.state = Layer::Dead;
217  break;
218 
219  case TrackingRecHit::bad:
220  matchedHit.state = Layer::Bad;
221  break;
222 
223  default:
224  matchedHit.state = Layer::Missed;
225  }
226  matchedHit.recHitId = hit - tr->recHitsBegin();
227  matchedHits.push_back(matchedHit);
228  continue;
229  }
230 
231  // find simulated tracks causing hit that was reconstructed
232  std::vector<SimHitIdpr> simIds = associator_->associateHitId(**hit);
233 
234  // must be noise or so
235  if (simIds.empty())
236  {
237  MatchedHit matchedHit;
238  matchedHit.detId = detId;
239  matchedHit.simTrackId = NonMatchedTrackId;
240  matchedHit.state = Layer::Noise;
241  matchedHit.recHitId = hit - tr->recHitsBegin();
242  matchedHits.push_back(matchedHit);
243  continue;
244  }
245 
246  // register all simulated tracks contributing
247  for (std::vector<SimHitIdpr>::const_iterator i = simIds.begin();
248  i != simIds.end(); ++i)
249  {
250  MatchedHit matchedHit;
251  matchedHit.detId = detId;
252  matchedHit.simTrackId = i->first;
253  matchedHit.collision = i->second;
254  // RecHit <-> SimHit matcher currently doesn't support muon system
255  if (detId.det() == DetId::Muon)
256  matchedHit.state = Layer::Unknown;
257  else
258  // assume hit was mismatched (until possible confirmation)
259  matchedHit.state = Layer::Misassoc;
260  matchedHit.recHitId = hit - tr->recHitsBegin();
261  matchedHits.push_back(matchedHit);
262  }
263  }
264 
265  // sort hits found so far by module id
266  std::stable_sort(matchedHits.begin(), matchedHits.end());
267 
268  std::vector<MatchedHit>::size_type size = matchedHits.size();
269 
270  // now iterate over simulated hits and compare (tracks in chain first)
271  for (SimParticleTrail::const_iterator track = spt.begin();
272  track != spt.end(); ++track)
273  {
274  // iterate over all hits in track
275  for (std::vector<PSimHit>::const_iterator hit =
276  (*track)->pSimHit_begin();
277  hit != (*track)->pSimHit_end(); ++hit)
278  {
279  MatchedHit matchedHit;
280  matchedHit.detId = DetId(hit->detUnitId());
281  matchedHit.simTrackId = hit->trackId();
282  matchedHit.collision = hit->eventId();
283 
284  // find range of reconstructed hits belonging to this module
285  std::pair<std::vector<MatchedHit>::iterator,
286  std::vector<MatchedHit>::iterator>
287  range = std::equal_range(
288  matchedHits.begin(),
289  matchedHits.begin() + size,
290  matchedHit.detId);
291 
292  // no reconstructed hit found, remember this as a missed module
293  if (range.first == range.second)
294  {
295  matchedHit.state = Layer::Missed;
296  matchedHit.recHitId = -1;
297  matchedHits.push_back(matchedHit);
298  continue;
299  }
300 
301  // now find if the hit belongs to the correct simulated track
302  std::vector<MatchedHit>::iterator pos =
303  std::lower_bound(range.first,
304  range.second,
305  matchedHit);
306 
307  // if it does, check for being a shared hit (was Misassoc before)
308  if (pos != range.second)
309  {
310  if (range.second - range.first > 1) // more than one SimHit
311  pos->state = Layer::Shared;
312  else
313  pos->state = Layer::Good; // only hit -> good hit
314  }
315  }
316  }
317 
318  // in case we added missed modules, re-sort
319  std::stable_sort(matchedHits.begin(), matchedHits.end());
320 
321  // prepare for ordering results by layer enum and layer/disk number
322  typedef std::multimap<DetLayer, const MatchedHit*> LayerHitMap;
323  LayerHitMap layerHitMap;
324 
325  // iterate over all simulated/reconstructed hits again
326  for (std::vector<MatchedHit>::const_iterator hit = matchedHits.begin();
327  hit != matchedHits.end();)
328  {
329  // we can have multiple reco-to-sim matches per module, find best one
330  const MatchedHit *best = 0;
331 
332  // this loop iterates over all subsequent hits in the same module
333  do
334  {
335  // update our best hit pointer
336  if (!best ||
337  statePriorities[hit->state] > statePriorities[best->state] ||
338  best->simTrackId == NonMatchedTrackId)
339  {
340  best = &*hit;
341  }
342  ++hit;
343  }
344  while (hit != matchedHits.end() &&
345  hit->detId == best->detId);
346 
347  // ignore hit in case track reco was looking at the wrong module
348  if (best->simTrackId != NonMatchedTrackId ||
349  best->state != Layer::Missed)
350  {
351  layerHitMap.insert(std::make_pair(
352  getDetLayer(best->detId), best));
353  }
354  }
355 
356  layers_.clear();
357 
358 #ifdef DEBUG_TRACK_QUALITY
359  std::cout << "---------------------" << std::endl;
360 #endif
361  // now prepare final collection
362  for (LayerHitMap::const_iterator hit = layerHitMap.begin();
363  hit != layerHitMap.end(); ++hit)
364  {
365 #ifdef DEBUG_TRACK_QUALITY
366  std::cout
367  << "detLayer (" << hit->first.first << ", " << hit->first.second << ")"
368  << " [" << (uint32_t)hit->second->detId << "] sim(" << (int)hit->second->simTrackId << ")"
369  << " hit(" << hit->second->recHitId << ") -> " << hit->second->state << std::endl;
370 #endif
371 
372  // find out if we need to start a new layer
373  Layer *layer = layers_.empty() ? 0 : &layers_.back();
374  if (!layer ||
375  hit->first.first != layer->subDet ||
376  hit->first.second != layer->layer)
377  {
378  Layer newLayer;
379  newLayer.subDet = hit->first.first;
380  newLayer.layer = hit->first.second;
381  layers_.push_back(newLayer);
382  layer = &layers_.back();
383  }
384 
385  // create hit and add it to layer
386  Layer::Hit newHit;
387  newHit.recHitId = hit->second->recHitId;
388  newHit.state = hit->second->state;
389 
390  layer->hits.push_back(newHit);
391  }
392 }
int i
Definition: DBlmapReader.cc:9
static DetLayer getDetLayer(DetId detId)
std::auto_ptr< TrackerHitAssociator > associator_
Definition: TrackQuality.h:92
uint16_t size_type
const Layer & layer(unsigned int index) const
Return information about the given layer by index.
Definition: TrackQuality.h:85
std::vector< Layer > layers_
Definition: TrackQuality.h:94
Definition: DetId.h:20
TransientTrackingRecHit::ConstRecHitPointer Hit
tuple cout
Definition: gather_cfg.py:121
Detector det() const
get the detector field from this detid
Definition: DetId.h:37
tuple size
Write out results.
static const int statePriorities[]
Definition: TrackQuality.cc:89
const Layer& TrackQuality::layer ( unsigned int  index) const
inline

Return information about the given layer by index.

Definition at line 85 of file TrackQuality.h.

References getHLTprescales::index, and layers_.

Referenced by evaluate(), geometryXMLparser.DTAlignable::index(), geometryXMLparser.CSCAlignable::index(), and TrackClassifier::qualityInformation().

86  {
87  return layers_[index];
88  }
std::vector< Layer > layers_
Definition: TrackQuality.h:94
void TrackQuality::newEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Pre-process event information (for accessing reconstruction information)

Definition at line 187 of file TrackQuality.cc.

References associator_, and associatorPSet_.

Referenced by TrackClassifier::newEvent().

188 {
190 }
std::auto_ptr< TrackerHitAssociator > associator_
Definition: TrackQuality.h:92
const edm::ParameterSet associatorPSet_
Definition: TrackQuality.h:91
unsigned int TrackQuality::numberOfLayers ( ) const
inline

Return the number of layers with simulated and/or reconstructed hits.

Definition at line 79 of file TrackQuality.h.

References layers_.

Referenced by TrackClassifier::qualityInformation().

80  {
81  return layers_.size();
82  }
std::vector< Layer > layers_
Definition: TrackQuality.h:94

Member Data Documentation

std::auto_ptr<TrackerHitAssociator> TrackQuality::associator_
private

Definition at line 92 of file TrackQuality.h.

Referenced by evaluate(), and newEvent().

const edm::ParameterSet TrackQuality::associatorPSet_
private

Definition at line 91 of file TrackQuality.h.

Referenced by newEvent().

std::vector<Layer> TrackQuality::layers_
private

Definition at line 94 of file TrackQuality.h.

Referenced by evaluate(), layer(), and numberOfLayers().