CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QuickTrackAssociatorByHitsImpl.h
Go to the documentation of this file.
1 #ifndef QuickTrackAssociatorByHitsImpl_h
2 #define QuickTrackAssociatorByHitsImpl_h
3 
5 
7 
8 // Forward declarations
10 
59  inline bool clusterTPAssociationListGreater(std::pair<OmniClusterRef, TrackingParticleRef> i,std::pair<OmniClusterRef, TrackingParticleRef> j) { return (i.first.rawIndex()>j.first.rawIndex()); }
60 
62 {
63 public:
64  typedef std::vector<std::pair<OmniClusterRef, TrackingParticleRef> > ClusterTPAssociationList;
66 
67  QuickTrackAssociatorByHitsImpl(std::shared_ptr<const TrackerHitAssociator> hitAssoc,
68  std::shared_ptr<const ClusterTPAssociationList> clusterToTPMap,
69  bool absoluteNumberOfHits,
70  double qualitySimToReco,
71  double puritySimToReco,
72  double cutRecoToSim,
73  bool threeHitTracksAreSpecial,
74  SimToRecoDenomType simToRecoDenominator);
75 
76  virtual
78  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
79  virtual
81  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
82  virtual
84  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
85 
86  virtual
88  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
89 
90  //seed
91  virtual
93  const edm::Handle<TrackingParticleCollection>&) const override;
94 
95  virtual
97  const edm::Handle<TrackingParticleCollection>&) const override;
98 
99 
100  private:
101  typedef std::pair<uint32_t,EncodedEventId> SimTrackIdentifiers;
102 
103  // - added by S. Sarkar
104  static bool tpIntPairGreater(std::pair<edm::Ref<TrackingParticleCollection>,size_t> i, std::pair<edm::Ref<TrackingParticleCollection>,size_t> j) { return (i.first.key()>j.first.key()); }
105 
112  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
113  reco::RecoToSimCollection associateRecoToSimImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
114 
121  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
122  reco::SimToRecoCollection associateSimToRecoImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
123 
124 
130  template<typename T_TPCollection,typename iter> std::vector< std::pair<edm::Ref<TrackingParticleCollection>,size_t> > associateTrack( const TrackerHitAssociator& hitAssociator, T_TPCollection trackingParticles, iter begin, iter end ) const;
138  template<typename T_TPCollection,typename iter> std::vector< std::pair<edm::Ref<TrackingParticleCollection>,size_t> > associateTrack( const ClusterTPAssociationList& clusterToTPMap, T_TPCollection trackingParticles, iter begin, iter end ) const;
139 
140 
142  bool trackingParticleContainsIdentifier( const TrackingParticle* pTrackingParticle, const SimTrackIdentifiers& identifier ) const;
143 
148  template<typename iter> int getDoubleCount( const TrackerHitAssociator& hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
151  template<typename iter> int getDoubleCount( const ClusterTPAssociationList& clusterToTPList, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
152 
159  template<typename iter> std::vector< std::pair<SimTrackIdentifiers,size_t> > getAllSimTrackIdentifiers( const TrackerHitAssociator& hitAssociator, iter begin, iter end ) const;
160 
161  // Added by S. Sarkar
162  template<typename iter> std::vector< OmniClusterRef> getMatchedClusters( iter begin, iter end ) const;
163 
165  return &(**iter);
166  }
167 
168  const TrackingRecHit* getHitFromIter(TrackingRecHitCollection::const_iterator iter) const {
169  return &(*iter);
170  }
171 
185  //void prepareEitherHitAssociatorOrClusterToTPMap( const edm::Event* pEvent, std::unique_ptr<ClusterTPAssociationList>& pClusterToTPMap, std::unique_ptr<TrackerHitAssociator>& pHitAssociator ) const;
186 
187  std::shared_ptr<const TrackerHitAssociator> hitAssociator_;
188  std::shared_ptr<const ClusterTPAssociationList> clusterToTPMap_;
189 
196 
197  // Added by S. Sarkar
198  //bool useClusterTPAssociation_;
199 }; // end of the QuickTrackAssociatorByHitsImpl class
200 
201 #endif // end of ifndef QuickTrackAssociatorByHitsImpl_h
int i
Definition: DBlmapReader.cc:9
static bool tpIntPairGreater(std::pair< edm::Ref< TrackingParticleCollection >, size_t > i, std::pair< edm::Ref< TrackingParticleCollection >, size_t > j)
std::vector< std::pair< OmniClusterRef, TrackingParticleRef > > ClusterTPAssociationList
std::shared_ptr< const ClusterTPAssociationList > clusterToTPMap_
virtual reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &trackCollectionHandle, const edm::Handle< TrackingParticleCollection > &trackingParticleCollectionHandle) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
virtual reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track > > &trackCollectionHandle, const edm::Handle< TrackingParticleCollection > &trackingParticleCollectionHandle) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
reco::RecoToSimCollection associateRecoToSimImplementation(T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateRecoToSim.
const TrackingRecHit * getHitFromIter(trackingRecHit_iterator iter) const
std::vector< std::pair< SimTrackIdentifiers, size_t > > getAllSimTrackIdentifiers(const TrackerHitAssociator &hitAssociator, iter begin, iter end) const
Returns a vector of pairs where first is a SimTrackIdentifiers (see typedef above) and second is the ...
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
bool clusterTPAssociationListGreater(std::pair< OmniClusterRef, TrackingParticleRef > i, std::pair< OmniClusterRef, TrackingParticleRef > j)
TrackToTrackingParticleAssociator that associates by hits a bit quicker than the normal TrackAssociat...
std::shared_ptr< const TrackerHitAssociator > hitAssociator_
creates either a ClusterTPAssociationList OR a TrackerHitAssociator and stores it in the provided uni...
std::pair< uint32_t, EncodedEventId > SimTrackIdentifiers
This is enough information to uniquely identify a sim track.
std::vector< std::pair< edm::Ref< TrackingParticleCollection >, size_t > > associateTrack(const TrackerHitAssociator &hitAssociator, T_TPCollection trackingParticles, iter begin, iter end) const
Returns the TrackingParticle that has the most associated hits to the given track.
int getDoubleCount(const TrackerHitAssociator &hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle) const
This method was copied almost verbatim from the standard TrackAssociatorByHits.
QuickTrackAssociatorByHitsImpl(std::shared_ptr< const TrackerHitAssociator > hitAssoc, std::shared_ptr< const ClusterTPAssociationList > clusterToTPMap, bool absoluteNumberOfHits, double qualitySimToReco, double puritySimToReco, double cutRecoToSim, bool threeHitTracksAreSpecial, SimToRecoDenomType simToRecoDenominator)
#define begin
Definition: vmac.h:30
QuickTrackAssociatorByHitsImpl::ClusterTPAssociationList ClusterTPAssociationList
reco::SimToRecoCollection associateSimToRecoImplementation(T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateSimToReco.
Monte Carlo truth information used for tracking validation.
const TrackingRecHit * getHitFromIter(TrackingRecHitCollection::const_iterator iter) const
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
bool trackingParticleContainsIdentifier(const TrackingParticle *pTrackingParticle, const SimTrackIdentifiers &identifier) const
Returns true if the supplied TrackingParticle has the supplied g4 track identifiers.
std::vector< OmniClusterRef > getMatchedClusters(iter begin, iter end) const