Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef L1Trigger_TrackFindingTMTT_L1track3D_h
0002 #define L1Trigger_TrackFindingTMTT_L1track3D_h
0003 
0004 #include "DataFormats/Math/interface/deltaPhi.h"
0005 #include "L1Trigger/TrackFindingTMTT/interface/L1trackBase.h"
0006 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0007 #include "L1Trigger/TrackFindingTMTT/interface/Utility.h"
0008 #include "L1Trigger/TrackFindingTMTT/interface/TP.h"
0009 #include "L1Trigger/TrackFindingTMTT/interface/Stub.h"
0010 #include "L1Trigger/TrackFindingTMTT/interface/Sector.h"
0011 #include "L1Trigger/TrackFindingTMTT/interface/HTrphi.h"
0012 
0013 #include <vector>
0014 #include <string>
0015 #include <unordered_set>
0016 #include <utility>
0017 
0018 //=== L1 track candidate found in 3 dimensions.
0019 //=== Gives access to all stubs on track and to its 3D helix parameters.
0020 //=== Also calculates & gives access to associated truth particle (Tracking Particle) if any.
0021 
0022 namespace tmtt {
0023 
0024   class L1track3D : public L1trackBase {
0025   public:
0026     // Seeding layers of tracklet pattern reco.
0027     enum TrackletSeedType { L1L2, L2L3, L3L4, L5L6, D1D2, D3D4, L1D1, L2D1, L3L4L2, L5L6L4, L2L3D1, D1D2L2, NONE };
0028 
0029   public:
0030     L1track3D(const Settings* settings,
0031               const std::vector<Stub*>& stubs,
0032               std::pair<unsigned int, unsigned int> cellLocationHT,
0033               std::pair<float, float> helixRphi,
0034               std::pair<float, float> helixRz,
0035               float helixD0,
0036               unsigned int iPhiSec,
0037               unsigned int iEtaReg,
0038               unsigned int optoLinkID,
0039               bool mergedHTcell)
0040         : L1trackBase(),
0041           settings_(settings),
0042           stubs_(stubs),
0043           stubsConst_(stubs_.begin(), stubs_.end()),
0044           cellLocationHT_(cellLocationHT),
0045           helixRphi_(helixRphi),
0046           helixRz_(helixRz),
0047           helixD0_(helixD0),
0048           iPhiSec_(iPhiSec),
0049           iEtaReg_(iEtaReg),
0050           optoLinkID_(optoLinkID),
0051           mergedHTcell_(mergedHTcell),
0052           seedLayerType_(TrackletSeedType::NONE),
0053           seedPS_(999) {
0054       nLayers_ = Utility::countLayers(settings, stubs_);  // Count tracker layers these stubs are in
0055       matchedTP_ = Utility::matchingTP(settings,
0056                                        stubs_,
0057                                        nMatchedLayers_,
0058                                        matchedStubs_);  // Find associated truth particle & calculate info about match.
0059     }
0060 
0061     // TMTT tracking: constructor
0062 
0063     L1track3D(const Settings* settings,
0064               const std::vector<Stub*>& stubs,
0065               std::pair<unsigned int, unsigned int> cellLocationHT,
0066               std::pair<float, float> helixRphi,
0067               std::pair<float, float> helixRz,
0068               unsigned int iPhiSec,
0069               unsigned int iEtaReg,
0070               unsigned int optoLinkID,
0071               bool mergedHTcell)
0072         : L1track3D(
0073               settings, stubs, cellLocationHT, helixRphi, helixRz, 0.0, iPhiSec, iEtaReg, optoLinkID, mergedHTcell) {}
0074 
0075     ~L1track3D() override = default;
0076 
0077     //--- Set/get optional info for tracklet tracks.
0078 
0079     // Tracklet seeding layer pair (from Tracklet::calcSeedIndex())
0080     void setSeedLayerType(unsigned int seedLayerType) { seedLayerType_ = static_cast<TrackletSeedType>(seedLayerType); }
0081     TrackletSeedType seedLayerType() const { return seedLayerType_; }
0082 
0083     // Tracklet seed stub pair uses PS modules (from FPGATracket::PSseed())
0084     void setSeedPS(unsigned int seedPS) { seedPS_ = seedPS; }
0085     unsigned int seedPS() const { return seedPS_; }
0086 
0087     // Best stub (stub with smallest Phi residual in each layer/disk)
0088     void setBestStubs(std::unordered_set<const Stub*> bestStubs) { bestStubs_ = bestStubs; }
0089     std::unordered_set<const Stub*> bestStubs() const { return bestStubs_; }
0090 
0091     //--- Get information about the reconstructed track.
0092 
0093     // Get stubs on track candidate.
0094     const std::vector<const Stub*>& stubsConst() const override { return stubsConst_; }
0095     const std::vector<Stub*>& stubs() const override { return stubs_; }
0096     // Get number of stubs on track candidate.
0097     unsigned int numStubs() const override { return stubs_.size(); }
0098     // Get number of tracker layers these stubs are in.
0099     unsigned int numLayers() const override { return nLayers_; }
0100     // Get cell location of track candidate in r-phi Hough Transform array in units of bin number.
0101     std::pair<unsigned int, unsigned int> cellLocationHT() const override { return cellLocationHT_; }
0102     // The two conventionally agreed track helix parameters relevant in r-phi plane. i.e. (q/Pt, phi0)
0103     std::pair<float, float> helixRphi() const { return helixRphi_; }
0104     // The two conventionally agreed track helix parameters relevant in r-z plane. i.e. (z0, tan_lambda)
0105     std::pair<float, float> helixRz() const { return helixRz_; }
0106 
0107     //--- Return chi variables, (both digitized & undigitized), which are the stub coords. relative to track.
0108 
0109     std::vector<float> chiPhi() {
0110       std::vector<float> result;
0111       for (const Stub* s : stubs_) {
0112         float chi_phi = reco::deltaPhi(s->phi(), this->phi0() - s->r() * this->qOverPt() * settings_->invPtToDphi());
0113         result.push_back(chi_phi);
0114       }
0115       return result;
0116     }
0117 
0118     std::vector<int> chiPhiDigi() {
0119       std::vector<int> result;
0120       const float phiMult = pow(2, settings_->phiSBits()) / settings_->phiSRange();
0121       for (const float& chi_phi : this->chiPhi()) {
0122         int iDigi_chi_phi = floor(chi_phi * phiMult);
0123         result.push_back(iDigi_chi_phi);
0124       }
0125       return result;
0126     }
0127 
0128     std::vector<float> chiZ() {
0129       std::vector<float> result;
0130       for (const Stub* s : stubs_) {
0131         float chi_z = s->z() - (this->z0() + s->r() * this->tanLambda());
0132         result.push_back(chi_z);
0133       }
0134       return result;
0135     }
0136 
0137     std::vector<int> chiZDigi() {
0138       std::vector<int> result;
0139       const float zMult = pow(2, settings_->zBits()) / settings_->zRange();
0140       for (const float& chi_z : this->chiZ()) {
0141         int iDigi_chi_z = floor(chi_z * zMult);
0142         result.push_back(iDigi_chi_z);
0143       }
0144       return result;
0145     }
0146 
0147     //--- User-friendlier access to the helix parameters.
0148 
0149     float qOverPt() const override { return helixRphi_.first; }
0150     float charge() const { return (this->qOverPt() > 0 ? 1 : -1); }
0151     float invPt() const { return std::abs(this->qOverPt()); }
0152     // Protect pt against 1/pt = 0.
0153     float pt() const {
0154       constexpr float small = 1.0e-6;
0155       return 1. / (small + this->invPt());
0156     }
0157     float d0() const { return helixD0_; }  // Hough transform assumes d0 = 0.
0158     float phi0() const override { return helixRphi_.second; }
0159     float z0() const { return helixRz_.first; }
0160     float tanLambda() const { return helixRz_.second; }
0161     float theta() const { return atan2(1., this->tanLambda()); }  // Use atan2 to ensure 0 < theta < pi.
0162     float eta() const { return -log(tan(0.5 * this->theta())); }
0163 
0164     // Phi and z coordinates at which track crosses "chosenR" values used by r-phi HT and rapidity sectors respectively.
0165     float phiAtChosenR() const {
0166       return reco::deltaPhi(this->phi0() - (settings_->invPtToDphi() * settings_->chosenRofPhi()) * this->qOverPt(),
0167                             0.);
0168     }
0169     float zAtChosenR() const {
0170       return (this->z0() + (settings_->chosenRofZ()) * this->tanLambda());
0171     }  // neglects transverse impact parameter & track curvature.
0172 
0173     //--- Get phi sector and eta region used by track finding code that this track is in.
0174     unsigned int iPhiSec() const override { return iPhiSec_; }
0175     unsigned int iEtaReg() const override { return iEtaReg_; }
0176 
0177     //--- Opto-link ID used to send this track from HT to Track Fitter
0178     unsigned int optoLinkID() const override { return optoLinkID_; }
0179 
0180     //--- Was this track produced from a marged HT cell (e.g. 2x2)?
0181     bool mergedHTcell() const { return mergedHTcell_; }
0182 
0183     //--- Get information about its association (if any) to a truth Tracking Particle.
0184 
0185     // Get best matching tracking particle (=nullptr if none).
0186     const TP* matchedTP() const override { return matchedTP_; }
0187     // Get the matched stubs with this Tracking Particle
0188     const std::vector<const Stub*>& matchedStubs() const override { return matchedStubs_; }
0189     // Get number of matched stubs with this Tracking Particle
0190     unsigned int numMatchedStubs() const override { return matchedStubs_.size(); }
0191     // Get number of tracker layers with matched stubs with this Tracking Particle
0192     unsigned int numMatchedLayers() const override { return nMatchedLayers_; }
0193     // Get purity of stubs on track candidate (i.e. fraction matching best Tracking Particle)
0194     float purity() const { return numMatchedStubs() / float(numStubs()); }
0195 
0196     //--- For debugging purposes.
0197 
0198     // Remove incorrect stubs from the track using truth information.
0199     // Also veto tracks where the HT cell estimated from the true helix parameters is inconsistent with the cell the HT found the track in, (since probable duplicates).
0200     // Also veto tracks that match a truth particle not used for the algo efficiency measurement.
0201     // Return a boolean indicating if the track should be kept. (i.e. Is genuine & non-duplicate).
0202     bool cheat() {
0203       bool keep = false;
0204 
0205       std::vector<Stub*> stubsSel;
0206       if (matchedTP_ != nullptr) {  // Genuine track
0207         for (Stub* s : stubs_) {
0208           const TP* tp = s->assocTP();
0209           if (tp != nullptr) {
0210             if (matchedTP_->index() == tp->index()) {
0211               stubsSel.push_back(s);  // This stub was produced by same truth particle as rest of track, so keep it.
0212             }
0213           }
0214         }
0215       }
0216       stubs_ = stubsSel;
0217       stubsConst_ = std::vector<const Stub*>(stubs_.begin(), stubs_.end());
0218 
0219       nLayers_ = Utility::countLayers(settings_, stubs_);  // Count tracker layers these stubs are in
0220       matchedTP_ = Utility::matchingTP(settings_,
0221                                        stubs_,
0222                                        nMatchedLayers_,
0223                                        matchedStubs_);  // Find associated truth particle & calculate info about match.
0224 
0225       bool genuine = (matchedTP_ != nullptr);
0226 
0227       if (genuine && matchedTP_->useForAlgEff()) {
0228         Sector secTmp(settings_, iPhiSec_, iEtaReg_);
0229         HTrphi htRphiTmp(settings_, iPhiSec_, iEtaReg_, secTmp.etaMin(), secTmp.etaMax(), secTmp.phiCentre());
0230         std::pair<unsigned int, unsigned int> trueCell = htRphiTmp.trueCell(matchedTP_);
0231 
0232         std::pair<unsigned int, unsigned int> htCell = this->cellLocationHT();
0233         bool consistent = (htCell == trueCell);  // If true, track is probably not a duplicate.
0234         if (mergedHTcell_) {
0235           // If this is a merged cell, check other elements of merged cell.
0236           std::pair<unsigned int, unsigned int> htCell10(htCell.first + 1, htCell.second);
0237           std::pair<unsigned int, unsigned int> htCell01(htCell.first, htCell.second + 1);
0238           std::pair<unsigned int, unsigned int> htCell11(htCell.first + 1, htCell.second + 1);
0239           if (htCell10 == trueCell)
0240             consistent = true;
0241           if (htCell01 == trueCell)
0242             consistent = true;
0243           if (htCell11 == trueCell)
0244             consistent = true;
0245         }
0246         if (consistent)
0247           keep = true;
0248       }
0249 
0250       return keep;  // Indicate if track should be kept.
0251     }
0252 
0253   private:
0254     //--- Configuration parameters
0255     const Settings* settings_;
0256 
0257     //--- Information about the reconstructed track.
0258     std::vector<Stub*> stubs_;
0259     std::vector<const Stub*> stubsConst_;
0260     std::unordered_set<const Stub*> bestStubs_;
0261     unsigned int nLayers_;
0262     std::pair<unsigned int, unsigned int> cellLocationHT_;
0263     std::pair<float, float> helixRphi_;
0264     std::pair<float, float> helixRz_;
0265     float helixD0_;
0266     unsigned int iPhiSec_;
0267     unsigned int iEtaReg_;
0268     unsigned int optoLinkID_;
0269     bool mergedHTcell_;
0270 
0271     //--- Optional info used for tracklet tracks.
0272     TrackletSeedType seedLayerType_;
0273     unsigned int seedPS_;
0274 
0275     //--- Information about its association (if any) to a truth Tracking Particle.
0276     const TP* matchedTP_;
0277     std::vector<const Stub*> matchedStubs_;
0278     unsigned int nMatchedLayers_;
0279   };
0280 
0281 }  // namespace tmtt
0282 
0283 #endif