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 
9 
10 #include <unordered_set>
11 
12 // Forward declarations
14 
15 namespace edm {
16  class EDProductGetter;
17 }
18 
67 {
68 public:
70 
72  std::unique_ptr<const TrackerHitAssociator> hitAssoc,
73  const ClusterTPAssociation *clusterToTPMap,
74  bool absoluteNumberOfHits,
75  double qualitySimToReco,
76  double puritySimToReco,
77  double cutRecoToSim,
78  bool threeHitTracksAreSpecial,
79  SimToRecoDenomType simToRecoDenominator);
80 
81  virtual
83  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
84  virtual
86  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
87  virtual
89  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
90 
91  virtual
93  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
94 
95  //seed
96  virtual
98  const edm::Handle<TrackingParticleCollection>&) const override;
99 
100  virtual
102  const edm::Handle<TrackingParticleCollection>&) const override;
103 
104 
105  private:
106  typedef std::pair<uint32_t,EncodedEventId> SimTrackIdentifiers;
107 
108  typedef std::unordered_set<reco::RecoToSimCollection::index_type> TrackingParticleRefKeySet;
109 
110  // - added by S. Sarkar
111  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()); }
112 
119  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
120  reco::RecoToSimCollection associateRecoToSimImplementation( const T_TrackCollection& trackCollection, const T_TrackingParticleCollection& trackingParticleCollection, const TrackingParticleRefKeySet *trackingParticleKeys, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
121 
128  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
129  reco::SimToRecoCollection associateSimToRecoImplementation( const T_TrackCollection& trackCollection, const T_TrackingParticleCollection& trackingParticleCollection, const TrackingParticleRefKeySet *trackingParticleKeys, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
130 
131 
137  template<typename T_TPCollection,typename iter> std::vector< std::pair<edm::Ref<TrackingParticleCollection>,size_t> > associateTrack( const TrackerHitAssociator& hitAssociator, const T_TPCollection& trackingParticles, const TrackingParticleRefKeySet *trackingParticleKeys, iter begin, iter end ) const;
145  template<typename T_TPCollection,typename iter> std::vector< std::pair<edm::Ref<TrackingParticleCollection>,size_t> > associateTrack( const ClusterTPAssociation& clusterToTPMap, const T_TPCollection& trackingParticles, const TrackingParticleRefKeySet *trackingParticleKeys, iter begin, iter end ) const;
146 
147 
149  bool trackingParticleContainsIdentifier( const TrackingParticle* pTrackingParticle, const SimTrackIdentifiers& identifier ) const;
150 
155  template<typename iter> int getDoubleCount( const TrackerHitAssociator& hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
158  template<typename iter> int getDoubleCount( const ClusterTPAssociation& clusterToTPList, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
159 
166  template<typename iter> std::vector< std::pair<SimTrackIdentifiers,size_t> > getAllSimTrackIdentifiers( const TrackerHitAssociator& hitAssociator, iter begin, iter end ) const;
167 
168  // Added by S. Sarkar
169  template<typename iter> std::vector< OmniClusterRef> getMatchedClusters( iter begin, iter end ) const;
170 
172  return &(**iter);
173  }
174 
175  const TrackingRecHit* getHitFromIter(TrackingRecHitCollection::const_iterator iter) const {
176  return &(*iter);
177  }
178 
192  //void prepareEitherHitAssociatorOrClusterToTPMap( const edm::Event* pEvent, std::unique_ptr<ClusterTPAssociation>& pClusterToTPMap, std::unique_ptr<TrackerHitAssociator>& pHitAssociator ) const;
193 
195  std::unique_ptr<const TrackerHitAssociator> hitAssociator_;
197 
204 
205  // Added by S. Sarkar
206  //bool useClusterTPAssociation_;
207 }; // end of the QuickTrackAssociatorByHitsImpl class
208 
209 #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)
reco::RecoToSimCollection associateRecoToSimImplementation(const T_TrackCollection &trackCollection, const T_TrackingParticleCollection &trackingParticleCollection, const TrackingParticleRefKeySet *trackingParticleKeys, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateRecoToSim.
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
edm::EDProductGetter const * productGetter_
creates either a ClusterTPAssociation OR a TrackerHitAssociator and stores it in the provided unique_...
const ClusterTPAssociation * clusterToTPMap_
EDProductGetter const * productGetter(std::atomic< void const * > const &iCache)
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::SimToRecoCollection associateSimToRecoImplementation(const T_TrackCollection &trackCollection, const T_TrackingParticleCollection &trackingParticleCollection, const TrackingParticleRefKeySet *trackingParticleKeys, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateSimToReco.
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 ...
std::vector< std::pair< edm::Ref< TrackingParticleCollection >, size_t > > associateTrack(const TrackerHitAssociator &hitAssociator, const T_TPCollection &trackingParticles, const TrackingParticleRefKeySet *trackingParticleKeys, iter begin, iter end) const
Returns the TrackingParticle that has the most associated hits to the given track.
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
tuple trackCollection
std::pair< uint32_t, EncodedEventId > SimTrackIdentifiers
This is enough information to uniquely identify a sim track.
int getDoubleCount(const TrackerHitAssociator &hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle) const
This method was copied almost verbatim from the standard TrackAssociatorByHits.
std::unordered_set< reco::RecoToSimCollection::index_type > TrackingParticleRefKeySet
Set for TrackingParticleRef keys.
#define begin
Definition: vmac.h:30
QuickTrackAssociatorByHitsImpl(edm::EDProductGetter const &productGetter, std::unique_ptr< const TrackerHitAssociator > hitAssoc, const ClusterTPAssociation *clusterToTPMap, bool absoluteNumberOfHits, double qualitySimToReco, double puritySimToReco, double cutRecoToSim, bool threeHitTracksAreSpecial, SimToRecoDenomType simToRecoDenominator)
Monte Carlo truth information used for tracking validation.
TrackToTrackingParticleAssociator that associates by hits a bit quicker than the normal TrackAssociat...
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
std::unique_ptr< const TrackerHitAssociator > hitAssociator_