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 &, const TrackerTopology *tTopo)
 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 28 of file TrackQuality.h.

Member Typedef Documentation

Definition at line 31 of file TrackQuality.h.

Constructor & Destructor Documentation

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

Constructor by pset.

Definition at line 143 of file TrackQuality.cc.

143  :
144  associatorPSet_(config.getParameter<edm::ParameterSet>("hitAssociator"))
145 {
146 }
T getParameter(std::string const &) const
const edm::ParameterSet associatorPSet_
Definition: TrackQuality.h:93

Member Function Documentation

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

Compute information about the track reconstruction quality.

Definition at line 153 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().

156 {
157  std::vector<MatchedHit> matchedHits;
158 
159  // iterate over reconstructed hits
160  for (trackingRecHit_iterator hit = tr->recHitsBegin();
161  hit != tr->recHitsEnd(); ++hit)
162  {
163  // on which module the hit lies
164  DetId detId = (*hit)->geographicalId();
165 
166  // FIXME: check for double-sided modules?
167 
168  // didn't find a hit on that module
169  if (!(*hit)->isValid())
170  {
171  MatchedHit matchedHit;
172  matchedHit.detId = detId;
173  matchedHit.simTrackId = NonMatchedTrackId;
174  // check why hit wasn't valid and propagate information
175  switch ((*hit)->getType())
176  {
178  matchedHit.state = Layer::Dead;
179  break;
180 
181  case TrackingRecHit::bad:
182  matchedHit.state = Layer::Bad;
183  break;
184 
185  default:
186  matchedHit.state = Layer::Missed;
187  }
188  matchedHit.recHitId = hit - tr->recHitsBegin();
189  matchedHits.push_back(matchedHit);
190  continue;
191  }
192 
193  // find simulated tracks causing hit that was reconstructed
194  std::vector<SimHitIdpr> simIds = associator_->associateHitId(**hit);
195 
196  // must be noise or so
197  if (simIds.empty())
198  {
199  MatchedHit matchedHit;
200  matchedHit.detId = detId;
201  matchedHit.simTrackId = NonMatchedTrackId;
202  matchedHit.state = Layer::Noise;
203  matchedHit.recHitId = hit - tr->recHitsBegin();
204  matchedHits.push_back(matchedHit);
205  continue;
206  }
207 
208  // register all simulated tracks contributing
209  for (std::vector<SimHitIdpr>::const_iterator i = simIds.begin();
210  i != simIds.end(); ++i)
211  {
212  MatchedHit matchedHit;
213  matchedHit.detId = detId;
214  matchedHit.simTrackId = i->first;
215  matchedHit.collision = i->second;
216  // RecHit <-> SimHit matcher currently doesn't support muon system
217  if (detId.det() == DetId::Muon)
218  matchedHit.state = Layer::Unknown;
219  else
220  // assume hit was mismatched (until possible confirmation)
221  matchedHit.state = Layer::Misassoc;
222  matchedHit.recHitId = hit - tr->recHitsBegin();
223  matchedHits.push_back(matchedHit);
224  }
225  }
226 
227  // sort hits found so far by module id
228  std::stable_sort(matchedHits.begin(), matchedHits.end());
229 
230  std::vector<MatchedHit>::size_type size = matchedHits.size();
231 
232 #warning "This file has been modified just to get it to compile without any regard as to whether it still functions as intended"
233 #ifdef REMOVED_JUST_TO_GET_IT_TO_COMPILE__THIS_CODE_NEEDS_TO_BE_CHECKED
234  // now iterate over simulated hits and compare (tracks in chain first)
235  for (SimParticleTrail::const_iterator track = spt.begin();
236  track != spt.end(); ++track)
237  {
238  // iterate over all hits in track
239  for (std::vector<PSimHit>::const_iterator hit =
240  (*track)->pSimHit_begin();
241  hit != (*track)->pSimHit_end(); ++hit)
242  {
243  MatchedHit matchedHit;
244  matchedHit.detId = DetId(hit->detUnitId());
245  matchedHit.simTrackId = hit->trackId();
246  matchedHit.collision = hit->eventId();
247 
248  // find range of reconstructed hits belonging to this module
249  std::pair<std::vector<MatchedHit>::iterator,
250  std::vector<MatchedHit>::iterator>
251  range = std::equal_range(
252  matchedHits.begin(),
253  matchedHits.begin() + size,
254  matchedHit.detId);
255 
256  // no reconstructed hit found, remember this as a missed module
257  if (range.first == range.second)
258  {
259  matchedHit.state = Layer::Missed;
260  matchedHit.recHitId = -1;
261  matchedHits.push_back(matchedHit);
262  continue;
263  }
264 
265  // now find if the hit belongs to the correct simulated track
266  std::vector<MatchedHit>::iterator pos =
267  std::lower_bound(range.first,
268  range.second,
269  matchedHit);
270 
271  // if it does, check for being a shared hit (was Misassoc before)
272  if (pos != range.second)
273  {
274  if (range.second - range.first > 1) // more than one SimHit
275  pos->state = Layer::Shared;
276  else
277  pos->state = Layer::Good; // only hit -> good hit
278  }
279  }
280  }
281 #endif
282 
283  // in case we added missed modules, re-sort
284  std::stable_sort(matchedHits.begin(), matchedHits.end());
285 
286  // prepare for ordering results by layer enum and layer/disk number
287  typedef std::multimap<DetLayer, const MatchedHit*> LayerHitMap;
288  LayerHitMap layerHitMap;
289 
290  // iterate over all simulated/reconstructed hits again
291  for (std::vector<MatchedHit>::const_iterator hit = matchedHits.begin();
292  hit != matchedHits.end();)
293  {
294  // we can have multiple reco-to-sim matches per module, find best one
295  const MatchedHit *best = 0;
296 
297  // this loop iterates over all subsequent hits in the same module
298  do
299  {
300  // update our best hit pointer
301  if (!best ||
302  statePriorities[hit->state] > statePriorities[best->state] ||
303  best->simTrackId == NonMatchedTrackId)
304  {
305  best = &*hit;
306  }
307  ++hit;
308  }
309  while (hit != matchedHits.end() &&
310  hit->detId == best->detId);
311 
312  // ignore hit in case track reco was looking at the wrong module
313  if (best->simTrackId != NonMatchedTrackId ||
314  best->state != Layer::Missed)
315  {
316  layerHitMap.insert(std::make_pair(
317  getDetLayer(best->detId,tTopo), best));
318  }
319  }
320 
321  layers_.clear();
322 
323 #ifdef DEBUG_TRACK_QUALITY
324  std::cout << "---------------------" << std::endl;
325 #endif
326  // now prepare final collection
327  for (LayerHitMap::const_iterator hit = layerHitMap.begin();
328  hit != layerHitMap.end(); ++hit)
329  {
330 #ifdef DEBUG_TRACK_QUALITY
331  std::cout
332  << "detLayer (" << hit->first.first << ", " << hit->first.second << ")"
333  << " [" << (uint32_t)hit->second->detId << "] sim(" << (int)hit->second->simTrackId << ")"
334  << " hit(" << hit->second->recHitId << ") -> " << hit->second->state << std::endl;
335 #endif
336 
337  // find out if we need to start a new layer
338  Layer *layer = layers_.empty() ? 0 : &layers_.back();
339  if (!layer ||
340  hit->first.first != layer->subDet ||
341  hit->first.second != layer->layer)
342  {
343  Layer newLayer;
344  newLayer.subDet = hit->first.first;
345  newLayer.layer = hit->first.second;
346  layers_.push_back(newLayer);
347  layer = &layers_.back();
348  }
349 
350  // create hit and add it to layer
351  Layer::Hit newHit;
352  newHit.recHitId = hit->second->recHitId;
353  newHit.state = hit->second->state;
354 
355  layer->hits.push_back(newHit);
356  }
357 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< TrackerHitAssociator > associator_
Definition: TrackQuality.h:94
uint16_t size_type
DetLayer getDetLayer(DetId detId, const TrackerTopology *tTopo)
Definition: TrackQuality.cc:97
const Layer & layer(unsigned int index) const
Return information about the given layer by index.
Definition: TrackQuality.h:87
std::vector< Layer > layers_
Definition: TrackQuality.h:96
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:85
const Layer& TrackQuality::layer ( unsigned int  index) const
inline

Return information about the given layer by index.

Definition at line 87 of file TrackQuality.h.

References getHLTprescales::index, and layers_.

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

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

Pre-process event information (for accessing reconstruction information)

Definition at line 148 of file TrackQuality.cc.

References associator_, and associatorPSet_.

Referenced by TrackClassifier::newEvent().

149 {
151 }
std::auto_ptr< TrackerHitAssociator > associator_
Definition: TrackQuality.h:94
const edm::ParameterSet associatorPSet_
Definition: TrackQuality.h:93
unsigned int TrackQuality::numberOfLayers ( ) const
inline

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

Definition at line 81 of file TrackQuality.h.

References layers_.

Referenced by TrackClassifier::qualityInformation().

82  {
83  return layers_.size();
84  }
std::vector< Layer > layers_
Definition: TrackQuality.h:96

Member Data Documentation

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

Definition at line 94 of file TrackQuality.h.

Referenced by evaluate(), and newEvent().

const edm::ParameterSet TrackQuality::associatorPSet_
private

Definition at line 93 of file TrackQuality.h.

Referenced by newEvent().

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

Definition at line 96 of file TrackQuality.h.

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