Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef L1T_OmtfP1_OMTFConfiguration_H
0002 #define L1T_OmtfP1_OMTFConfiguration_H
0003 
0004 #include <map>
0005 #include <set>
0006 #include <vector>
0007 #include <ostream>
0008 #include <memory>
0009 
0010 #include "L1Trigger/L1TMuonOverlapPhase1/interface/ProcConfigurationBase.h"
0011 #include "CondFormats/L1TObjects/interface/L1TMuonOverlapParams.h"
0012 #include "DataFormats/L1TMuon/interface/RegionalMuonCandFwd.h"
0013 
0014 //typedef int omtfPdfValueType; //normal emulation is with int type
0015 typedef float PdfValueType;  //but floats are needed for the PatternOptimizer
0016 
0017 namespace edm {
0018   class ParameterSet;
0019 }
0020 
0021 class RefHitDef {
0022 public:
0023   //FIXME: default values should be sonnected to configuration values
0024   RefHitDef(unsigned int aInput = 15,
0025             int aPhiMin = 5760,
0026             int aPhiMax = 5760,
0027             unsigned int aRegion = 99,
0028             unsigned int aRefLayer = 99);
0029 
0030 public:
0031   bool fitsRange(int iPhi) const;
0032 
0033   ///Hit input number within a cone
0034   unsigned int iInput;
0035 
0036   ///Region number assigned to this referecne hit
0037   unsigned int iRegion;
0038 
0039   ///Reference layer logic number (0-7)
0040   unsigned int iRefLayer;
0041 
0042   ///Local to processor phi range.
0043   ///Hit has to fit into this range to be assigned to this iRegion;
0044   std::pair<int, int> range;
0045 
0046   friend std::ostream& operator<<(std::ostream& out, const RefHitDef& aRefHitDef);
0047 };
0048 
0049 class OMTFConfiguration : public ProcConfigurationBase {
0050 public:
0051   typedef std::vector<std::pair<unsigned int, unsigned int> > vector1D_pair;
0052   typedef std::vector<vector1D_pair> vector2D_pair;
0053   typedef std::vector<vector2D_pair> vector3D_pair;
0054 
0055   typedef std::vector<int> vector1D;
0056   typedef std::vector<vector1D> vector2D;
0057   typedef std::vector<vector2D> vector3D;
0058   typedef std::vector<vector3D> vector4D;
0059 
0060   OMTFConfiguration() { ; };
0061 
0062   virtual void configure(const L1TMuonOverlapParams* omtfParams);
0063 
0064   void initCounterMatrices();
0065 
0066   ///Find logic region number using first input number
0067   ///and then local phi value. The input and phi
0068   ///ranges are taken from DB.
0069   unsigned int getRegionNumberFromMap(unsigned int iInput, unsigned int iRefLayer, int iPhi) const;
0070 
0071   ///Check if given referecne hit is
0072   ///in phi range for some logic cone.
0073   ///Care is needed arounf +Pi and +2Pi points
0074   bool isInRegionRange(int iPhiStart, unsigned int coneSize, int iPhi) const;
0075 
0076   ///Return global phi for beggining of given processor
0077   ///Uses minim phi over all reference layers.
0078   int globalPhiStart(unsigned int iProcessor) const;
0079 
0080   ///Return layer number encoding subsystem, and
0081   ///station number in a simple formula:
0082   /// aLayer+100*detId.subdetId()
0083   ///where aLayer is a layer number counting from vertex
0084   uint32_t getLayerNumber(uint32_t rawId) const;
0085 
0086   unsigned int fwVersion() const { return (rawParams.fwVersion() >> 16) & 0xFFFF; };
0087   unsigned int patternsVersion() const { return rawParams.fwVersion() & 0xFFFF; };
0088 
0089   const L1TMuonOverlapParams* getRawParams() const { return &rawParams; };
0090 
0091   float minPdfVal() const { return 0.001; };
0092   unsigned int nLayers() const override { return rawParams.nLayers(); };
0093   unsigned int nHitsPerLayer() const { return rawParams.nHitsPerLayer(); };
0094   unsigned int nRefLayers() const { return rawParams.nRefLayers(); };
0095   unsigned int nPhiBits() const { return rawParams.nPhiBits(); };
0096   unsigned int nPdfAddrBits() const { return rawParams.nPdfAddrBits(); };
0097   unsigned int nPdfBins() const { return pdfBins; };
0098   unsigned int nPdfValBits() const { return rawParams.nPdfValBits(); };
0099   int pdfMaxValue() const { return pdfMaxVal; };
0100   unsigned int nPhiBins() const override { return rawParams.nPhiBins(); };
0101   double omtfPhiUnit() const { return 2 * M_PI / nPhiBins(); }
0102   unsigned int nRefHits() const { return rawParams.nRefHits(); };
0103   unsigned int nTestRefHits() const { return rawParams.nTestRefHits(); };
0104   //processors number per detector side
0105   unsigned int nProcessors() const override { return rawParams.nProcessors(); };
0106   //total number of processors in the system
0107   unsigned int processorCnt() const { return 2 * rawParams.nProcessors(); };
0108   unsigned int nLogicRegions() const { return rawParams.nLogicRegions(); };
0109   unsigned int nInputs() const { return rawParams.nInputs(); };
0110   unsigned int nGoldenPatterns() const { return rawParams.nGoldenPatterns(); };
0111 
0112   const std::map<int, int>& getHwToLogicLayer() const { return hwToLogicLayer; }
0113   const std::map<int, int>& getLogicToHwLayer() const { return logicToHwLayer; }
0114   const std::map<int, int>& getLogicToLogic() const { return logicToLogic; }
0115   const std::set<int>& getBendingLayers() const { return bendingLayers; }
0116   const std::vector<int>& getRefToLogicNumber() const { return refToLogicNumber; }
0117 
0118   const std::vector<unsigned int>& getBarrelMin() const { return barrelMin; }
0119   const std::vector<unsigned int>& getBarrelMax() const { return barrelMax; }
0120   const std::vector<unsigned int>& getEndcap10DegMin() const { return endcap10DegMin; }
0121   const std::vector<unsigned int>& getEndcap10DegMax() const { return endcap10DegMax; }
0122   const std::vector<unsigned int>& getEndcap20DegMin() const { return endcap20DegMin; }
0123   const std::vector<unsigned int>& getEndcap20DegMax() const { return endcap20DegMax; }
0124 
0125   const std::vector<std::vector<int> >& getProcessorPhiVsRefLayer() const { return processorPhiVsRefLayer; }
0126   const std::vector<std::vector<std::vector<std::pair<int, int> > > >& getRegionPhisVsRefLayerVsInput() const {
0127     return regionPhisVsRefLayerVsInput;
0128   }
0129   const std::vector<std::vector<RefHitDef> >& getRefHitsDefs() const { return refHitsDefs; }
0130 
0131   const vector3D_pair& getConnections() const { return connections; };
0132 
0133   vector4D& getMeasurements4D() { return measurements4D; }
0134   vector4D& getMeasurements4Dref() { return measurements4Dref; }
0135 
0136   const vector4D& getMeasurements4D() const { return measurements4D; }
0137   const vector4D& getMeasurements4Dref() const { return measurements4Dref; }
0138 
0139   double ptUnit = 0.5;  // GeV/unit
0140   ///uGMT pt scale conversion
0141   double hwPtToGev(int hwPt) const override { return (hwPt - 1.) * ptUnit; }
0142 
0143   double uptUnit = 1;  // GeV/unit
0144   double hwUPtToGev(int hwPt) const override { return (hwPt - 1.) * uptUnit; }
0145 
0146   ///uGMT pt scale conversion: [0GeV, 0.5GeV) = 1 [0.5GeV, 1 Gev) = 2
0147   int ptGevToHw(double ptGev) const override { return (ptGev / ptUnit + 1); }
0148 
0149   int calcGlobalPhi(int locPhi, int proc) const;
0150 
0151   double etaUnit = 0.010875;  //=2.61/240 TODO value taken from the interface note, should be defined somewhere globally
0152   ///center of eta bin
0153   virtual double hwEtaToEta(int hwEta) const { return (hwEta * etaUnit); }
0154 
0155   int etaToHwEta(double eta) const override { return (eta / etaUnit); }
0156 
0157   static unsigned int eta2Bits(unsigned int eta);
0158 
0159   static unsigned int etaBits2HwEta(unsigned int eta);
0160 
0161   static int etaBit2Code(unsigned int bit);
0162 
0163   double phiGmtUnit = 2. * M_PI / 576;  //TODO from the interface note, should be defined somewhere globally
0164   //phi in radians
0165   virtual int phiToGlobalHwPhi(double phi) const { return std::floor(phi / phiGmtUnit); }
0166 
0167   //phi in radians
0168   virtual double hwPhiToGlobalPhi(int phi) const { return phi * phiGmtUnit; }
0169 
0170   ///iProcessor - 0...5
0171   ///phiRad [-pi,pi]
0172   ///return phi inside the processor
0173   int getProcScalePhi(unsigned int iProcessor, double phiRad) const;
0174 
0175   int getProcScalePhi(double phiRad, double procPhiZeroRad = 0) const override {
0176     return 0;  // TODO replace getProcScalePhi(unsigned int iProcessor, double phiRad) with this function
0177   }
0178 
0179   double procHwPhiToGlobalPhi(int procHwPhi, int procHwPhi0) const;
0180 
0181   int procPhiToGmtPhi(int procPhi) const {
0182     ///conversion factor from OMTF to uGMT scale is  5400/576 i.e. phiValue/=9.375;
0183     return floor(procPhi * 437. / (1 << 12));  // ie. use as in hw: 9.3729977
0184     //cannot be (procPhi * 437) >> 12, because this floor is needed
0185   }
0186 
0187   ///input phi should be in the hardware scale (nPhiBins units for 2pi), can be in range  -nPhiBins ... nPhiBins,
0188   //is converted to range -nPhiBins/2 +1 ... nPhiBins/2, pi = nPhiBins/2
0189   //virtual int foldPhi(int phi) const;
0190 
0191   ///Continuous processor number [0...12], to be used in as array index,
0192   unsigned int getProcIndx(unsigned int iProcessor, l1t::tftype mtfType) const {
0193     return ((mtfType - l1t::tftype::omtf_neg) * rawParams.nProcessors() + iProcessor);
0194   };
0195 
0196   //TODO implement more efficient
0197   bool isBendingLayer(unsigned int iLayer) const override { return getBendingLayers().count(iLayer); }
0198 
0199   ///pattern pt range in Gev
0200   struct PatternPt {
0201     double ptFrom = 0;
0202     double ptTo = 0;
0203     int charge = 0;
0204   };
0205 
0206   PatternPt getPatternPtRange(unsigned int patNum) const;
0207 
0208   //call it when the patterns are read directly from the xml, without using the LUTs
0209   void setPatternPtRange(const std::vector<PatternPt>& patternPts) { this->patternPts = patternPts; }
0210 
0211   ///charge: -1 - negative, +1 - positive
0212   unsigned int getPatternNum(double pt, int charge) const;
0213 
0214   static const unsigned int patternsInGroup = 4;
0215 
0216   //takes the groups from the key, it should be set during xml reading, or creating the goldenPats
0217   template <class GoldenPatternType>
0218   vector2D getPatternGroups(const std::vector<std::unique_ptr<GoldenPatternType> >& goldenPats) const {
0219     //unsigned int mergedCnt = 4;
0220     vector2D mergedPatterns;
0221     for (unsigned int iPat = 0; iPat < goldenPats.size(); iPat++) {
0222       unsigned int group = goldenPats.at(iPat)->key().theGroup;
0223 
0224       if (mergedPatterns.size() == group) {
0225         mergedPatterns.push_back(vector1D());
0226       }
0227 
0228       if (group < mergedPatterns.size()) {
0229         //theIndexInGroup starts from 1, as in xml
0230         if (mergedPatterns[group].size() == (goldenPats.at(iPat)->key().theIndexInGroup - 1))
0231           mergedPatterns[group].push_back(iPat);
0232         else
0233           return mergedPatterns;  //TODO should throw error
0234       } else
0235         return mergedPatterns;  //TODO should throw error
0236     }
0237     return mergedPatterns;
0238   }
0239 
0240   /**configuration from the edm::ParameterSet
0241    * the parameters are set (i.e. overwritten) only if their exist in the edmParameterSet
0242    */
0243   void configureFromEdmParameterSet(const edm::ParameterSet& edmParameterSet) override;
0244 
0245   int getGoldenPatternResultFinalizeFunction() const { return goldenPatternResultFinalizeFunction; }
0246 
0247   void setGoldenPatternResultFinalizeFunction(int goldenPatternResultFinalizeFunction = 0) {
0248     this->goldenPatternResultFinalizeFunction = goldenPatternResultFinalizeFunction;
0249   }
0250 
0251   friend std::ostream& operator<<(std::ostream& out, const OMTFConfiguration& aConfig);
0252 
0253   bool isNoHitValueInPdf() const { return noHitValueInPdf; }
0254 
0255   void setNoHitValueInPdf(bool noHitValueInPdf = false) { this->noHitValueInPdf = noHitValueInPdf; }
0256 
0257   int getSorterType() const { return sorterType; }
0258 
0259   void setSorterType(int sorterType = 0) { this->sorterType = sorterType; }
0260 
0261   const std::string& getGhostBusterType() const { return ghostBusterType; }
0262 
0263   void setGhostBusterType(const std::string& ghostBusterType = "") { this->ghostBusterType = ghostBusterType; }
0264 
0265   bool usePhiBExtrapolationMB1() const { return this->usePhiBExtrapolationFromMB1_; }
0266 
0267   bool usePhiBExtrapolationMB2() const { return this->usePhiBExtrapolationFromMB2_; }
0268 
0269   int getDtRefHitMinQuality() const { return dtRefHitMinQuality; }
0270 
0271   void setDtRefHitMinQuality(int dtRefHitMinQuality = 2) { this->dtRefHitMinQuality = dtRefHitMinQuality; }
0272 
0273   bool getDumpResultToXML() const { return dumpResultToXML; }
0274 
0275   void printConfig() const;
0276 
0277   bool useEndcapStubsRInExtr() const { return useEndcapStubsRInExtr_; }
0278 
0279   bool useStubQualInExtr() const { return useStubQualInExtr_; }
0280 
0281   //[cm]
0282   int minCSCStubRME12() const { return minCSCStubRME12_; }
0283   //[cm], othere CSC station than ME1/2
0284   int minCscStubR() const { return minCSCStubR_; }
0285 
0286   bool cleanStubs() const { return cleanStubs_; }
0287 
0288 private:
0289   L1TMuonOverlapParams rawParams;
0290 
0291   std::map<int, int> hwToLogicLayer;
0292   std::map<int, int> logicToHwLayer;
0293   std::map<int, int> logicToLogic;
0294   std::set<int> bendingLayers;
0295   std::vector<int> refToLogicNumber;
0296 
0297   ///Starting and final sectors connected to
0298   ///processors.
0299   ///Index: processor number
0300   std::vector<unsigned int> barrelMin;
0301   std::vector<unsigned int> barrelMax;
0302   std::vector<unsigned int> endcap10DegMin;
0303   std::vector<unsigned int> endcap10DegMax;
0304   std::vector<unsigned int> endcap20DegMin;
0305   std::vector<unsigned int> endcap20DegMax;
0306 
0307   ///Starting iPhi for each processor and each referecne layer
0308   ///Global phi scale is used
0309   ///First index: processor number
0310   ///Second index: referecne layer number
0311   std::vector<std::vector<int> > processorPhiVsRefLayer;
0312 
0313   ///Begin and end local phi for each logis region
0314   ///First index: input number
0315   ///Second index: reference layer number
0316   ///Third index: region
0317   ///pair.first: starting phi of region (inclusive)
0318   ///pair.second: ending phi of region (inclusive)
0319   std::vector<std::vector<std::vector<std::pair<int, int> > > > regionPhisVsRefLayerVsInput;
0320 
0321   ///Vector with definitions of reference hits
0322   ///Vector has fixed size of nRefHits
0323   ///Order of elements defines priority order
0324   ///First index: processor number (0-5)
0325   ///Second index: ref hit number (0-127)
0326   std::vector<std::vector<RefHitDef> > refHitsDefs;
0327 
0328   ///Map of connections
0329   vector3D_pair connections;
0330 
0331   ///4D matrices used during creation of the connections tables.
0332   vector4D measurements4D;
0333   vector4D measurements4Dref;
0334 
0335   std::vector<PatternPt> patternPts;
0336 
0337   int pdfMaxVal = 0;
0338   unsigned int pdfBins = 0;
0339 
0340   int goldenPatternResultFinalizeFunction = 0;
0341 
0342   //likelihood of "no hit" in the pdf
0343   bool noHitValueInPdf = false;
0344 
0345   int sorterType = 0;
0346 
0347   std::string ghostBusterType = "";
0348 
0349   //if true, in the OMTFProcessor<GoldenPatternType>::processInput the phiB extrapolation is used for the refHit of the MB1, i.e. logicLayer 0 and 1
0350   bool usePhiBExtrapolationFromMB1_ = false;
0351 
0352   //if true, in the OMTFProcessor<GoldenPatternType>::processInput the phiB extrapolation is used for the refHit of the MB2, i.e. logicLayer 2 and 3
0353   bool usePhiBExtrapolationFromMB2_ = false;
0354 
0355   bool useStubQualInExtr_ = false;
0356 
0357   bool useEndcapStubsRInExtr_ = false;
0358 
0359   //min quality of the DT phi hit used as the reference hit
0360   //Remember that it is on the top of the minDtPhiQuality
0361   int dtRefHitMinQuality = 2;
0362 
0363   int minCSCStubRME12_ = 0;
0364   int minCSCStubR_ = 0;
0365 
0366   bool dumpResultToXML = false;
0367 
0368   bool cleanStubs_ = false;
0369 };
0370 
0371 #endif