CMS 3D CMS Logo

TrackAssociatorByPositionImpl.h
Go to the documentation of this file.
1 #ifndef TrackAssociatorByPositionImpl_h
2 #define TrackAssociatorByPositionImpl_h
3 
17 
20 
22 
24 
25 #include <map>
26 
27 //Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater:
28 //the track with the lowest association chi2 will be the first in the output map.
29 
31 public:
32  typedef std::pair<TrackingParticleRef, TrackPSimHitRef> SimHitTPPair;
33  typedef std::vector<SimHitTPPair> SimHitTPAssociationList;
34  enum class Method { chi2, dist, momdr, posdr };
35 
37  const Propagator* prop,
38  const SimHitTPAssociationList* assocList,
39  double qMinCut,
40  double qCut,
42  Method method,
43  bool minIfNoMatch,
44  bool considerAllSimHits)
45  : theGeometry(geo),
46  thePropagator(prop),
47  theSimHitsTPAssoc(assocList),
48  theQminCut(qMinCut),
49  theQCut(qCut),
50  thePositionMinimumDistance(positionMinimumDistance),
51  theMethod(method),
52  theMinIfNoMatch(minIfNoMatch),
53  theConsiderAllSimHits(considerAllSimHits) {}
54 
56 
58  const edm::RefVector<TrackingParticleCollection>&) const override;
59 
61 
63  const edm::RefVector<TrackingParticleCollection>&) const override;
64 
65 private:
66  double quality(const TrajectoryStateOnSurface&, const TrajectoryStateOnSurface&) const;
67 
70  const SimHitTPAssociationList* theSimHitsTPAssoc;
71  double theQminCut;
72  double theQCut;
77 
79  TrajectoryStateOnSurface getState(const TrackingParticleRef&, const SimHitTPAssociationList& simHitsTPAssoc) const;
80  //edm::InputTag _simHitTpMapTag;
81 };
82 
83 #endif
double quality(const TrajectoryStateOnSurface &, const TrajectoryStateOnSurface &) const
FreeTrajectoryState getState(const reco::Track &) const
const SimHitTPAssociationList * theSimHitsTPAssoc
reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
std::pair< TrackingParticleRef, TrackPSimHitRef > SimHitTPPair
TrackAssociatorByPositionImpl(const TrackingGeometry *geo, const Propagator *prop, const SimHitTPAssociationList *assocList, double qMinCut, double qCut, double positionMinimumDistance, Method method, bool minIfNoMatch, bool considerAllSimHits)
reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
std::vector< SimHitTPPair > SimHitTPAssociationList