Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:47

0001 #ifndef MuonAnalysis_MuonAssociators_MatcherUsingTracksAlgorithm_h
0002 #define MuonAnalysis_MuonAssociators_MatcherUsingTracksAlgorithm_h
0003 //
0004 //
0005 
0006 /**
0007   \class    pat::MatcherUsingTracksAlgorithm MatcherUsingTracksAlgorithm.h "MuonAnalysis/MuonAssociators/interface/MatcherUsingTracksAlgorithm.h"
0008   \brief    Matcher of reconstructed objects to other reconstructed objects using the tracks inside them 
0009             
0010   \author   Giovanni Petrucciani
0011   \version  $Id: MatcherUsingTracksAlgorithm.h,v 1.7 2010/07/12 20:56:10 gpetrucc Exp $
0012 */
0013 
0014 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0015 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0016 
0017 #include "FWCore/Framework/interface/ESHandle.h"
0018 #include "FWCore/Framework/interface/EventSetup.h"
0019 #include "FWCore/Framework/interface/ConsumesCollector.h"
0020 #include "FWCore/Utilities/interface/InputTag.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022 
0023 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
0024 #include "MagneticField/Engine/interface/MagneticField.h"
0025 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
0026 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0027 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToPoint.h"
0028 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0029 
0030 class IdealMagneticFieldRecord;
0031 class TrackingComponentsRecord;
0032 class GlobalTrackingGeometryRecord;
0033 
0034 class MatcherUsingTracksAlgorithm {
0035 public:
0036   explicit MatcherUsingTracksAlgorithm(const edm::ParameterSet &iConfig, edm::ConsumesCollector);
0037   virtual ~MatcherUsingTracksAlgorithm() {}
0038 
0039   /// Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators
0040   void init(const edm::EventSetup &iSetup);
0041 
0042   /// Try to match one track to another one. Return true if succeeded.
0043   /// For matches not by ref, it will update deltaR, deltaLocalPos and deltaPtRel if the match suceeded
0044   bool match(const reco::Candidate &c1,
0045              const reco::Candidate &c2,
0046              float &deltaR,
0047              float &deltaEta,
0048              float &deltaPhi,
0049              float &deltaLocalPos,
0050              float &deltaPtRel,
0051              float &chi2) const;
0052 
0053   /// Find the best match to another candidate, and return its index in the vector
0054   /// For matches not by ref, it will update deltaR, deltaLocalPos and deltaPtRel if the match suceeded
0055   /// Returns -1 if the match fails
0056   int match(const reco::Candidate &tk,
0057             const edm::View<reco::Candidate> &c2s,
0058             float &deltaR,
0059             float &deltaEta,
0060             float &deltaPhi,
0061             float &deltaLocalPos,
0062             float &deltaPtRel,
0063             float &chi2) const;
0064 
0065   /// Return 'true' if the matcher will produce meaningful deltaR, deltaLocalPos, deltaPtRel values
0066   bool hasMetrics() const { return algo_ != ByTrackRef; }
0067 
0068   /// Return 'true' if the matcher will produce also chi2
0069   bool hasChi2() const { return useChi2_; }
0070 
0071   /// Compute the chi2 of two free trajectory states, in the curvilinear frame (q/p, theta, phi, dxy, dsz)
0072   /// At least one must have errors
0073   /// diagonalOnly: don't use off-diagonal terms of covariance matrix
0074   /// useVertex   : use dxy, dsz in the chi2 (if false, use only q/p, theta, phi)
0075   /// useFirstMomentum : use the 'start' state momentum to compute dxy, dsx (if false, use 'other')
0076   static double getChi2(const FreeTrajectoryState &start,
0077                         const FreeTrajectoryState &other,
0078                         bool diagonalOnly,
0079                         bool useVertex,
0080                         bool useFirstMomentum);
0081 
0082   /// Compute the chi2 of one free trajectory state and a TrajectoryStateClosestToPoint closest to it, in the perigee frame
0083   /// At least one must have errors
0084   /// diagonalOnly: don't use off-diagonal terms of covariance matrix
0085   /// useVertex   : use dxy, dsz in the chi2 (if false, use only q/p, theta, phi)
0086   static double getChi2(const FreeTrajectoryState &start,
0087                         const TrajectoryStateClosestToPoint &other,
0088                         bool diagonalOnly,
0089                         bool useVertex);
0090 
0091   /// Compute the chi2 of two free trajectory states, in the local frame (q/p, dx, dy, dxdz, dydz)
0092   /// At least one must have errors
0093   /// diagonalOnly: don't use off-diagonal terms of covariance matrix
0094   /// useVertex   : use dx, dy in the chi2 (if false, use only direction and q/p)
0095   static double getChi2(const TrajectoryStateOnSurface &start,
0096                         const TrajectoryStateOnSurface &other,
0097                         bool diagonalOnly,
0098                         bool usePosition);
0099 
0100   /// Possibly crop the 3x3 part of the matrix or remove off-diagonal terms, then invert.
0101   static void cropAndInvert(AlgebraicSymMatrix55 &cov, bool diagonalOnly, bool top3by3only);
0102 
0103 private:
0104   enum AlgoType {
0105     ByTrackRef,
0106     ByDirectComparison,
0107     ByPropagatingSrc,
0108     ByPropagatingMatched,
0109     ByPropagatingSrcTSCP,
0110     ByPropagatingMatchedTSCP
0111   };  // propagate closest to point
0112   enum WhichTrack { None, TrackerTk, MuonTk, GlobalTk };
0113   enum WhichState { AtVertex, Innermost, Outermost };
0114 
0115   AlgoType algo_;
0116   WhichTrack whichTrack1_, whichTrack2_;
0117   WhichState whichState1_, whichState2_;
0118 
0119   // Preselection cuts on both sides
0120   StringCutObjectSelector<reco::Candidate, true> srcCut_, matchedCut_;
0121 
0122   // Matching cuts
0123   float maxLocalPosDiff_;
0124   float maxGlobalMomDeltaR_;
0125   float maxGlobalMomDeltaEta_;
0126   float maxGlobalMomDeltaPhi_;
0127   float maxGlobalDPtRel_;
0128   float maxChi2_;
0129   bool requireSameCharge_;
0130   bool useChi2_, chi2UseVertex_, chi2DiagonalOnly_, chi2FirstMomentum_;
0131   enum SortBy { LocalPosDiff, GlobalMomDeltaR, GlobalMomDeltaEta, GlobalMomDeltaPhi, GlobalDPtRel, Chi2 };
0132   SortBy sortBy_;
0133 
0134   //- Tools
0135   edm::ESHandle<MagneticField> magfield_;
0136   edm::ESHandle<Propagator> propagator_;
0137   edm::ESHandle<GlobalTrackingGeometry> geometry_;
0138 
0139   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magfieldToken_;
0140   edm::ESGetToken<Propagator, TrackingComponentsRecord> propagatorToken_;
0141   edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> geometryToken_;
0142 
0143   /// Get track reference out of a Candidate (via dynamic_cast to reco::RecoCandidate)
0144   reco::TrackRef getTrack(const reco::Candidate &reco, WhichTrack which) const;
0145 
0146   /// Starting state for the propagation
0147   FreeTrajectoryState startingState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const;
0148 
0149   /// End state for the propagation
0150   TrajectoryStateOnSurface targetState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const;
0151 
0152   /// Propagate and match. return true if current pair is the new best match (in that case, update also deltaR and deltaLocalPos)
0153   /// Uses TrajectoryStateClosestToPointBuilder
0154   bool matchWithPropagation(const FreeTrajectoryState &start,
0155                             const FreeTrajectoryState &target,
0156                             float &lastDeltaR,
0157                             float &lastDeltaEta,
0158                             float &lastDeltaPhi,
0159                             float &lastDeltaLocalPos,
0160                             float &lastGlobalDPtRel,
0161                             float &lastChi2) const;
0162 
0163   /// Propagate and match. return true if current pair is the new best match (in that case, update also deltaR and deltaLocalPos)
0164   /// Uses standard propagator to reach target's surface
0165   bool matchWithPropagation(const FreeTrajectoryState &start,
0166                             const TrajectoryStateOnSurface &target,
0167                             float &lastDeltaR,
0168                             float &lastDeltaEta,
0169                             float &lastDeltaPhi,
0170                             float &lastDeltaLocalPos,
0171                             float &lastGlobalDPtRel,
0172                             float &lastChi2) const;
0173 
0174   /// Compare directly two states. return true if current pair is the new best match (in that case, update also deltaR and deltaLocalPos)
0175   bool matchByDirectComparison(const FreeTrajectoryState &start,
0176                                const FreeTrajectoryState &other,
0177                                float &lastDeltaR,
0178                                float &lastDeltaEta,
0179                                float &lastDeltaPhi,
0180                                float &lastDeltaLocalPos,
0181                                float &lastGlobalDPtRel,
0182                                float &lastChi2) const;
0183 
0184   /// Parse some configuration
0185   void getConf(const edm::ParameterSet &iConfig,
0186                const std::string &whatFor,
0187                WhichTrack &whichTrack,
0188                WhichState &whichState);
0189 };
0190 
0191 #endif