Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-16 02:42:43

0001 #ifndef DataFormats_L1TMuonPhase2_KMTFTrack_h
0002 #define DataFormats_L1TMuonPhase2_KMTFTrack_h
0003 
0004 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0005 #include "DataFormats/L1TMuonPhase2/interface/MuonStub.h"
0006 #include "DataFormats/L1Trigger/interface/BXVector.h"
0007 
0008 namespace l1t {
0009 
0010   class KMTFTrack;
0011   typedef std::vector<KMTFTrack> KMTFTrackCollection;
0012   typedef BXVector<KMTFTrack> KMTFTrackBxCollection;
0013 
0014   class KMTFTrack : public reco::LeafCandidate {
0015   public:
0016     KMTFTrack()
0017         : reco::LeafCandidate(-1, reco::LeafCandidate::PolarLorentzVector(0.1, 0.0, 0.0, 0.105)),
0018           unconstrainedP4_(reco::LeafCandidate::PolarLorentzVector(0.1, 0.0, 0.0, 0.105)),
0019           covariance_(std::vector<double>(6, 0.0)),
0020           curvVertex_(0),
0021           ptC_(0),
0022           phiVertex_(0),
0023           dxy_(0),
0024           curvMuon_(0),
0025           ptU_(0),
0026           phiMuon_(0),
0027           phiBMuon_(0),
0028           curv_(0),
0029           phi_(0),
0030           phiB_(0),
0031           coarseEta_(0),
0032           approxPromptChi2_(0),
0033           approxPromptErrChi2_(0),
0034           approxDispChi2_(0),
0035           approxDispErrChi2_(0),
0036           hitPattern_(0),
0037           step_(1),
0038           sector_(0),
0039           wheel_(0),
0040           quality_(0),
0041           hasFineEta_(false),
0042           bx_(0),
0043           rankPrompt_(0),
0044           rankDisp_(0),
0045           idFlag_(0) {}
0046 
0047     ~KMTFTrack() override = default;
0048 
0049     KMTFTrack(const l1t::MuonStubRef& seed, int phi, int phiB)
0050         : reco::LeafCandidate(-1, reco::LeafCandidate::PolarLorentzVector(0.1, 0.0, 0.0, 0.105)),
0051           unconstrainedP4_(reco::LeafCandidate::PolarLorentzVector(0.1, 0.0, 0.0, 0.105)),
0052           covariance_(std::vector<double>(6, 0.0)),
0053           curvVertex_(0),
0054           ptC_(0),
0055           phiVertex_(0),
0056           dxy_(0),
0057           curvMuon_(0),
0058           ptU_(0),
0059           phiMuon_(0),
0060           phiBMuon_(0),
0061           curv_(0),
0062           phi_(phi),
0063           phiB_(phiB),
0064           coarseEta_(0),
0065           approxPromptChi2_(0),
0066           approxPromptErrChi2_(0),
0067           approxDispChi2_(0),
0068           approxDispErrChi2_(0),
0069           hitPattern_(0),
0070           step_(seed->depthRegion()),
0071           sector_(seed->phiRegion()),
0072           wheel_(seed->etaRegion()),
0073           quality_(seed->quality()),
0074           hasFineEta_(false),
0075           bx_(seed->bxNum()),
0076           rankPrompt_(0),
0077           rankDisp_(0),
0078           idFlag_(0) {
0079       stubs_.push_back(seed);
0080       residuals_.push_back(0);
0081       residuals_.push_back(0);
0082       residuals_.push_back(0);
0083     }
0084 
0085     reco::LeafCandidate::PolarLorentzVector displacedP4() const { return unconstrainedP4_; }
0086 
0087     //unconstrained pt
0088     int ptDisplaced() const { return ptU_; }
0089     //unconstrained curvature at station 1
0090     int curvatureAtMuon() const { return curvMuon_; }
0091     //unconstrained phi at station 1
0092     int phiAtMuon() const { return phiMuon_; }
0093     //unconstrained phiB at station 1
0094     int phiBAtMuon() const { return phiBMuon_; }
0095 
0096     //constrained pt
0097     int ptPrompt() const { return ptC_; }
0098     //Constrained curvature at vertex
0099     int curvatureAtVertex() const { return curvVertex_; }
0100     //constrained phi at the vertex
0101     int phiAtVertex() const { return phiVertex_; }
0102     //Impact parameter as calculated from the muon track
0103     int dxy() const { return dxy_; }
0104     //Unconstrained curvature at the Muon systen
0105     int curvature() const { return curv_; }
0106     //Unconstrained phi at the Muon systen
0107     int positionAngle() const { return phi_; }
0108     //Unconstrained bending angle at the Muon systen
0109     int bendingAngle() const { return phiB_; }
0110     //Coarse eta caluclated only using phi segments
0111     int coarseEta() const { return coarseEta_; }
0112     //Approximate Chi2 metrics
0113     int approxPromptChi2() const { return approxPromptChi2_; }
0114     int approxPromptErrChi2() const { return approxPromptErrChi2_; }
0115     int approxDispChi2() const { return approxDispChi2_; }
0116     int approxDispErrChi2() const { return approxDispErrChi2_; }
0117 
0118     int hitPattern() const { return hitPattern_; }
0119     //step;
0120     int step() const { return step_; }
0121     //sector;
0122     int sector() const { return sector_; }
0123     //wheel
0124     int wheel() const { return wheel_; }
0125     //quality
0126     int quality() const { return quality_; }
0127 
0128     //fine eta
0129     int fineEta() const { return fineEta_; }
0130     bool hasFineEta() const { return hasFineEta_; }
0131 
0132     //BX
0133     int bx() const { return bx_; }
0134 
0135     //rank
0136     int rankPrompt() const { return rankPrompt_; }
0137     int rankDisp() const { return rankDisp_; }
0138 
0139     int id() const { return idFlag_; }
0140 
0141     //Associated stubs
0142     const l1t::MuonStubRefVector& stubs() const { return stubs_; }
0143 
0144     //get Kalman gain
0145     const std::vector<float>& kalmanGain(unsigned int step) const {
0146       switch (step) {
0147         case 3:
0148           return kalmanGain3_;
0149         case 2:
0150           return kalmanGain2_;
0151         case 1:
0152           return kalmanGain1_;
0153         case 0:
0154           return kalmanGain0_;
0155       }
0156       return kalmanGain0_;
0157     }
0158 
0159     //get covariance
0160     const std::vector<double>& covariance() const { return covariance_; }
0161 
0162     //get residual
0163     int residual(uint i) const { return residuals_[i]; }
0164 
0165     //check overlap
0166     bool overlapTrack(const KMTFTrack& other) const {
0167       for (const auto& s1 : stubs_) {
0168         for (const auto& s2 : other.stubs()) {
0169           if (s1->phiRegion() == s2->phiRegion() && s1->etaRegion() == s2->etaRegion() &&
0170               s1->depthRegion() == s2->depthRegion() && s1->id() == s2->id())
0171             return true;
0172         }
0173       }
0174       return false;
0175     }
0176 
0177     bool operator==(const KMTFTrack& t2) const {
0178       if (this->stubs().size() != t2.stubs().size())
0179         return false;
0180       for (unsigned int i = 0; i < this->stubs().size(); ++i) {
0181         const l1t::MuonStubRef& s1 = this->stubs()[i];
0182         const l1t::MuonStubRef& s2 = t2.stubs()[i];
0183         if (s1->phiRegion() != s2->phiRegion() || s1->etaRegion() != s2->etaRegion() ||
0184             s1->depthRegion() != s2->depthRegion() || s1->id() != s2->id() || s1->tfLayer() != s2->tfLayer())
0185           return false;
0186       }
0187       return true;
0188     }
0189 
0190     //Set coordinates general
0191     void setCoordinates(int step, int curv, int phi, int phiB) {
0192       step_ = step;
0193       curv_ = curv;
0194       phiB_ = phiB;
0195       phi_ = phi;
0196     }
0197 
0198     void setCoordinatesAtVertex(int curv, int phi, int dxy) {
0199       curvVertex_ = curv;
0200       phiVertex_ = phi;
0201       dxy_ = dxy;
0202     }
0203 
0204     void setCoordinatesAtMuon(int curv, int phi, int phiB) {
0205       curvMuon_ = curv;
0206       phiMuon_ = phi;
0207       phiBMuon_ = phiB;
0208     }
0209 
0210     void setPt(int ptC, int ptU) {
0211       ptC_ = ptC;
0212       ptU_ = ptU;
0213     }
0214 
0215     void setCoarseEta(int eta) { coarseEta_ = eta; }
0216 
0217     void setHitPattern(int pattern) { hitPattern_ = pattern; }
0218 
0219     void setApproxChi2(int chi, int chiErr, bool prompt) {
0220       if (prompt) {
0221         approxPromptChi2_ = chi;
0222         approxPromptErrChi2_ = chiErr;
0223       } else {
0224         approxDispChi2_ = chi;
0225         approxDispErrChi2_ = chiErr;
0226       }
0227     }
0228 
0229     void setPtEtaPhi(double pt, double eta, double phi) {
0230       PolarLorentzVector v(pt, eta, phi, 0.105);
0231       setP4(v);
0232     }
0233     void setPtEtaPhiDisplaced(double pt, double eta, double phi) {
0234       unconstrainedP4_.SetPt(pt);
0235       unconstrainedP4_.SetEta(eta);
0236       unconstrainedP4_.SetPhi(phi);
0237     }
0238 
0239     void addStub(const l1t::MuonStubRef& stub) {
0240       if (stub->quality() < quality_)
0241         quality_ = stub->quality();
0242       stubs_.push_back(stub);
0243     }
0244 
0245     void setStubs(const l1t::MuonStubRefVector& stubs) { stubs_ = stubs; }
0246 
0247     void setRank(int rank, bool vertex) {
0248       if (vertex)
0249         rankPrompt_ = rank;
0250       else
0251         rankDisp_ = rank;
0252     }
0253 
0254     void setIDFlag(bool passPrompt, bool passDisp) {
0255       unsigned p0 = 0;
0256       unsigned p1 = 0;
0257 
0258       if (passPrompt)
0259         p0 = 1;
0260       if (passDisp)
0261         p1 = 2;
0262 
0263       idFlag_ = p0 | p1;
0264     }
0265 
0266     void setKalmanGain(
0267         unsigned int step, unsigned int K, float a1, float a2, float a3 = 0, float a4 = 0, float a5 = 0, float a6 = 0) {
0268       switch (step) {
0269         case 3:
0270           kalmanGain3_.push_back(K);
0271           kalmanGain3_.push_back(a1);
0272           kalmanGain3_.push_back(a2);
0273           kalmanGain3_.push_back(a3);
0274           kalmanGain3_.push_back(a4);
0275           kalmanGain3_.push_back(a5);
0276           kalmanGain3_.push_back(a6);
0277           break;
0278         case 2:
0279           kalmanGain2_.push_back(K);
0280           kalmanGain2_.push_back(a1);
0281           kalmanGain2_.push_back(a2);
0282           kalmanGain2_.push_back(a3);
0283           kalmanGain2_.push_back(a4);
0284           kalmanGain2_.push_back(a5);
0285           kalmanGain2_.push_back(a6);
0286           break;
0287         case 1:
0288           kalmanGain1_.push_back(K);
0289           kalmanGain1_.push_back(a1);
0290           kalmanGain1_.push_back(a2);
0291           kalmanGain1_.push_back(a3);
0292           kalmanGain1_.push_back(a4);
0293           kalmanGain1_.push_back(a5);
0294           kalmanGain1_.push_back(a6);
0295           break;
0296         case 0:
0297           kalmanGain0_.push_back(K);
0298           kalmanGain0_.push_back(a1);
0299           kalmanGain0_.push_back(a2);
0300           kalmanGain0_.push_back(a3);
0301           break;
0302 
0303         default:
0304           throw cms::Exception("WrongCondition") << "Critical ERROR on setting the Kalman gain\n";
0305       }
0306     }
0307 
0308     //set covariance
0309     void setCovariance(const CovarianceMatrix& c) {
0310       covariance_[0] = c(0, 0);
0311       covariance_[1] = c(0, 1);
0312       covariance_[2] = c(1, 1);
0313       covariance_[3] = c(0, 2);
0314       covariance_[4] = c(1, 2);
0315       covariance_[5] = c(2, 2);
0316     }
0317 
0318     //set fine eta
0319     void setFineEta(int eta) {
0320       fineEta_ = eta;
0321       hasFineEta_ = true;
0322     }
0323 
0324     //set residual
0325     void setResidual(uint i, int val) { residuals_[i] = val; }
0326 
0327   private:
0328     reco::LeafCandidate::PolarLorentzVector unconstrainedP4_;
0329 
0330     //Covariance matrix for studies
0331     std::vector<double> covariance_;
0332     l1t::MuonStubRefVector stubs_;
0333 
0334     //vertex coordinates
0335     int curvVertex_;
0336     int ptC_;
0337     int phiVertex_;
0338     int dxy_;
0339 
0340     //muon coordinates
0341     int curvMuon_;
0342     int ptU_;
0343     int phiMuon_;
0344     int phiBMuon_;
0345 
0346     //generic coordinates
0347     int curv_;
0348     int phi_;
0349     int phiB_;
0350     //common coordinates
0351     int coarseEta_;
0352 
0353     //Approximate Chi2 metric
0354     int approxPromptChi2_;
0355     int approxPromptErrChi2_;
0356     int approxDispChi2_;
0357     int approxDispErrChi2_;
0358 
0359     //phi bitmask
0360     int hitPattern_;
0361 
0362     //propagation step
0363     int step_;
0364 
0365     //sector
0366     int sector_;
0367     //wheel
0368     int wheel_;
0369 
0370     //quality
0371     int quality_;
0372 
0373     //Fine eta
0374     int fineEta_;
0375 
0376     //has fine eta?
0377     bool hasFineEta_;
0378 
0379     //BX
0380     int bx_;
0381 
0382     //rank
0383     int rankPrompt_;
0384     int rankDisp_;
0385 
0386     //flag
0387     int idFlag_;
0388 
0389     //Kalman Gain for making LUTs
0390     std::vector<float> kalmanGain0_;
0391     std::vector<float> kalmanGain1_;
0392     std::vector<float> kalmanGain2_;
0393     std::vector<float> kalmanGain3_;
0394 
0395     std::vector<int> residuals_;
0396   };
0397 
0398 }  // namespace l1t
0399 #endif