Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:01

0001 #ifndef QuickTrackAssociatorByHitsImpl_h
0002 #define QuickTrackAssociatorByHitsImpl_h
0003 
0004 #include "DataFormats/TrackerRecHit2D/interface/OmniClusterRef.h"
0005 
0006 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
0007 
0008 #include "SimTracker/TrackerHitAssociation/interface/ClusterTPAssociation.h"
0009 
0010 #include "FWCore/Utilities/interface/IndexSet.h"
0011 
0012 // Forward declarations
0013 class TrackerHitAssociator;
0014 
0015 namespace edm {
0016   class EDProductGetter;
0017 }
0018 
0019 /** @brief TrackToTrackingParticleAssociator that associates by hits a bit quicker than the normal TrackAssociatorByHitsImpl class.
0020  *
0021  * NOTE - Doesn't implement the TrackCandidate association methods (from TrackAssociatorBase) so will always
0022  * return empty associations for those.
0023  *
0024  * This track associator (mostly) does the same as TrackAssociatorByHitsImpl, but faster. I've tested it a fair bit and can't find
0025  * any differences between the results of this and the standard TrackAssociatorByHitsImpl.
0026  *
0027  * Configuration parameters:
0028  *
0029  * AbsoluteNumberOfHits - bool - if true, Quality_SimToReco and Cut_RecoToSim are the absolute number of shared hits required for
0030  * association, not the percentage.
0031  *
0032  * Quality_SimToReco - double - The minimum amount of shared hits required, as a percentage of either the reconstructed hits or
0033  * simulated hits (see SimToRecoDenominator), for the track to be considered associated during a call to associateSimToReco. See
0034  * also AbsoluteNumberOfHits.
0035  *
0036  * Purity_SimToReco - double - The minimum amount of shared hits required, as a percentage of the reconstructed hits, for the
0037  * track to be considered associated during a call to associateSimToReco. Has no effect if AbsoluteNumberOfHits is true.
0038  *
0039  * Cut_RecoToSim - double - The minimum amount of shared hits required, as a percentage of the reconstructed hits, for the track
0040  * to be considered associated during a call to associateRecoToSim. See also AbsoluteNumberOfHits.
0041  *
0042  * ThreeHitTracksAreSpecial - bool - If true, tracks with 3 hits must have all their hits associated.
0043  *
0044  * SimToRecoDenominator - string - Must be either "sim" or "reco". If "sim" Quality_SimToReco is the percentage of simulated hits
0045  * that need to be shared. If "reco" then it's the percentage of reconstructed hits (i.e. same as Purity_SimToReco).
0046  *
0047  * associatePixel - bool - Passed on to the hit associator.
0048  *
0049  * associateStrip - bool - Passed on to the hit associator.
0050  *
0051  * requireStoredHits - bool - Whether or not to insist all TrackingParticles have at least one PSimHit. The PSimHits are not required
0052  * for the association, but the old TrackAssociatorByHitsImpl still had this requirement. Storing PSimHits in the TrackingParticle is now
0053  * optional (see TrackingTruthAccumulator which replaces TrackingTruthProducer). Having requireStoredHits set to true will mean no
0054  * TrackingParticles will be associated if you have chosen not to store the hits. The flag is only kept in order to retain the old
0055  * behaviour which can give very slightly different results.
0056  *
0057  * Note that the TrackAssociatorByHitsImpl parameters UseGrouped and UseSplitting are not used.
0058  *
0059  * @author Mark Grimes (mark.grimes@cern.ch)
0060  * @date 09/Nov/2010
0061  * Significant changes to remove any differences to the standard TrackAssociatorByHitsImpl results 07/Jul/2011.
0062  * Association for TrajectorySeeds added by Giuseppe Cerati sometime between 2011 and 2013.
0063  * Functionality to associate using pre calculated cluster to TrackingParticle maps added by Subir Sarker sometime in 2013.
0064  * Overhauled to remove mutables to make it thread safe by Mark Grimes 01/May/2014.
0065  */
0066 class QuickTrackAssociatorByHitsImpl : public reco::TrackToTrackingParticleAssociatorBaseImpl {
0067 public:
0068   enum SimToRecoDenomType { denomnone, denomsim, denomreco };
0069 
0070   QuickTrackAssociatorByHitsImpl(edm::EDProductGetter const& productGetter,
0071                                  std::unique_ptr<const TrackerHitAssociator> hitAssoc,
0072                                  const ClusterTPAssociation* clusterToTPMap,
0073                                  bool absoluteNumberOfHits,
0074                                  double qualitySimToReco,
0075                                  double puritySimToReco,
0076                                  double cutRecoToSim,
0077                                  double pixelHitWeight,
0078                                  bool threeHitTracksAreSpecial,
0079                                  SimToRecoDenomType simToRecoDenominator);
0080 
0081   reco::RecoToSimCollection associateRecoToSim(
0082       const edm::Handle<edm::View<reco::Track> >& trackCollectionHandle,
0083       const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
0084   reco::SimToRecoCollection associateSimToReco(
0085       const edm::Handle<edm::View<reco::Track> >& trackCollectionHandle,
0086       const edm::Handle<TrackingParticleCollection>& trackingParticleCollectionHandle) const override;
0087   reco::RecoToSimCollection associateRecoToSim(
0088       const edm::RefToBaseVector<reco::Track>& trackCollection,
0089       const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
0090 
0091   reco::SimToRecoCollection associateSimToReco(
0092       const edm::RefToBaseVector<reco::Track>& trackCollection,
0093       const edm::RefVector<TrackingParticleCollection>& trackingParticleCollection) const override;
0094 
0095   //seed
0096   reco::RecoToSimCollectionSeed associateRecoToSim(const edm::Handle<edm::View<TrajectorySeed> >&,
0097                                                    const edm::Handle<TrackingParticleCollection>&) const override;
0098 
0099   reco::SimToRecoCollectionSeed associateSimToReco(const edm::Handle<edm::View<TrajectorySeed> >&,
0100                                                    const edm::Handle<TrackingParticleCollection>&) const override;
0101 
0102   //candidate
0103   reco::RecoToSimCollectionTCandidate associateRecoToSim(
0104       const edm::Handle<TrackCandidateCollection>&, const edm::RefVector<TrackingParticleCollection>&) const override;
0105 
0106   reco::SimToRecoCollectionTCandidate associateSimToReco(
0107       const edm::Handle<TrackCandidateCollection>&, const edm::RefVector<TrackingParticleCollection>&) const override;
0108 
0109 private:
0110   typedef std::pair<uint32_t, EncodedEventId>
0111       SimTrackIdentifiers;  ///< @brief This is enough information to uniquely identify a sim track
0112 
0113   typedef edm::IndexSet TrackingParticleRefKeySet;  ///< @brief Set for TrackingParticleRef keys
0114 
0115   // - added by S. Sarkar
0116   static bool tpIntPairGreater(std::pair<edm::Ref<TrackingParticleCollection>, size_t> i,
0117                                std::pair<edm::Ref<TrackingParticleCollection>, size_t> j) {
0118     return (i.first.key() > j.first.key());
0119   }
0120 
0121   /** @brief The method that does the work for both overloads of associateRecoToSim.
0122    *
0123    * Parts that actually rely on the type of the collections are delegated out to overloaded functions
0124    * in the unnamed namespace of the .cc file. Parts that rely on the type of T_hitOrClusterAssociator
0125    * are delegated out to overloaded methods.
0126    */
0127   template <class T_RecoToSimCollection,
0128             class T_TrackCollection,
0129             class T_TrackingParticleCollection,
0130             class T_hitOrClusterAssociator>
0131   T_RecoToSimCollection associateRecoToSimImplementation(const T_TrackCollection& trackCollection,
0132                                                          const T_TrackingParticleCollection& trackingParticleCollection,
0133                                                          const TrackingParticleRefKeySet* trackingParticleKeys,
0134                                                          T_hitOrClusterAssociator hitOrClusterAssociator) const;
0135 
0136   /** @brief The method that does the work for both overloads of associateSimToReco.
0137    *
0138    * Parts that actually rely on the type of the collections are delegated out to overloaded functions
0139    * in the unnamed namespace of the .cc file. Parts that rely on the type of T_hitOrClusterAssociator
0140    * are delegated out to overloaded methods.
0141    */
0142   template <class T_SimToRecoCollection,
0143             class T_TrackCollection,
0144             class T_TrackingParticleCollection,
0145             class T_hitOrClusterAssociator>
0146   T_SimToRecoCollection associateSimToRecoImplementation(const T_TrackCollection& trackCollection,
0147                                                          const T_TrackingParticleCollection& trackingParticleCollection,
0148                                                          const TrackingParticleRefKeySet* trackingParticleKeys,
0149                                                          T_hitOrClusterAssociator hitOrClusterAssociator) const;
0150 
0151   /** @brief Returns the TrackingParticle that has the most associated hits to the given track.
0152    *
0153    * Return value is a vector of pairs, where first is an edm::Ref to the associated TrackingParticle, and second is
0154    * the number of associated hits.
0155    */
0156   template <typename T_TPCollection, typename iter>
0157   std::vector<std::pair<edm::Ref<TrackingParticleCollection>, double> > associateTrack(
0158       const TrackerHitAssociator& hitAssociator,
0159       const T_TPCollection& trackingParticles,
0160       const TrackingParticleRefKeySet* trackingParticleKeys,
0161       iter begin,
0162       iter end) const;
0163   /** @brief Returns the TrackingParticle that has the most associated hits to the given track.
0164    *
0165    * See the notes for the other overload for the return type.
0166    *
0167    * Note that the trackingParticles parameter is not actually required since all the information is in clusterToTPMap,
0168    * but the method signature has to match the other overload because it is called from a templated method.
0169    */
0170   template <typename T_TPCollection, typename iter>
0171   std::vector<std::pair<edm::Ref<TrackingParticleCollection>, double> > associateTrack(
0172       const ClusterTPAssociation& clusterToTPMap,
0173       const T_TPCollection& trackingParticles,
0174       const TrackingParticleRefKeySet* trackingParticleKeys,
0175       iter begin,
0176       iter end) const;
0177 
0178   /** @brief Returns true if the supplied TrackingParticle has the supplied g4 track identifiers. */
0179   bool trackingParticleContainsIdentifier(const TrackingParticle* pTrackingParticle,
0180                                           const SimTrackIdentifiers& identifier) const;
0181 
0182   /** @brief This method was copied almost verbatim from the standard TrackAssociatorByHits.
0183    *
0184    * Modified 01/May/2014 to take the TrackerHitAssociator as a parameter rather than using a member.
0185    */
0186   template <typename iter>
0187   double getDoubleCount(const TrackerHitAssociator& hitAssociator,
0188                         iter begin,
0189                         iter end,
0190                         TrackingParticleRef associatedTrackingParticle) const;
0191   /** @brief Overload for when using cluster to TrackingParticle association list.
0192    */
0193   template <typename iter>
0194   double getDoubleCount(const ClusterTPAssociation& clusterToTPList,
0195                         iter begin,
0196                         iter end,
0197                         TrackingParticleRef associatedTrackingParticle) const;
0198 
0199   /** @brief Returns a vector of pairs where first is a SimTrackIdentifiers (see typedef above) and second is the number of hits that came from that sim track.
0200    *
0201    * This is used so that the TrackingParticle collection only has to be looped over once to search for each sim track, rather than once per hit.
0202    * E.g. If all the hits in the reco track come from the same sim track, then there will only be one entry with second as the number of hits in
0203    * the track.
0204    */
0205   template <typename iter>
0206   std::vector<std::pair<SimTrackIdentifiers, double> > getAllSimTrackIdentifiers(
0207       const TrackerHitAssociator& hitAssociator, iter begin, iter end) const;
0208 
0209   // The last parameter is used to decide whether we cound hits or clusters
0210   template <typename T_Track>
0211   double weightedNumberOfTrackClusters(const T_Track& track, const TrackerHitAssociator&) const;
0212   template <typename T_Track>
0213   double weightedNumberOfTrackClusters(const T_Track& track, const ClusterTPAssociation&) const;
0214 
0215   // called only by weightedNumberOfTrackClusters(..., ClusterTPAssociation)
0216   template <typename iter>
0217   double weightedNumberOfTrackClusters(iter begin, iter end) const;
0218 
0219   /** @brief creates either a ClusterTPAssociation OR a TrackerHitAssociator and stores it in the provided unique_ptr. The other will be null.
0220    *
0221    * A decision is made whether to create a ClusterTPAssociation or a TrackerHitAssociator depending on how this
0222    * track associator was configured. If the ClusterTPAssociation couldn't be fetched from the event then it
0223    * falls back to creating a TrackerHitAssociator.
0224    *
0225    * Only one type will be created, never both. The other unique_ptr reference will be null so check for that
0226    * and decide which to use.
0227    *
0228    * N.B. The value of useClusterTPAssociation_ should not be used to decide which of the two pointers to use. If
0229    * the cluster to TrackingParticle couldn't be retrieved from the event then pClusterToTPMap will be null but
0230    * useClusterTPAssociation_ is no longer changed to false.
0231    */
0232   //void prepareEitherHitAssociatorOrClusterToTPMap( const edm::Event* pEvent, std::unique_ptr<ClusterTPAssociation>& pClusterToTPMap, std::unique_ptr<TrackerHitAssociator>& pHitAssociator ) const;
0233 
0234   edm::EDProductGetter const* productGetter_;
0235   std::unique_ptr<const TrackerHitAssociator> hitAssociator_;
0236   const ClusterTPAssociation* clusterToTPMap_;
0237 
0238   double qualitySimToReco_;
0239   double puritySimToReco_;
0240   double pixelHitWeight_;
0241   double cutRecoToSim_;
0242   SimToRecoDenomType simToRecoDenominator_;
0243   bool threeHitTracksAreSpecial_;
0244   bool absoluteNumberOfHits_;
0245 
0246   // Added by S. Sarkar
0247   //bool useClusterTPAssociation_;
0248 };  // end of the QuickTrackAssociatorByHitsImpl class
0249 
0250 #endif  // end of ifndef QuickTrackAssociatorByHitsImpl_h