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 144 of file TrackQuality.cc.

144  :
145  associatorPSet_(config.getParameter<edm::ParameterSet>("hitAssociator"))
146 {
147 }
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 154 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, 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().

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

References associator_, and associatorPSet_.

Referenced by TrackClassifier::newEvent().

150 {
152 }
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().