Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:48

0001 #ifndef L1Trigger_TrackFindingTMTT_KFTrackletTrack_h
0002 #define L1Trigger_TrackFindingTMTT_KFTrackletTrack_h
0003 
0004 #include "FWCore/Utilities/interface/Exception.h"
0005 #include "DataFormats/Math/interface/deltaPhi.h"
0006 #include "L1Trigger/TrackFindingTMTT/interface/L1trackBase.h"
0007 #include "L1Trigger/TrackFindingTMTT/interface/L1track3D.h"
0008 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0009 #include "L1Trigger/TrackFindingTMTT/interface/TP.h"
0010 #include "L1Trigger/TrackFindingTMTT/interface/Stub.h"
0011 #include "L1Trigger/TrackFindingTMTT/interface/DigitalTrack.h"
0012 
0013 #include <vector>
0014 #include <utility>
0015 #include <string>
0016 
0017 //=== This is used uniquely for HYBRID TRACKING.
0018 //=== It is the equivalent of class L1fittedTrack.
0019 //===
0020 //=== This represents a fitted L1 track candidate found in 3 dimensions.
0021 //=== It gives access to the fitted helix parameters & chi2 etc.
0022 //=== It also calculates & gives access to associated truth particle (Tracking Particle) if any.
0023 //=== It also gives access to the 3D hough-transform track candidate (L1track3D) on which the fit was run.
0024 
0025 namespace tmtt {
0026 
0027   class KFTrackletTrack {
0028   public:
0029     // Store a new fitted track, specifying the input Hough transform track, the stubs used for the fit,
0030     // bit-encoded hit layers,
0031     // the fitted helix parameters & chi2,
0032     // and the number of helix parameters being fitted (=5 if d0 is fitted, or =4 if d0 is not fitted).
0033     // Also specify phi sector and eta region used by track-finding code that this track was in.
0034     // And if track fit declared this to be a valid track (enough stubs left on track after fit etc.).
0035     KFTrackletTrack(const L1track3D* l1track3D,
0036                     const std::vector<const Stub*>& stubs,
0037                     unsigned int hitPattern,
0038                     float qOverPt,
0039                     float d0,
0040                     float phi0,
0041                     float z0,
0042                     float tanLambda,
0043                     float chi2rphi,
0044                     float chi2rz,
0045                     unsigned int nHelixParam,
0046                     unsigned int iPhiSec,
0047                     unsigned int iEtaReg,
0048                     bool accepted = true,
0049                     bool done_bcon = false,
0050                     float qOverPt_bcon = 0.,
0051                     float d0_bcon = 0.,
0052                     float phi0_bcon = 0.,
0053                     float chi2rphi_bcon = 0.)
0054         : l1track3D_(l1track3D),
0055           stubs_(stubs),
0056           hitPattern_(hitPattern),
0057           qOverPt_(qOverPt),
0058           d0_(d0),
0059           phi0_(phi0),
0060           z0_(z0),
0061           tanLambda_(tanLambda),
0062           chi2rphi_(chi2rphi),
0063           chi2rz_(chi2rz),
0064           done_bcon_(done_bcon),
0065           qOverPt_bcon_(qOverPt_bcon),
0066           d0_bcon_(d0_bcon),
0067           phi0_bcon_(phi0_bcon),
0068           chi2rphi_bcon_(chi2rphi_bcon),
0069           nHelixParam_(nHelixParam),
0070           iPhiSec_(iPhiSec),
0071           iEtaReg_(iEtaReg),
0072           optoLinkID_(l1track3D->optoLinkID()),
0073           nSkippedLayers_(0),
0074           numUpdateCalls_(0),
0075           numIterations_(0),
0076           accepted_(accepted) {}
0077 
0078     //--- Set/get additional info about fitted track that is specific to individual track fit algorithms (KF, LR, chi2)
0079     //--- and is used for debugging/histogramming purposes.
0080 
0081     void setInfoKF(unsigned int nSkippedLayers, unsigned int numUpdateCalls) {
0082       nSkippedLayers_ = nSkippedLayers;
0083       numUpdateCalls_ = numUpdateCalls;
0084     }
0085 
0086     void infoKF(unsigned int& nSkippedLayers, unsigned int& numUpdateCalls) const {
0087       nSkippedLayers = nSkippedLayers_;
0088       numUpdateCalls = numUpdateCalls_;
0089     }
0090 
0091     const L1track3D* l1track3D() const { return l1track3D_; }
0092 
0093     // Get stubs on fitted track (can differ from those on HT track if track fit kicked out stubs with bad residuals)
0094     const std::vector<const Stub*>& stubs() const { return stubs_; }
0095     // Get number of stubs on fitted track.
0096     unsigned int numStubs() const { return stubs_.size(); }
0097     // Get number of tracker layers these stubs are in.
0098     unsigned int numLayers() const { return nLayers_; }
0099     // Get number of stubs deleted from track candidate by fitter (because they had large residuals)
0100     unsigned int numKilledStubs() const { return l1track3D_->numStubs() - this->numStubs(); }
0101 
0102     // Get bit-encoded hit pattern (where layer number assigned by increasing distance from origin, according to layers track expected to cross).
0103     unsigned int hitPattern() const { return hitPattern_; }
0104 
0105     //--- Get the fitted track helix parameters.
0106 
0107     float qOverPt() const { return qOverPt_; }
0108     float charge() const { return (qOverPt_ > 0 ? 1 : -1); }
0109     float invPt() const { return std::abs(qOverPt_); }
0110     // Protect pt against 1/pt = 0.
0111     float pt() const {
0112       constexpr float small = 1.0e-6;
0113       return 1. / (small + this->invPt());
0114     }
0115     float d0() const { return d0_; }
0116     float phi0() const { return phi0_; }
0117     float z0() const { return z0_; }
0118     float tanLambda() const { return tanLambda_; }
0119     float theta() const { return atan2(1., tanLambda_); }  // Use atan2 to ensure 0 < theta < pi.
0120     float eta() const { return -log(tan(0.5 * this->theta())); }
0121 
0122     //--- Get the fitted helix parameters with beam-spot constraint.
0123     //--- If constraint not applied (e.g. 4 param fit) then these are identical to unconstrained values.
0124 
0125     bool done_bcon() const { return done_bcon_; }  // Was beam-spot constraint aplied?
0126     float qOverPt_bcon() const { return qOverPt_bcon_; }
0127     float charge_bcon() const { return (qOverPt_bcon_ > 0 ? 1 : -1); }
0128     float invPt_bcon() const { return std::abs(qOverPt_bcon_); }
0129     float pt_bcon() const { return 1. / (1.0e-6 + this->invPt_bcon()); }
0130     float phi0_bcon() const { return phi0_bcon_; }
0131     float d0_bcon() const { return d0_bcon_; }
0132 
0133     // Phi and z coordinates at which track crosses "chosenR" values used by r-phi HT and rapidity sectors respectively.
0134     // (Optionally with beam-spot constraint applied).
0135     float phiAtChosenR(bool beamConstraint) const {
0136       if (beamConstraint) {
0137         return reco::deltaPhi(phi0_bcon_ -
0138                                   asin((settings_->invPtToDphi() * settings_->chosenRofPhi()) * qOverPt_bcon_) -
0139                                   d0_bcon_ / (settings_->chosenRofPhi()),
0140                               0.);
0141       } else {
0142         return reco::deltaPhi(phi0_ - asin((settings_->invPtToDphi() * settings_->chosenRofPhi()) * qOverPt_) -
0143                                   d0_ / (settings_->chosenRofPhi()),
0144                               0.);
0145       }
0146     }
0147     float zAtChosenR() const {
0148       return (z0_ + (settings_->chosenRofZ()) * tanLambda_);
0149     }  // neglects transverse impact parameter & track curvature.
0150 
0151     // Get the number of helix parameters being fitted (=5 if d0 is fitted or =4 if d0 is not fitted).
0152     float nHelixParam() const { return nHelixParam_; }
0153 
0154     // Get the fit degrees of freedom, chi2 & chi2/DOF
0155     unsigned int numDOF() const { return 2 * this->numStubs() - nHelixParam_; }
0156     unsigned int numDOFrphi() const { return this->numStubs() - (nHelixParam_ - 2); }
0157     unsigned int numDOFrz() const { return this->numStubs() - 2; }
0158     float chi2rphi() const { return chi2rphi_; }
0159     float chi2rz() const { return chi2rz_; }
0160     float chi2() const { return chi2rphi_ + chi2rz_; }
0161     float chi2dof() const { return (this->chi2()) / this->numDOF(); }
0162 
0163     //--- Ditto, but if beam-spot constraint is applied.
0164     //--- If constraint not applied (e.g. 4 param fit) then these are identical to unconstrained values.
0165     unsigned int numDOF_bcon() const { return (this->numDOF() - 1); }
0166     unsigned int numDOFrphi_bcon() const { return (this->numDOFrphi() - 1); }
0167     float chi2rphi_bcon() const { return chi2rphi_bcon_; }
0168     float chi2_bcon() const { return chi2rphi_bcon_ + chi2rz_; }
0169     float chi2dof_bcon() const { return (this->chi2_bcon()) / this->numDOF_bcon(); }
0170 
0171     //--- Get phi sector and eta region used by track finding code that this track is in.
0172     unsigned int iPhiSec() const { return iPhiSec_; }
0173     unsigned int iEtaReg() const { return iEtaReg_; }
0174 
0175     //--- Opto-link ID used to send this track from HT to Track Fitter
0176     unsigned int optoLinkID() const { return optoLinkID_; }
0177 
0178     //--- Get whether the track has been rejected or accepted by the fit
0179 
0180     bool accepted() const { return accepted_; }
0181 
0182     // Digitize track and degrade helix parameter resolution according to effect of digitisation.
0183     void digitizeTrack(const std::string& fitterName);
0184 
0185     // Access to detailed info about digitized track
0186     const DigitalTrack* digitaltrack() const { return digitalTrack_.get(); }
0187 
0188   private:
0189     //--- Configuration parameters
0190     const Settings* settings_;
0191 
0192     //--- The 3D hough-transform track candidate which was fitted.
0193     const L1track3D* l1track3D_;
0194 
0195     //--- The stubs on the fitted track (can differ from those on HT track if fit kicked off stubs with bad residuals)
0196     std::vector<const Stub*> stubs_;
0197     unsigned int nLayers_;
0198 
0199     //--- Bit-encoded hit pattern (where layer number assigned by increasing distance from origin, according to layers track expected to cross).
0200     unsigned int hitPattern_;
0201 
0202     //--- The fitted helix parameters and fit chi-squared.
0203     float qOverPt_;
0204     float d0_;
0205     float phi0_;
0206     float z0_;
0207     float tanLambda_;
0208     float chi2rphi_;
0209     float chi2rz_;
0210 
0211     //--- Ditto with beam-spot constraint applied in case of 5-parameter fit, plus boolean to indicate
0212     bool done_bcon_;
0213     float qOverPt_bcon_;
0214     float d0_bcon_;
0215     float phi0_bcon_;
0216     float chi2rphi_bcon_;
0217 
0218     //--- The number of helix parameters being fitted (=5 if d0 is fitted or =4 if d0 is not fitted).
0219     unsigned int nHelixParam_;
0220 
0221     //--- Phi sector and eta region used track finding code that this track was in.
0222     unsigned int iPhiSec_;
0223     unsigned int iEtaReg_;
0224     //--- Opto-link ID from HT to Track Fitter.
0225     unsigned int optoLinkID_;
0226 
0227     //--- Information about its association (if any) to a truth Tracking Particle.
0228     const TP* matchedTP_;
0229     std::vector<const Stub*> matchedStubs_;
0230     unsigned int nMatchedLayers_;
0231 
0232     //--- Info specific to KF fitter.
0233     unsigned int nSkippedLayers_;
0234     unsigned int numUpdateCalls_;
0235     //--- Info specific to LR fitter.
0236     int numIterations_;
0237     std::string lostMatchingState_;
0238     std::unordered_map<std::string, int> stateCalls_;
0239 
0240     std::shared_ptr<DigitalTrack> digitalTrack_;  // Class used to digitize track if required.
0241 
0242     //--- Has the track fit declared this to be a valid track?
0243     bool accepted_;
0244   };
0245 
0246 }  // namespace tmtt
0247 
0248 #endif