CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QuickTrackAssociatorByHits.h
Go to the documentation of this file.
1 #ifndef QuickTrackAssociatorByHits_h
2 #define QuickTrackAssociatorByHits_h
3 
7 
9 //#include "SimTracker/TrackerHitAssociation/plugins/Cluster2TPMapProducer.h"
10 //#include "SimTracker/TrackerHitAssociation/plugins/TP2ClusterMapProducer.h"
11 
12 // Forward declarations
14 
63 {
64 public:
67  virtual
69  edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle,
70  const edm::Event* pEvent=0,
71  const edm::EventSetup* pSetup=0 ) const override;
72  virtual
74  edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle,
75  const edm::Event* pEvent=0,
76  const edm::EventSetup* pSetup=0 ) const override;
77  virtual
79  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection,
80  const edm::Event* pEvent=0,
81  const edm::EventSetup* pSetup=0 ) const override;
82  virtual
84  const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection,
85  const edm::Event* pEvent=0,
86  const edm::EventSetup* pSetup=0 ) const override;
87 
88  //seed
89  virtual
92  const edm::Event * event ,
93  const edm::EventSetup * setup ) const override;
94 
95  virtual
98  const edm::Event * event ,
99  const edm::EventSetup * setup ) const override;
100 
101 
102 private:
103  typedef std::pair<uint32_t,EncodedEventId> SimTrackIdentifiers;
105 
106  // - added by S. Sarkar
107  typedef std::vector<std::pair<OmniClusterRef, TrackingParticleRef> > ClusterTPAssociationList;
108  static bool clusterTPAssociationListGreater(std::pair<OmniClusterRef, TrackingParticleRef> i,std::pair<OmniClusterRef, TrackingParticleRef> j) { return (i.first.rawIndex()>j.first.rawIndex()); }
109  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()); }
110 
117  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
118  reco::RecoToSimCollection associateRecoToSimImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
119 
126  template<class T_TrackCollection, class T_TrackingParticleCollection, class T_hitOrClusterAssociator>
127  reco::SimToRecoCollection associateSimToRecoImplementation( T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator ) const;
128 
129 
135  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;
143  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;
144 
145 
147  bool trackingParticleContainsIdentifier( const TrackingParticle* pTrackingParticle, const SimTrackIdentifiers& identifier ) const;
148 
153  template<typename iter> int getDoubleCount( const TrackerHitAssociator& hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
156  template<typename iter> int getDoubleCount( const ClusterTPAssociationList& clusterToTPList, iter begin, iter end, TrackingParticleRef associatedTrackingParticle ) const;
157 
164  template<typename iter> std::vector< std::pair<SimTrackIdentifiers,size_t> > getAllSimTrackIdentifiers( const TrackerHitAssociator& hitAssociator, iter begin, iter end ) const;
165 
166  // Added by S. Sarkar
167  template<typename iter> std::vector< OmniClusterRef> getMatchedClusters( iter begin, iter end ) const;
168 
170  return &(**iter);
171  }
172 
173  const TrackingRecHit* getHitFromIter(TrackingRecHitCollection::const_iterator iter) const {
174  return &(*iter);
175  }
176 
190  void prepareEitherHitAssociatorOrClusterToTPMap( const edm::Event* pEvent, std::unique_ptr<ClusterTPAssociationList>& pClusterToTPMap, std::unique_ptr<TrackerHitAssociator>& pHitAssociator ) const;
191 
193 
200 
201  // Added by S. Sarkar
204 }; // end of the QuickTrackAssociatorByHits class
205 
206 #endif // end of ifndef QuickTrackAssociatorByHits_h
reco::SimToRecoCollection associateSimToRecoImplementation(T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateSimToReco.
const TrackingRecHit * getHitFromIter(trackingRecHit_iterator iter) const
int i
Definition: DBlmapReader.cc:9
TrackAssociator that associates by hits a bit quicker than the normal TrackAssociatorByHits class...
static bool clusterTPAssociationListGreater(std::pair< OmniClusterRef, TrackingParticleRef > i, std::pair< OmniClusterRef, TrackingParticleRef > j)
const TrackingRecHit * getHitFromIter(TrackingRecHitCollection::const_iterator iter) const
bool trackingParticleContainsIdentifier(const TrackingParticle *pTrackingParticle, const SimTrackIdentifiers &identifier) const
Returns true if the supplied TrackingParticle has the supplied g4 track identifiers.
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &trackCollectionHandle, edm::Handle< TrackingParticleCollection > &trackingParticleCollectionHandle, const edm::Event *pEvent=0, const edm::EventSetup *pSetup=0) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
QuickTrackAssociatorByHits(const edm::ParameterSet &config)
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.
static bool tpIntPairGreater(std::pair< edm::Ref< TrackingParticleCollection >, size_t > i, std::pair< edm::Ref< TrackingParticleCollection >, size_t > j)
std::pair< uint32_t, EncodedEventId > SimTrackIdentifiers
This is enough information to uniquely identify a sim track.
std::vector< std::pair< OmniClusterRef, TrackingParticleRef > > ClusterTPAssociationList
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &trackCollectionHandle, edm::Handle< TrackingParticleCollection > &trackingParticleCollectionHandle, const edm::Event *pEvent=0, const edm::EventSetup *pSetup=0) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void prepareEitherHitAssociatorOrClusterToTPMap(const edm::Event *pEvent, std::unique_ptr< ClusterTPAssociationList > &pClusterToTPMap, std::unique_ptr< TrackerHitAssociator > &pHitAssociator) const
creates either a ClusterTPAssociationList OR a TrackerHitAssociator and stores it in the provided uni...
int getDoubleCount(const TrackerHitAssociator &hitAssociator, iter begin, iter end, TrackingParticleRef associatedTrackingParticle) const
This method was copied almost verbatim from the standard TrackAssociatorByHits.
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 ...
#define begin
Definition: vmac.h:30
QuickTrackAssociatorByHitsImpl::ClusterTPAssociationList ClusterTPAssociationList
std::vector< OmniClusterRef > getMatchedClusters(iter begin, iter end) const
reco::RecoToSimCollection associateRecoToSimImplementation(T_TrackCollection trackCollection, T_TrackingParticleCollection trackingParticleCollection, T_hitOrClusterAssociator hitOrClusterAssociator) const
The method that does the work for both overloads of associateRecoToSim.
Monte Carlo truth information used for tracking validation.
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection