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 // Forward declarations
12 
13 namespace edm {
14  class EDProductGetter;
15 }
16 
65 {
66 public:
68 
70  std::unique_ptr<const TrackerHitAssociator> hitAssoc,
71  const ClusterTPAssociationList *clusterToTPMap,
72  bool absoluteNumberOfHits,
73  double qualitySimToReco,
74  double puritySimToReco,
75  double cutRecoToSim,
76  bool threeHitTracksAreSpecial,
77  SimToRecoDenomType simToRecoDenominator);
78 
79  virtual
81  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
82  virtual
84  const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
85  virtual
87  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
88 
89  virtual
91  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
92 
93  //seed
94  virtual
96  const edm::Handle<TrackingParticleCollection>&) const override;
97 
98  virtual
100  const edm::Handle<TrackingParticleCollection>&) const override;
101 
102 
103  private:
104  typedef std::pair<uint32_t,EncodedEventId> SimTrackIdentifiers;
105 
106  // - added by S. Sarkar
107  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()); }
108 
115  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
116  reco::RecoToSimCollection associateRecoToSimImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
117 
124  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
125  reco::SimToRecoCollection associateSimToRecoImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
126 
127 
133  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;
141  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;
142 
143 
145  bool trackingParticleContainsIdentifier( const TrackingParticle* pTrackingParticle, const SimTrackIdentifiers& identifier ) const;
146 
151  template<typename iter> int getDoubleCount( const TrackerHitAssociator& hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
154  template<typename iter> int getDoubleCount( const ClusterTPAssociationList& clusterToTPList, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
155 
162  template<typename iter> std::vector< std::pair<SimTrackIdentifiers,size_t> > getAllSimTrackIdentifiers( const TrackerHitAssociator& hitAssociator, iter begin, iter end ) const;
163 
164  // Added by S. Sarkar
165  template<typename iter> std::vector< OmniClusterRef> getMatchedClusters( iter begin, iter end ) const;
166 
168  return &(**iter);
169  }
170 
171  const TrackingRecHit* getHitFromIter(TrackingRecHitCollection::const_iterator iter) const {
172  return &(*iter);
173  }
174 
188  //void prepareEitherHitAssociatorOrClusterToTPMap( const edm::Event* pEvent, std::unique_ptr<ClusterTPAssociationList>& pClusterToTPMap, std::unique_ptr<TrackerHitAssociator>& pHitAssociator ) const;
189 
191  std::unique_ptr<const TrackerHitAssociator> hitAssociator_;
193 
200 
201  // Added by S. Sarkar
202  //bool useClusterTPAssociation_;
203 }; // end of the QuickTrackAssociatorByHitsImpl class
204 
205 #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)
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 ClusterTPAssociationList OR a TrackerHitAssociator and stores it in the provided uni...
EDProductGetter const * productGetter(std::atomic< void const * > const &iCache)
QuickTrackAssociatorByHitsImpl(edm::EDProductGetter const &productGetter, std::unique_ptr< const TrackerHitAssociator > hitAssoc, const ClusterTPAssociationList *clusterToTPMap, bool absoluteNumberOfHits, double qualitySimToReco, double puritySimToReco, double cutRecoToSim, bool threeHitTracksAreSpecial, SimToRecoDenomType simToRecoDenominator)
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
std::vector< std::pair< OmniClusterRef, TrackingParticleRef > > ClusterTPAssociationList
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.
const ClusterTPAssociationList * clusterToTPMap_
#define begin
Definition: vmac.h:30
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.
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_