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_, TrackQuality::Layer::Bad, TrackingRecHit::bad, gather_cfg::cout, TrackQuality::Layer::Dead, DetId::det(), getDetLayer(), 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::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  // in case we added missed modules, re-sort
232  std::stable_sort(matchedHits.begin(), matchedHits.end());
233 
234  // prepare for ordering results by layer enum and layer/disk number
235  typedef std::multimap<DetLayer, const MatchedHit*> LayerHitMap;
236  LayerHitMap layerHitMap;
237 
238  // iterate over all simulated/reconstructed hits again
239  for (std::vector<MatchedHit>::const_iterator hit = matchedHits.begin();
240  hit != matchedHits.end();)
241  {
242  // we can have multiple reco-to-sim matches per module, find best one
243  const MatchedHit *best = 0;
244 
245  // this loop iterates over all subsequent hits in the same module
246  do
247  {
248  // update our best hit pointer
249  if (!best ||
250  statePriorities[hit->state] > statePriorities[best->state] ||
251  best->simTrackId == NonMatchedTrackId)
252  {
253  best = &*hit;
254  }
255  ++hit;
256  }
257  while (hit != matchedHits.end() &&
258  hit->detId == best->detId);
259 
260  // ignore hit in case track reco was looking at the wrong module
261  if (best->simTrackId != NonMatchedTrackId ||
262  best->state != Layer::Missed)
263  {
264  layerHitMap.insert(std::make_pair(
265  getDetLayer(best->detId,tTopo), best));
266  }
267  }
268 
269  layers_.clear();
270 
271 #ifdef DEBUG_TRACK_QUALITY
272  std::cout << "---------------------" << std::endl;
273 #endif
274  // now prepare final collection
275  for (LayerHitMap::const_iterator hit = layerHitMap.begin();
276  hit != layerHitMap.end(); ++hit)
277  {
278 #ifdef DEBUG_TRACK_QUALITY
279  std::cout
280  << "detLayer (" << hit->first.first << ", " << hit->first.second << ")"
281  << " [" << (uint32_t)hit->second->detId << "] sim(" << (int)hit->second->simTrackId << ")"
282  << " hit(" << hit->second->recHitId << ") -> " << hit->second->state << std::endl;
283 #endif
284 
285  // find out if we need to start a new layer
286  Layer *layer = layers_.empty() ? 0 : &layers_.back();
287  if (!layer ||
288  hit->first.first != layer->subDet ||
289  hit->first.second != layer->layer)
290  {
291  Layer newLayer;
292  newLayer.subDet = hit->first.first;
293  newLayer.layer = hit->first.second;
294  layers_.push_back(newLayer);
295  layer = &layers_.back();
296  }
297 
298  // create hit and add it to layer
299  Layer::Hit newHit;
300  newHit.recHitId = hit->second->recHitId;
301  newHit.state = hit->second->state;
302 
303  layer->hits.push_back(newHit);
304  }
305 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< TrackerHitAssociator > associator_
Definition: TrackQuality.h:94
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
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
SeedingHitSet::ConstRecHitPointer Hit
Definition: DetId.h:18
tuple cout
Definition: gather_cfg.py:121
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
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 cmsHarvester::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().