Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-12 23:00:01

0001 /*
0002 Kalman Filter L1 Muon algorithm
0003 Michalis Bachtis (UCLA)
0004 Sep. 2017
0005 
0006 */
0007 
0008 #ifndef L1TMuonBarrelKalmanAlgo_H
0009 #define L1TMuonBarrelKalmanAlgo_H
0010 
0011 #include "DataFormats/L1TMuon/interface/L1MuKBMTrack.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "L1Trigger/L1TMuonBarrel/interface/L1TMuonBarrelKalmanLUTs.h"
0014 #include "DataFormats/L1TMuon/interface/RegionalMuonCand.h"
0015 
0016 class L1TMuonBarrelKalmanAlgo {
0017 public:
0018   typedef ROOT::Math::SVector<double, 2> Vector2;
0019   typedef ROOT::Math::SMatrix<double, 2, 2, ROOT::Math::MatRepSym<double, 2> > CovarianceMatrix2;
0020   typedef ROOT::Math::SMatrix<double, 3, 2> Matrix32;
0021   typedef ROOT::Math::SMatrix<double, 2, 3> Matrix23;
0022   typedef ROOT::Math::SMatrix<double, 1, 3> Matrix13;
0023   typedef ROOT::Math::SMatrix<double, 3, 1> Matrix31;
0024   typedef ROOT::Math::SMatrix<double, 3, 3> Matrix33;
0025 
0026   L1TMuonBarrelKalmanAlgo(const edm::ParameterSet& settings);
0027   std::pair<bool, L1MuKBMTrack> chain(const L1MuKBMTCombinedStubRef&, const L1MuKBMTCombinedStubRefVector&);
0028 
0029   L1MuKBMTrackCollection clean(const L1MuKBMTrackCollection&, uint);
0030 
0031   void addBMTFMuon(int, const L1MuKBMTrack&, std::unique_ptr<l1t::RegionalMuonCandBxCollection>&);
0032   l1t::RegionalMuonCand convertToBMTF(const L1MuKBMTrack& track);
0033 
0034 private:
0035   bool verbose_;
0036   std::pair<bool, uint> match(const L1MuKBMTCombinedStubRef&, const L1MuKBMTCombinedStubRefVector&, int);
0037   int correctedPhi(const L1MuKBMTCombinedStubRef&, int);
0038   int correctedPhiB(const L1MuKBMTCombinedStubRef&);
0039   void propagate(L1MuKBMTrack&);
0040   void updateEta(L1MuKBMTrack&, const L1MuKBMTCombinedStubRef&);
0041   bool update(L1MuKBMTrack&, const L1MuKBMTCombinedStubRef&, int, int);
0042   bool updateOffline(L1MuKBMTrack&, const L1MuKBMTCombinedStubRef&);
0043   bool updateOffline1D(L1MuKBMTrack&, const L1MuKBMTCombinedStubRef&);
0044   bool updateLUT(L1MuKBMTrack&, const L1MuKBMTCombinedStubRef&, int, int);
0045   void vertexConstraint(L1MuKBMTrack&);
0046   void vertexConstraintOffline(L1MuKBMTrack&);
0047   void vertexConstraintLUT(L1MuKBMTrack&);
0048   int hitPattern(const L1MuKBMTrack&);
0049   int customBitmask(unsigned int, unsigned int, unsigned int, unsigned int);
0050   bool getBit(int, int);
0051   void setFloatingPointValues(L1MuKBMTrack&, bool);
0052   int phiAt2(const L1MuKBMTrack& track);
0053   bool estimateChiSquare(L1MuKBMTrack&);
0054   void estimateCompatibility(L1MuKBMTrack&);
0055   int rank(const L1MuKBMTrack&);
0056   int wrapAround(int, int);
0057   std::pair<bool, uint> getByCode(const L1MuKBMTrackCollection& tracks, int mask);
0058   std::map<int, int> trackAddress(const L1MuKBMTrack&, int&);
0059   int encode(bool ownwheel, int sector, bool tag);
0060   uint twosCompToBits(int);
0061   int fp_product(float, int, uint);
0062 
0063   uint etaStubRank(const L1MuKBMTCombinedStubRef&);
0064 
0065   void calculateEta(L1MuKBMTrack& track);
0066   uint matchAbs(std::map<uint, uint>&, uint, uint);
0067 
0068   //LUT service
0069   L1TMuonBarrelKalmanLUTs* lutService_;
0070   int ptLUT(int K);
0071 
0072   //Initial Curvature
0073   std::vector<double> initK_;
0074   std::vector<double> initK2_;
0075 
0076   //propagation coefficients
0077   std::vector<double> eLoss_;
0078   std::vector<double> aPhi_;
0079   std::vector<double> aPhiB_;
0080   std::vector<double> aPhiBNLO_;
0081   std::vector<double> bPhi_;
0082   std::vector<double> bPhiB_;
0083   double phiAt2_;
0084   std::vector<double> etaLUT0_;
0085   std::vector<double> etaLUT1_;
0086 
0087   //Chi Square estimator input
0088   uint globalChi2Cut_;
0089   uint globalChi2CutLimit_;
0090   std::vector<double> chiSquare_;
0091   std::vector<int> chiSquareCutPattern_;
0092   std::vector<int> chiSquareCutCurv_;
0093   std::vector<int> chiSquareCut_;
0094 
0095   std::vector<double> trackComp_;
0096   std::vector<double> trackCompErr1_;
0097   std::vector<double> trackCompErr2_;
0098   std::vector<int> trackCompPattern_;
0099   std::vector<int> trackCompCutCurv_;
0100   std::vector<int> trackCompCut_;
0101   std::vector<int> chiSquareCutTight_;
0102 
0103   //bitmasks to run== diferent combinations for a given seed in a given station
0104   std::vector<int> combos4_;
0105   std::vector<int> combos3_;
0106   std::vector<int> combos2_;
0107   std::vector<int> combos1_;
0108 
0109   //bits for fixed point precision
0110   static const int BITSCURV = 14;
0111   static const int BITSPHI = 12;
0112   static const int BITSPHIB = 13;
0113   static const int BITSPARAM = 14;
0114   static const int GAIN_0 = 9;
0115   static const int GAIN_0INT = 6;
0116   static const int GAIN_4 = 9;
0117   static const int GAIN_4INT = 4;
0118   static const int GAIN_V0 = 9;
0119   static const int GAIN_V0INT = 3;
0120 
0121   static const int GAIN2_0 = 12;
0122   static const int GAIN2_0INT = 8;
0123   static const int GAIN2_1 = 12;
0124   static const int GAIN2_1INT = 4;
0125   static const int GAIN2_4 = 12;
0126   static const int GAIN2_4INT = 4;
0127   static const int GAIN2_5 = 12;
0128   static const int GAIN2_5INT = 0;
0129   //STUFF NOT USED IN THE FIRMWARE BUT ONLY FOR DEBUGGING
0130   ///////////////////////////////////////////////////////
0131 
0132   bool useOfflineAlgo_;
0133   std::vector<double> mScatteringPhi_;
0134   std::vector<double> mScatteringPhiB_;
0135   //point resolution for phi
0136   double pointResolutionPhi_;
0137   //point resolution for phiB
0138   double pointResolutionPhiB_;
0139   std::vector<double> pointResolutionPhiBH_;
0140   std::vector<double> pointResolutionPhiBL_;
0141   //double pointResolutionPhiB_;
0142   //point resolution for vertex
0143   double pointResolutionVertex_;
0144   //Toggle for the new quality calculation in the emulator
0145   bool useNewQualityCalculation_;
0146 
0147   //Sorter
0148   class StubSorter {
0149   public:
0150     StubSorter(uint sector) { sec_ = sector; }
0151 
0152     bool operator()(const L1MuKBMTCombinedStubRef& a, const L1MuKBMTCombinedStubRef& b) {
0153       if (correctedPhi(a) < correctedPhi(b))
0154         return true;
0155       return false;
0156     }
0157 
0158   private:
0159     int sec_;
0160     int correctedPhi(const L1MuKBMTCombinedStubRef& stub) {
0161       if (stub->scNum() == sec_)
0162         return stub->phi();
0163       else if ((stub->scNum() == sec_ - 1) || (stub->scNum() == 11 && sec_ == 0))
0164         return stub->phi() - 2144;
0165       else if ((stub->scNum() == sec_ + 1) || (stub->scNum() == 0 && sec_ == 11))
0166         return stub->phi() + 2144;
0167       return 0;
0168     }
0169   };
0170 };
0171 #endif