Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-13 03:23:53

0001 /*
0002  * OMTFProcessor.cpp
0003  *
0004  *  Created on: Oct 7, 2017
0005  *      Author: kbunkow
0006  */
0007 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/OMTFProcessor.h"
0008 #include "L1Trigger/L1TMuonOverlapPhase1/interface/MuonStub.h"
0009 #include "L1Trigger/L1TMuonOverlapPhase1/interface/MuonStubsInput.h"
0010 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/GhostBuster.h"
0011 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/GhostBusterPreferRefDt.h"
0012 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/GoldenPatternWithStat.h"
0013 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/IOMTFEmulationObserver.h"
0014 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/OMTFinput.h"
0015 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Omtf/OMTFSorter.h"
0016 #include "L1Trigger/L1TMuonOverlapPhase1/interface/StubResult.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 #include <bitset>
0021 #include <cmath>
0022 #include <cstdlib>
0023 #include <iostream>
0024 #include <iomanip>
0025 #include <map>
0026 #include <string>
0027 #include <vector>
0028 
0029 #include <boost/property_tree/ptree.hpp>
0030 #include <boost/property_tree/xml_parser.hpp>
0031 
0032 ///////////////////////////////////////////////
0033 ///////////////////////////////////////////////
0034 template <class GoldenPatternType>
0035 OMTFProcessor<GoldenPatternType>::OMTFProcessor(OMTFConfiguration* omtfConfig,
0036                                                 const edm::ParameterSet& edmCfg,
0037                                                 edm::EventSetup const& evSetup,
0038                                                 const L1TMuonOverlapParams* omtfPatterns)
0039     : ProcessorBase<GoldenPatternType>(omtfConfig, omtfPatterns) {
0040   init(edmCfg, evSetup);
0041 };
0042 
0043 template <class GoldenPatternType>
0044 OMTFProcessor<GoldenPatternType>::OMTFProcessor(OMTFConfiguration* omtfConfig,
0045                                                 const edm::ParameterSet& edmCfg,
0046                                                 edm::EventSetup const& evSetup,
0047                                                 GoldenPatternVec<GoldenPatternType>&& gps)
0048     : ProcessorBase<GoldenPatternType>(omtfConfig, std::forward<GoldenPatternVec<GoldenPatternType> >(gps)) {
0049   init(edmCfg, evSetup);
0050 };
0051 
0052 template <class GoldenPatternType>
0053 OMTFProcessor<GoldenPatternType>::~OMTFProcessor() {
0054   if (useFloatingPointExtrapolation)
0055     saveExtrapolFactors();
0056 }
0057 
0058 template <class GoldenPatternType>
0059 void OMTFProcessor<GoldenPatternType>::init(const edm::ParameterSet& edmCfg, edm::EventSetup const& evSetup) {
0060   setSorter(new OMTFSorter<GoldenPatternType>(this->myOmtfConfig->getSorterType()));
0061   //initialize with the default sorter
0062 
0063   if (this->myOmtfConfig->getGhostBusterType() == "GhostBusterPreferRefDt" ||
0064       this->myOmtfConfig->getGhostBusterType() == "byLLH" || this->myOmtfConfig->getGhostBusterType() == "byFPLLH" ||
0065       this->myOmtfConfig->getGhostBusterType() == "byRefLayer") {
0066     setGhostBuster(new GhostBusterPreferRefDt(this->myOmtfConfig));
0067     edm::LogVerbatim("OMTFReconstruction") << "setting " << this->myOmtfConfig->getGhostBusterType() << std::endl;
0068   } else {
0069     setGhostBuster(new GhostBuster(this->myOmtfConfig));  //initialize with the default sorter
0070     edm::LogVerbatim("OMTFReconstruction") << "setting GhostBuster" << std::endl;
0071   }
0072 
0073   edm::LogVerbatim("OMTFReconstruction") << "fwVersion 0x" << hex << this->myOmtfConfig->fwVersion() << std::endl;
0074 
0075   useStubQualInExtr = this->myOmtfConfig->useStubQualInExtr();
0076   useEndcapStubsRInExtr = this->myOmtfConfig->useEndcapStubsRInExtr();
0077 
0078   if (edmCfg.exists("useFloatingPointExtrapolation"))
0079     useFloatingPointExtrapolation = edmCfg.getParameter<bool>("useFloatingPointExtrapolation");
0080 
0081   std::string extrapolFactorsFilename;
0082   if (edmCfg.exists("extrapolFactorsFilename")) {
0083     extrapolFactorsFilename = edmCfg.getParameter<edm::FileInPath>("extrapolFactorsFilename").fullPath();
0084   }
0085 
0086   if (this->myOmtfConfig->usePhiBExtrapolationMB1() || this->myOmtfConfig->usePhiBExtrapolationMB2()) {
0087     extrapolFactors.resize(2, std::vector<std::map<int, double> >(this->myOmtfConfig->nLayers()));
0088     extrapolFactorsNorm.resize(2, std::vector<std::map<int, int> >(this->myOmtfConfig->nLayers()));
0089 
0090     //when useFloatingPointExtrapolation is true the extrapolFactors are not used,
0091     //all calculations are done in the extrapolateDtPhiBFloatPoint
0092     if (!extrapolFactorsFilename.empty() && !useFloatingPointExtrapolation)
0093       loadExtrapolFactors(extrapolFactorsFilename);
0094   }
0095 
0096   edm::LogVerbatim("OMTFReconstruction") << "useFloatingPointExtrapolation " << useFloatingPointExtrapolation
0097                                          << std::endl;
0098 }
0099 
0100 template <class GoldenPatternType>
0101 std::vector<l1t::RegionalMuonCand> OMTFProcessor<GoldenPatternType>::getFinalcandidates(unsigned int iProcessor,
0102                                                                                         l1t::tftype mtfType,
0103                                                                                         const AlgoMuons& algoCands) {
0104   std::vector<l1t::RegionalMuonCand> result;
0105 
0106   for (auto& myCand : algoCands) {
0107     l1t::RegionalMuonCand candidate;
0108 
0109     //the charge is only for the constrained measurement. The constrained measurement is always defined for a valid candidate
0110     if (ptAssignment) {
0111       if (myCand->getPdfSumConstr() > 0 && myCand->getFiredLayerCntConstr() >= 3)
0112         candidate.setHwPt(myCand->getPtNNConstr());
0113       else if (myCand->getPtUnconstr() > 0)
0114         candidate.setHwPt(1);
0115       else
0116         candidate.setHwPt(0);
0117 
0118       candidate.setHwSign(myCand->getChargeNNConstr() < 0 ? 1 : 0);
0119     } else {
0120       if (myCand->getPdfSumConstr() > 0 && myCand->getFiredLayerCntConstr() >= 3)
0121         candidate.setHwPt(myCand->getPtConstr());
0122       else if (myCand->getPtUnconstr() > 0)
0123         //if myCand->getPdfSumConstr() == 0, the myCand->getPtConstr() might not be 0, see the end of GhostBusterPreferRefDt::select
0124         //but 0 means empty candidate, 1 means pt=0, therefore here we set HwPt to 1, as the PtUnconstr > 0
0125         candidate.setHwPt(1);
0126       else
0127         candidate.setHwPt(0);
0128 
0129       candidate.setHwSign(myCand->getChargeConstr() < 0 ? 1 : 0);
0130     }
0131 
0132     if (mtfType == l1t::omtf_pos)
0133       candidate.setHwEta(myCand->getEtaHw());
0134     else
0135       candidate.setHwEta((-1) * myCand->getEtaHw());
0136 
0137     int phiValue = myCand->getPhi();
0138     if (phiValue >= int(this->myOmtfConfig->nPhiBins()))
0139       phiValue -= this->myOmtfConfig->nPhiBins();
0140     phiValue = this->myOmtfConfig->procPhiToGmtPhi(phiValue);
0141     candidate.setHwPhi(phiValue);
0142 
0143     candidate.setHwSignValid(1);
0144 
0145     if (myCand->getPtUnconstr() >= 0) {  //empty PtUnconstrained is -1, maybe should be corrected on the source
0146       //the upt has different hardware scale than the pt, the upt unit is 1 GeV
0147       candidate.setHwPtUnconstrained(myCand->getPtUnconstr());
0148     } else
0149       candidate.setHwPtUnconstrained(0);
0150 
0151     unsigned int quality = 12;
0152     if (this->myOmtfConfig->fwVersion() <= 6)
0153       quality = checkHitPatternValidity(myCand->getFiredLayerBits()) ? 0 | (1 << 2) | (1 << 3) : 0 | (1 << 2);  //12 : 4
0154 
0155     if (abs(myCand->getEtaHw()) == 115 &&  //115 is eta 1.25                        rrrrrrrrccccdddddd
0156         (static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000001110000000").to_ulong() ||
0157          static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000001110000000").to_ulong() ||
0158          static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000000110000000").to_ulong() ||
0159          static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000001100000000").to_ulong() ||
0160          static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000001010000000").to_ulong())) {
0161       if (this->myOmtfConfig->fwVersion() <= 6)
0162         quality = 4;
0163       else
0164         quality = 1;
0165     }
0166 
0167     if (this->myOmtfConfig->fwVersion() >= 5 && this->myOmtfConfig->fwVersion() <= 6) {
0168       if (static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000000011").to_ulong() ||
0169           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000000011").to_ulong() ||
0170           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000000011").to_ulong() ||
0171           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000000011").to_ulong() ||
0172           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000100000000000011").to_ulong() ||
0173           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000000000000011").to_ulong() ||
0174           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000000000000011").to_ulong() ||
0175           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000000000000011").to_ulong() ||
0176           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000001100").to_ulong() ||
0177           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000001100").to_ulong() ||
0178           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000001100").to_ulong() ||
0179           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000001100").to_ulong() ||
0180           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000100000000001100").to_ulong() ||
0181           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000000000001100").to_ulong() ||
0182           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000000000001100").to_ulong() ||
0183           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000000000001100").to_ulong() ||
0184           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000110000").to_ulong() ||
0185           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000110000").to_ulong() ||
0186           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000110000").to_ulong() ||
0187           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000110000").to_ulong() ||
0188           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000100000000110000").to_ulong() ||
0189           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000000000110000").to_ulong() ||
0190           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000000000110000").to_ulong() ||
0191           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000000000110000").to_ulong())
0192         quality = 1;
0193     } else if (this->myOmtfConfig->fwVersion() >= 8) {  //TODO fix the fwVersion     rrrrrrrrccccdddddd
0194       if (static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110000000011").to_ulong() ||
0195           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000000011").to_ulong() ||
0196           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000000011").to_ulong() ||
0197           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110000000001").to_ulong() ||
0198 
0199           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000001100").to_ulong() ||
0200           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011000000001100").to_ulong() ||
0201           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000001100").to_ulong() ||
0202           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011000000000100").to_ulong() ||
0203 
0204           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000011000000001").to_ulong() ||
0205           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000010000000001").to_ulong())
0206         quality = 1;
0207       else if (
0208           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000000101").to_ulong() ||
0209           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010001000001").to_ulong() ||
0210           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000011000000001").to_ulong() ||
0211           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000011000000011").to_ulong() ||
0212           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000011100000001").to_ulong() ||
0213           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000000011").to_ulong() ||
0214           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100001000100").to_ulong() ||
0215           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100100000001").to_ulong() ||
0216           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110100000001").to_ulong() ||
0217           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000111000000000").to_ulong() ||
0218           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000111000000001").to_ulong() ||
0219           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000111000000011").to_ulong() ||
0220           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000001000100").to_ulong() ||
0221           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001010000000001").to_ulong() ||
0222           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001010000000011").to_ulong() ||
0223           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001010000000100").to_ulong() ||
0224           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000000001").to_ulong() ||
0225           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000000100").to_ulong() ||
0226           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000000111").to_ulong() ||
0227           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100001000000").to_ulong() ||
0228           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001110000000100").to_ulong() ||
0229           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001110000000101").to_ulong() ||
0230           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000000101").to_ulong() ||
0231           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010010000000001").to_ulong() ||
0232           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010010000000100").to_ulong() ||
0233           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010010000000101").to_ulong() ||
0234           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010100000000001").to_ulong() ||
0235           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010100000000101").to_ulong() ||
0236           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011110000000100").to_ulong() ||
0237           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011110000000101").to_ulong() ||
0238           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000101000000010101").to_ulong() ||
0239           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000010000000001").to_ulong() ||
0240           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000011000000000").to_ulong() ||
0241           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000011000000001").to_ulong() ||
0242           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000100000000001").to_ulong() ||
0243           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000110000000000").to_ulong() ||
0244           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001001000000000100").to_ulong() ||
0245           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001001100000000100").to_ulong() ||
0246           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001010000000000100").to_ulong() ||
0247           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000000010000001").to_ulong() ||
0248           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000000011000100").to_ulong() ||
0249           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000010000000001").to_ulong() ||
0250           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000100000000001").to_ulong() ||
0251           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000011000000000").to_ulong() ||
0252           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110000000001").to_ulong() ||
0253           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010000000011").to_ulong() ||
0254           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110000000011").to_ulong() ||
0255           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011000000001100").to_ulong() ||
0256           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000011000000000100").to_ulong() ||
0257           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000010010000001").to_ulong() ||
0258           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010000000001100").to_ulong() ||
0259           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001001000001000100").to_ulong() ||
0260           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000000101").to_ulong() ||
0261           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000100000000101").to_ulong() ||
0262           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000000011").to_ulong() ||
0263           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001110000000111").to_ulong() ||
0264           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000110001000001").to_ulong() ||
0265           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001110000000011").to_ulong() ||
0266           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000000001000100").to_ulong() ||
0267           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000110001000001").to_ulong() ||
0268           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000000101").to_ulong() ||
0269           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001010000001000000").to_ulong() ||
0270           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001100000001000000").to_ulong() ||
0271           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("100000010000000001").to_ulong() ||
0272           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000010010000000").to_ulong() ||
0273           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000010100000001100").to_ulong() ||
0274           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000110000000011").to_ulong() ||
0275           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001000000001100").to_ulong() ||
0276           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000000000000111101").to_ulong() ||
0277           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000001100000110001").to_ulong() ||
0278           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000100000000010100").to_ulong() ||
0279           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000100000000011").to_ulong() ||
0280           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("001000110000000001").to_ulong() ||
0281           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("010000100010000001").to_ulong() ||
0282           static_cast<unsigned int>(myCand->getFiredLayerBits()) == std::bitset<18>("000100000000110000").to_ulong())
0283         quality = 8;
0284     }  //  if (abs(myCand->getEta()) == 121) quality = 4;
0285     if (abs(myCand->getEtaHw()) >= 121)
0286       quality = 0;  // changed from 4 on request from HI
0287 
0288     candidate.setHwQual(quality);
0289 
0290     std::map<int, int> trackAddr;
0291     trackAddr[0] = myCand->getFiredLayerBits();
0292     //TODO in the hardware, the uPt is sent to the uGMT at the trackAddr = (uPt << 18) + trackAddr;
0293     //check if it matters if it needs to be here as well
0294     trackAddr[1] = myCand->getRefLayer();
0295     trackAddr[2] = myCand->getDisc();
0296     if (candidate.hwPt() > 0 || candidate.hwPtUnconstrained() > 0) {
0297       candidate.setTrackAddress(trackAddr);
0298       candidate.setTFIdentifiers(iProcessor, mtfType);
0299       result.push_back(candidate);
0300     }
0301   }
0302   return result;
0303 }
0304 ///////////////////////////////////////////////////////
0305 ///////////////////////////////////////////////////////
0306 
0307 ///////////////////////////////////////////////////////
0308 ///////////////////////////////////////////////////////
0309 template <class GoldenPatternType>
0310 bool OMTFProcessor<GoldenPatternType>::checkHitPatternValidity(unsigned int hits) {
0311   ///FIXME: read the list from configuration so this can be controlled at runtime.
0312   std::vector<unsigned int> badPatterns = {
0313       99840, 34304, 3075, 36928, 12300, 98816, 98944, 33408, 66688, 66176, 7171, 20528, 33856, 35840, 4156, 34880};
0314 
0315   /*
0316 99840 01100001 1000 000000      011000011000000000
0317 34304 00100001 1000 000000      001000011000000000
0318  3075 00000011 0000 000011      000000110000000011
0319 36928 00100100 0001 000000      001001000001000000
0320 12300 00001100 0000 001100      000011000000001100
0321 98816 01100000 1000 000000      011000001000000000
0322 98944 01100000 1010 000000      011000001010000000
0323 33408 00100000 1010 000000      001000001010000000
0324 66688 01000001 0010 000000      010000010010000000
0325 66176 01000000 1010 000000      010000001010000000
0326  7171 00000111 0000 000011      000001110000000011
0327 20528 00010100 0000 110000      000101000000110000
0328 33856 00100001 0001 000000      001000010001000000
0329 35840 00100011 0000 000000      001000110000000000
0330  4156 00000100 0000 111100      000001000000111100
0331 34880 00100010 0001 000000      001000100001000000
0332    */
0333   for (auto aHitPattern : badPatterns) {
0334     if (hits == aHitPattern)
0335       return false;
0336   }
0337 
0338   return true;
0339 }
0340 ///////////////////////////////////////////////////////
0341 ///////////////////////////////////////////////////////
0342 template <class GoldenPatternType>
0343 AlgoMuons OMTFProcessor<GoldenPatternType>::sortResults(unsigned int iProcessor, l1t::tftype mtfType, int charge) {
0344   unsigned int procIndx = this->myOmtfConfig->getProcIndx(iProcessor, mtfType);
0345   return sorter->sortResults(procIndx, this->getPatterns(), charge);
0346 }
0347 
0348 template <class GoldenPatternType>
0349 int OMTFProcessor<GoldenPatternType>::extrapolateDtPhiBFloatPoint(const int& refLogicLayer,
0350                                                                   const int& refPhi,
0351                                                                   const int& refPhiB,
0352                                                                   unsigned int targetLayer,
0353                                                                   const int& targetStubPhi,
0354                                                                   const int& targetStubQuality,
0355                                                                   const int& targetStubEta,
0356                                                                   const int& targetStubR,
0357                                                                   const OMTFConfiguration* omtfConfig) {
0358   LogTrace("l1tOmtfEventPrint") << "\n"
0359                                 << __FUNCTION__ << ":" << __LINE__ << " refLogicLayer " << refLogicLayer
0360                                 << " targetLayer " << targetLayer << std::endl;
0361   LogTrace("l1tOmtfEventPrint") << "refPhi " << refPhi << " refPhiB " << refPhiB << " targetStubPhi " << targetStubPhi
0362                                 << " targetStubQuality " << targetStubQuality << std::endl;
0363 
0364   int phiExtr = 0;  //delta phi extrapolated
0365 
0366   float rRefLayer = 431.133;  //[cm], MB1 i.e. refLogicLayer = 0
0367   if (refLogicLayer == 2)
0368     rRefLayer = 512.401;  //MB2
0369   else if (refLogicLayer != 0) {
0370     return 0;
0371     //throw cms::Exception("OMTFProcessor<GoldenPatternType>::extrapolateDtPhiB: wrong refStubLogicLayer " + std::to_string(refLogicLayer) );
0372   }
0373 
0374   int reflLayerIndex = refLogicLayer == 0 ? 0 : 1;
0375 
0376   if (targetLayer == 0 || targetLayer == 2 || targetLayer == 4 || (targetLayer >= 10 && targetLayer <= 14)) {
0377     //all units are cm. Values from the CMS geometry
0378     float rTargetLayer = 512.401;  //MB2
0379 
0380     if (targetLayer == 0)
0381       rTargetLayer = 431.133;  //MB1
0382     else if (targetLayer == 4)
0383       rTargetLayer = 617.946;  //MB3
0384 
0385     else if (targetLayer == 10)
0386       rTargetLayer = 413.675;  //RB1in
0387     else if (targetLayer == 11)
0388       rTargetLayer = 448.675;  //RB1out
0389     else if (targetLayer == 12)
0390       rTargetLayer = 494.975;  //RB2in
0391     else if (targetLayer == 13)
0392       rTargetLayer = 529.975;  //RB2out
0393     else if (targetLayer == 14)
0394       rTargetLayer = 602.150;  //RB3
0395 
0396     if (useStubQualInExtr) {
0397       if (targetLayer == 0 || targetLayer == 2 || targetLayer == 4) {
0398         if (targetStubQuality == 2 || targetStubQuality == 0)
0399           rTargetLayer = rTargetLayer - 23.5 / 2;  //inner superlayer
0400         else if (targetStubQuality == 3 || targetStubQuality == 1)
0401           rTargetLayer = rTargetLayer + 23.5 / 2;  //outer superlayer
0402       }
0403     }
0404 
0405     float d = rTargetLayer - rRefLayer;
0406     //formula in the form as in the slides explaining the extrapolation algorithm
0407     //float deltaPhiExtr = d/rTargetLayer * refPhiB / omtfConfig->dtPhiBUnitsRad(); //[rad]
0408     //phiExtr = round(deltaPhiExtr / omtfConfig->omtfPhiUnit()); //[halfStrip]
0409 
0410     //formula with approximation, used to calculate extrFactor
0411     float extrFactor = d / rTargetLayer / omtfConfig->dtPhiBUnitsRad() / omtfConfig->omtfPhiUnit();
0412     phiExtr = extrFactor * (float)refPhiB;  //[halfStrip]
0413 
0414     //formula without approximation
0415     float deltaPhiExtr = atan(d / rTargetLayer * tan(refPhiB / omtfConfig->dtPhiBUnitsRad()));  //[rad]
0416     phiExtr = round(deltaPhiExtr / omtfConfig->omtfPhiUnit());                                  //[halfStrip]
0417 
0418     if (useStubQualInExtr & (targetLayer == 0 || targetLayer == 2 || targetLayer == 4)) {
0419       extrapolFactors[reflLayerIndex][targetLayer][targetStubQuality] = extrFactor;
0420       extrapolFactorsNorm[reflLayerIndex][targetLayer][targetStubQuality] = 1;
0421     } else {
0422       extrapolFactors[reflLayerIndex][targetLayer][0] = extrFactor;
0423       extrapolFactorsNorm[reflLayerIndex][targetLayer][0] = 1;
0424     }
0425 
0426     //LogTrace("l1tOmtfEventPrint") <<__FUNCTION__<<":"<<__LINE__<<" deltaPhiExtr "<<deltaPhiExtr<<" phiExtr "<<phiExtr<<std::endl;
0427 
0428     LogTrace("l1tOmtfEventPrint") << "\n"
0429                                   << __FUNCTION__ << ":" << __LINE__ << " refLogicLayer " << refLogicLayer
0430                                   << " targetLayer " << std::setw(2) << targetLayer << " targetStubQuality "
0431                                   << targetStubQuality << " extrFactor " << extrFactor << std::endl;
0432 
0433     LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << ":" << __LINE__ << " refPhiB " << refPhiB << " phiExtr " << phiExtr
0434                                   << std::endl;
0435 
0436   } else if (targetLayer == 1 || targetLayer == 3 || targetLayer == 5) {
0437     int deltaPhi = targetStubPhi - refPhi;  //[halfStrip]
0438 
0439     //deltaPhi is here in phi_b hw scale
0440     deltaPhi = round(deltaPhi * omtfConfig->omtfPhiUnit() * omtfConfig->dtPhiBUnitsRad());
0441 
0442     phiExtr = refPhiB - deltaPhi;  //phiExtr is also in phi_b hw scale
0443     LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << ":" << __LINE__ << " deltaPhi " << deltaPhi << " phiExtr "
0444                                   << phiExtr << std::endl;
0445   } else if ((targetLayer >= 6 && targetLayer <= 9) || (targetLayer >= 15 && targetLayer <= 17)) {
0446     //if true, for the CSC and endcap RPC the R is taken from the hit coordinates
0447 
0448     float rME = targetStubR;
0449     if (!useEndcapStubsRInExtr) {
0450       //all units are cm. This are the average R values for a given chamber (more or less middle of the chamber, but taking into account the OMTF eta range)
0451       if (targetLayer == 6 || targetLayer == 15)  //ME1/3, RE1/3,
0452         rME = 600.;
0453       else if (targetLayer == 7 || targetLayer == 15) {  //ME2/2, RE2/3,
0454         if (refLogicLayer == 0)
0455           rME = 600.;
0456         else
0457           rME = 640.;
0458       } else if (targetLayer == 8 || rME == 16) {  //ME3/2, RE3/3,
0459         if (refLogicLayer == 0)
0460           rME = 620.;
0461         else
0462           rME = 680.;
0463       } else if (targetLayer == 9) {
0464         rME = 460.;  //for the refLogicLayer = 1. refLogicLayer = 2 is impossible
0465       }
0466     }
0467 
0468     float d = rME - rRefLayer;
0469     //formula in the form as in the slides explaining the extrapolation algorithm
0470     //float deltaPhiExtr = d / rME * refPhiB / omtfConfig->dtPhiBUnitsRad();  //[rad]
0471     //phiExtr = round(deltaPhiExtr / omtfConfig->omtfPhiUnit()); //[halfStrip]
0472 
0473     //formula with approximation, used to calculate extrFactor
0474     float extrFactor = d / rME / omtfConfig->dtPhiBUnitsRad() / omtfConfig->omtfPhiUnit();
0475     phiExtr = extrFactor * refPhiB;  //[halfStrip]
0476 
0477     //formula without approximation
0478     float deltaPhiExtr = atan(d / rME * tan(refPhiB / omtfConfig->dtPhiBUnitsRad()));  //[rad]
0479     phiExtr = round(deltaPhiExtr / omtfConfig->omtfPhiUnit());                         //[halfStrip]
0480 
0481     if (useEndcapStubsRInExtr) {
0482       //extrapolFactors[reflLayerIndex][targetLayer][std::abs(targetStubEta)] += extrFactor;
0483       //extrapolFactorsNorm[reflLayerIndex][targetLayer][std::abs(targetStubEta)]++;
0484       extrapolFactors[reflLayerIndex][targetLayer][std::abs(rME)] += extrFactor;
0485       extrapolFactorsNorm[reflLayerIndex][targetLayer][std::abs(rME)]++;
0486       //extrapolFactors[reflLayerIndex][targetLayer][0] += extrFactor;
0487       //extrapolFactorsNorm[reflLayerIndex][targetLayer][0]++;
0488     } else {
0489       extrapolFactors[reflLayerIndex][targetLayer][0] = extrFactor;
0490       extrapolFactorsNorm[reflLayerIndex][targetLayer][0] = 1;
0491     }
0492     LogTrace("l1tOmtfEventPrint") << "\n"
0493                                   << __FUNCTION__ << ":" << __LINE__ << " refLogicLayer " << refLogicLayer
0494                                   << " targetLayer " << std::setw(2) << targetLayer << " targetStubR " << targetStubR
0495                                   << " targetStubEta " << targetStubEta << " extrFactor "
0496                                   << " rRefLayer " << rRefLayer << " d " << d << " deltaPhiExtr " << deltaPhiExtr
0497                                   << " phiExtr " << phiExtr << std::endl;
0498   }
0499   //TODO restrict the range of the phiExtr and refPhiB !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
0500 
0501   return phiExtr;
0502 }
0503 
0504 template <class GoldenPatternType>
0505 int OMTFProcessor<GoldenPatternType>::extrapolateDtPhiBFixedPoint(const int& refLogicLayer,
0506                                                                   const int& refPhi,
0507                                                                   const int& refPhiB,
0508                                                                   unsigned int targetLayer,
0509                                                                   const int& targetStubPhi,
0510                                                                   const int& targetStubQuality,
0511                                                                   const int& targetStubEta,
0512                                                                   const int& targetStubR,
0513                                                                   const OMTFConfiguration* omtfConfig) {
0514   int phiExtr = 0;  //delta phi extrapolated
0515 
0516   int reflLayerIndex = refLogicLayer == 0 ? 0 : 1;
0517   int extrFactor = 0;
0518 
0519   if (targetLayer == 0 || targetLayer == 2 || targetLayer == 4) {
0520     if (useStubQualInExtr)
0521       extrFactor = extrapolFactors[reflLayerIndex][targetLayer][targetStubQuality];
0522     else
0523       extrFactor = extrapolFactors[reflLayerIndex][targetLayer][0];
0524   } else if (targetLayer == 1 || targetLayer == 3 || targetLayer == 5) {
0525     int deltaPhi = targetStubPhi - refPhi;  //here targetStubPhi is phi, not phiB
0526 
0527     int scaleFactor = this->myOmtfConfig->omtfPhiUnit() * this->myOmtfConfig->dtPhiBUnitsRad() * 512;
0528     //= 305 for phase-1, 512 is multiplier so that scaleFactor is non-zero integer
0529 
0530     deltaPhi = (deltaPhi * scaleFactor) / 512;  //here deltaPhi is converted to the phi_b hw scale
0531 
0532     phiExtr = refPhiB - deltaPhi;  //phiExtr is also in phi_b hw scale
0533     //LogTrace("l1tOmtfEventPrint") <<__FUNCTION__<<":"<<__LINE__<<" deltaPhi "<<deltaPhi<<" phiExtr "<<phiExtr<<std::endl;
0534 
0535   } else if (targetLayer >= 10 && targetLayer <= 14) {
0536     extrFactor = extrapolFactors[reflLayerIndex][targetLayer][0];
0537   } else if ((targetLayer >= 6 && targetLayer <= 9) || (targetLayer >= 15 && targetLayer <= 17)) {
0538     if (useEndcapStubsRInExtr) {
0539       //if given abs(targetStubEta) value is not present in the map, it is added with default value of 0
0540       //so it should be good. The only problem is that the map can grow...
0541       //TODO change to targetStubR when it is implemented in the FW
0542       extrFactor = extrapolFactors[reflLayerIndex][targetLayer][abs(targetStubEta)];
0543       //extrFactor = extrapolFactors[reflLayerIndex][targetLayer][abs(targetStubR)];
0544     } else {
0545       extrFactor = extrapolFactors[reflLayerIndex][targetLayer][0];
0546     }
0547   }
0548 
0549   if (this->myOmtfConfig->isBendingLayer(targetLayer) == false) {
0550     phiExtr = extrFactor * refPhiB / extrapolMultiplier;
0551   }
0552 
0553   LogTrace("l1tOmtfEventPrint") << "\n"
0554                                 << __FUNCTION__ << ":" << __LINE__ << " refLogicLayer " << refLogicLayer
0555                                 << " targetLayer " << targetLayer << std::endl;
0556   LogTrace("l1tOmtfEventPrint") << "refPhi " << refPhi << " refPhiB " << refPhiB << " targetStubPhi " << targetStubPhi
0557                                 << " targetStubQuality " << targetStubQuality << " targetStubEta " << targetStubEta
0558                                 << " extrFactor " << extrFactor << " phiExtr " << phiExtr << std::endl;
0559 
0560   return phiExtr;
0561 }
0562 
0563 template <class GoldenPatternType>
0564 int OMTFProcessor<GoldenPatternType>::extrapolateDtPhiB(const MuonStubPtr& refStub,
0565                                                         const MuonStubPtr& targetStub,
0566                                                         unsigned int targetLayer,
0567                                                         const OMTFConfiguration* omtfConfig) {
0568   if (useFloatingPointExtrapolation)
0569     return OMTFProcessor<GoldenPatternType>::extrapolateDtPhiBFloatPoint(refStub->logicLayer,
0570                                                                          refStub->phiHw,
0571                                                                          refStub->phiBHw,
0572                                                                          targetLayer,
0573                                                                          targetStub->phiHw,
0574                                                                          targetStub->qualityHw,
0575                                                                          targetStub->etaHw,
0576                                                                          targetStub->r,
0577                                                                          omtfConfig);
0578   return OMTFProcessor<GoldenPatternType>::extrapolateDtPhiBFixedPoint(refStub->logicLayer,
0579                                                                        refStub->phiHw,
0580                                                                        refStub->phiBHw,
0581                                                                        targetLayer,
0582                                                                        targetStub->phiHw,
0583                                                                        targetStub->qualityHw,
0584                                                                        targetStub->etaHw,
0585                                                                        targetStub->r,
0586                                                                        omtfConfig);
0587 }
0588 ///////////////////////////////////////////////
0589 ///////////////////////////////////////////////
0590 //const std::vector<OMTFProcessor::resultsMap> &
0591 template <class GoldenPatternType>
0592 void OMTFProcessor<GoldenPatternType>::processInput(unsigned int iProcessor,
0593                                                     l1t::tftype mtfType,
0594                                                     const OMTFinput& aInput,
0595                                                     std::vector<std::unique_ptr<IOMTFEmulationObserver> >& observers) {
0596   unsigned int procIndx = this->myOmtfConfig->getProcIndx(iProcessor, mtfType);
0597   for (auto& itGP : this->theGPs) {
0598     for (auto& result : itGP->getResults()[procIndx]) {
0599       result.reset();
0600     }
0601   }
0602 
0603   LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << "\n"
0604                                 << __LINE__ << " iProcessor " << iProcessor << " mtfType " << mtfType << " procIndx "
0605                                 << procIndx << " ----------------------" << std::endl;
0606   //////////////////////////////////////
0607   //////////////////////////////////////
0608   std::vector<const RefHitDef*> refHitDefs;
0609 
0610   {
0611     auto refHitsBits = aInput.getRefHits(iProcessor);
0612     if (refHitsBits.none())
0613       return;  // myResults;
0614 
0615     //loop over all possible refHits, e.g. 128
0616     for (unsigned int iRefHit = 0; iRefHit < this->myOmtfConfig->nRefHits(); ++iRefHit) {
0617       if (!refHitsBits[iRefHit])
0618         continue;
0619 
0620       refHitDefs.push_back(&(this->myOmtfConfig->getRefHitsDefs()[iProcessor][iRefHit]));
0621 
0622       if (refHitDefs.size() == this->myOmtfConfig->nTestRefHits())
0623         break;
0624     }
0625   }
0626 
0627   boost::property_tree::ptree procDataTree;
0628   LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << " " << __LINE__;
0629   for (unsigned int iLayer = 0; iLayer < this->myOmtfConfig->nLayers(); ++iLayer) {
0630     //debug
0631     /*for(auto& h : layerHits) {
0632       if(h != 5400)
0633         LogTrace("l1tOmtfEventPrint")<<__FUNCTION__<<" "<<__LINE__<<" iLayer "<<iLayer<<" layerHit "<<h<<std::endl;
0634     }*/
0635 
0636     for (unsigned int iRefHit = 0; iRefHit < refHitDefs.size(); iRefHit++) {
0637       const RefHitDef& aRefHitDef = *(refHitDefs[iRefHit]);
0638 
0639       unsigned int refLayerLogicNum = this->myOmtfConfig->getRefToLogicNumber()[aRefHitDef.iRefLayer];
0640       const MuonStubPtr refStub = aInput.getMuonStub(refLayerLogicNum, aRefHitDef.iInput);
0641       //int etaRef = refStub->etaHw;
0642 
0643       unsigned int iRegion = aRefHitDef.iRegion;
0644 
0645       MuonStubPtrs1D restrictedLayerStubs = this->restrictInput(iProcessor, iRegion, iLayer, aInput);
0646 
0647       //LogTrace("l1tOmtfEventPrint")<<__FUNCTION__<<" "<<__LINE__<<" iLayer "<<iLayer<<" iRefLayer "<<aRefHitDef.iRefLayer<<std::endl;
0648       //LogTrace("l1tOmtfEventPrint")<<"iLayer "<<iLayer<<" iRefHit "<<iRefHit;
0649       //LogTrace("l1tOmtfEventPrint")<<" nTestedRefHits "<<nTestedRefHits<<" aRefHitDef "<<aRefHitDef<<std::endl;
0650 
0651       std::vector<int> extrapolatedPhi(restrictedLayerStubs.size(), 0);
0652 
0653       //TODO make sure the that the iRefLayer numbers used here corresponds to this in the hwToLogicLayer_0x000X.xml
0654       if ((this->myOmtfConfig->usePhiBExtrapolationMB1() && aRefHitDef.iRefLayer == 0) ||
0655           (this->myOmtfConfig->usePhiBExtrapolationMB2() && aRefHitDef.iRefLayer == 2)) {
0656         if ((iLayer != refLayerLogicNum) && (iLayer != refLayerLogicNum + 1)) {
0657           unsigned int iStub = 0;
0658           for (auto& targetStub : restrictedLayerStubs) {
0659             if (targetStub) {
0660               extrapolatedPhi[iStub] = extrapolateDtPhiB(refStub, targetStub, iLayer, this->myOmtfConfig);
0661 
0662               LogTrace("l1tOmtfEventPrint")
0663                   << "\n"
0664                   << __FUNCTION__ << ":" << __LINE__ << " extrapolating from layer " << refLayerLogicNum
0665                   << " - iRefLayer " << aRefHitDef.iRefLayer << " to layer " << iLayer << " stub " << targetStub
0666                   << " value " << extrapolatedPhi[iStub] << std::endl;
0667 
0668               if (this->myOmtfConfig->getDumpResultToXML()) {
0669                 auto& extrapolatedPhiTree = procDataTree.add_child("extrapolatedPhi", boost::property_tree::ptree());
0670                 extrapolatedPhiTree.add("<xmlattr>.refLayer", refLayerLogicNum);
0671                 extrapolatedPhiTree.add("<xmlattr>.layer", iLayer);
0672                 extrapolatedPhiTree.add("<xmlattr>.refPhiBHw", refStub->phiBHw);
0673                 extrapolatedPhiTree.add("<xmlattr>.iStub", iStub);
0674                 extrapolatedPhiTree.add("<xmlattr>.qualityHw", targetStub->qualityHw);
0675                 extrapolatedPhiTree.add("<xmlattr>.etaHw", targetStub->etaHw);
0676                 extrapolatedPhiTree.add("<xmlattr>.phiExtr", extrapolatedPhi[iStub]);
0677 
0678                 if (this->myOmtfConfig->isBendingLayer(iLayer))
0679                   extrapolatedPhiTree.add("<xmlattr>.dist_phi", targetStub->phiBHw - extrapolatedPhi[iStub]);
0680                 else
0681                   extrapolatedPhiTree.add("<xmlattr>.dist_phi", targetStub->phiHw - extrapolatedPhi[iStub]);
0682               }
0683             }
0684             iStub++;
0685           }
0686         }
0687       }
0688 
0689       for (auto& itGP : this->theGPs) {
0690         if (itGP->key().thePt == 0)  //empty pattern
0691           continue;
0692 
0693         StubResult stubResult =
0694             itGP->process1Layer1RefLayer(aRefHitDef.iRefLayer, iLayer, restrictedLayerStubs, extrapolatedPhi, refStub);
0695 
0696         /* LogTrace("l1tOmtfEventPrint")<<__FUNCTION__<<":"<<__LINE__
0697                                      <<" layerResult: valid"<<stubResult.getValid()
0698                                      <<" pdfVal "<<stubResult.getPdfVal()
0699                                      <<std::endl;*/
0700 
0701         itGP->getResults()[procIndx][iRefHit].setStubResult(iLayer, stubResult);
0702       }
0703     }
0704   }
0705 
0706   for (unsigned int iRefHit = 0; iRefHit < refHitDefs.size(); iRefHit++) {
0707     const RefHitDef& aRefHitDef = *(refHitDefs[iRefHit]);
0708 
0709     unsigned int refLayerLogicNum = this->myOmtfConfig->getRefToLogicNumber()[aRefHitDef.iRefLayer];
0710     const MuonStubPtr refStub = aInput.getMuonStub(refLayerLogicNum, aRefHitDef.iInput);
0711 
0712     int phiRef = refStub->phiHw;
0713     int etaRef = refStub->etaHw;
0714 
0715     //calculating the phiExtrp in the case the RefLayer is MB1, to include it in the  candidate phi of candidate
0716     int phiExtrp = 0;
0717     if ((this->myOmtfConfig->usePhiBExtrapolationMB1() && aRefHitDef.iRefLayer == 0)) {
0718       //||(this->myOmtfConfig->getUsePhiBExtrapolationMB2() && aRefHitDef.iRefLayer == 2) ) {  //the extrapolation from the layer 2 to the layer 2 has no sense, so phiExtrp is 0
0719       LogTrace("l1tOmtfEventPrint") << "\n"
0720                                     << __FUNCTION__ << ":" << __LINE__
0721                                     << "extrapolating ref hit to get the phi of the candidate" << std::endl;
0722       if (useFloatingPointExtrapolation)
0723         phiExtrp = extrapolateDtPhiBFloatPoint(
0724             aRefHitDef.iRefLayer, phiRef, refStub->phiBHw, 2, 0, 6, 0, 0, this->myOmtfConfig);
0725       else
0726         phiExtrp = extrapolateDtPhiBFixedPoint(
0727             aRefHitDef.iRefLayer, phiRef, refStub->phiBHw, 2, 0, 6, 0, 0, this->myOmtfConfig);
0728     }
0729 
0730     for (auto& itGP : this->theGPs) {
0731       if (itGP->key().thePt == 0)  //empty pattern
0732         continue;
0733 
0734       int phiRefSt2 = itGP->propagateRefPhi(phiRef + phiExtrp, etaRef, aRefHitDef.iRefLayer);
0735       itGP->getResults()[procIndx][iRefHit].set(aRefHitDef.iRefLayer, phiRefSt2, etaRef, phiRef);
0736     }
0737   }
0738 
0739   //////////////////////////////////////
0740   //////////////////////////////////////
0741   {
0742     for (auto& itGP : this->theGPs) {
0743       itGP->finalise(procIndx);
0744       //debug
0745       /*for(unsigned int iRefHit = 0; iRefHit < itGP->getResults()[procIndx].size(); ++iRefHit) {
0746         if(itGP->getResults()[procIndx][iRefHit].isValid()) {
0747           LogTrace("l1tOmtfEventPrint")<<__FUNCTION__<<":"<<"__LINE__"<<itGP->getResults()[procIndx][iRefHit]<<std::endl;
0748         }
0749       }*/
0750     }
0751   }
0752 
0753   for (auto& obs : observers)
0754     obs->addProcesorData("extrapolation", procDataTree);
0755 
0756   return;
0757 }
0758 ///////////////////////////////////////////////////////
0759 ///////////////////////////////////////////////////////
0760 
0761 template <class GoldenPatternType>
0762 std::vector<l1t::RegionalMuonCand> OMTFProcessor<GoldenPatternType>::run(
0763     unsigned int iProcessor,
0764     l1t::tftype mtfType,
0765     int bx,
0766     OMTFinputMaker* inputMaker,
0767     std::vector<std::unique_ptr<IOMTFEmulationObserver> >& observers) {
0768   //uncomment if you want to check execution time of each method
0769   //boost::timer::auto_cpu_timer t("%ws wall, %us user in getProcessorCandidates\n");
0770 
0771   for (auto& obs : observers)
0772     obs->observeProcesorBegin(iProcessor, mtfType);
0773 
0774   //input is shared_ptr because the observers may need them after the run() method execution is finished
0775   std::shared_ptr<OMTFinput> input = std::make_shared<OMTFinput>(this->myOmtfConfig);
0776   inputMaker->buildInputForProcessor(input->getMuonStubs(), iProcessor, mtfType, bx, bx, observers);
0777 
0778   if (this->myOmtfConfig->cleanStubs()) {
0779     //this has sense for the pattern generation from the tracks with the secondaries
0780     //if more than one stub is in a given layer, all stubs are removed from this layer
0781     for (unsigned int iLayer = 0; iLayer < input->getMuonStubs().size(); ++iLayer) {
0782       auto& layerStubs = input->getMuonStubs()[iLayer];
0783       int count = std::count_if(layerStubs.begin(), layerStubs.end(), [](auto& ptr) { return ptr != nullptr; });
0784       if (count > 1) {
0785         for (auto& ptr : layerStubs)
0786           ptr.reset();
0787 
0788         LogTrace("OMTFReconstruction") << __FUNCTION__ << ":" << __LINE__ << "cleaning stubs in the layer " << iLayer
0789                                        << " stubs count :" << count << std::endl;
0790       }
0791     }
0792   }
0793 
0794   //LogTrace("l1tOmtfEventPrint")<<"buildInputForProce "; t.report();
0795   processInput(iProcessor, mtfType, *(input.get()), observers);
0796 
0797   //LogTrace("l1tOmtfEventPrint")<<"processInput       "; t.report();
0798   AlgoMuons algoCandidates = sortResults(iProcessor, mtfType);
0799 
0800   if (ptAssignment) {
0801     for (auto& myCand : algoCandidates) {
0802       if (myCand->isValid()) {
0803         auto pts = ptAssignment->getPts(myCand, observers);
0804         /*for (unsigned int i = 0; i < pts.size(); i++) {
0805         trackAddr[10 + i] = this->myOmtfConfig->ptGevToHw(pts[i]);
0806       }*/
0807       }
0808     }
0809   }
0810 
0811   //LogTrace("l1tOmtfEventPrint")<<"sortResults        "; t.report();
0812   // perform GB
0813   //watch out: etaBits2HwEta is used in the ghostBust to convert the AlgoMuons eta, it affect algoCandidates as they are pointers
0814   AlgoMuons gbCandidates = ghostBust(algoCandidates);
0815 
0816   //LogTrace("l1tOmtfEventPrint")<<"ghostBust"; t.report();
0817   // fill RegionalMuonCand colleciton
0818   std::vector<l1t::RegionalMuonCand> candMuons = getFinalcandidates(iProcessor, mtfType, gbCandidates);
0819 
0820   //LogTrace("l1tOmtfEventPrint")<<"getFinalcandidates "; t.report();
0821   //fill outgoing collection
0822   for (auto& candMuon : candMuons) {
0823     candMuon.setHwQual(candMuon.hwQual());
0824   }
0825 
0826   for (auto& obs : observers) {
0827     obs->observeProcesorEmulation(iProcessor, mtfType, input, algoCandidates, gbCandidates, candMuons);
0828   }
0829 
0830   return candMuons;
0831 }
0832 
0833 template <class GoldenPatternType>
0834 void OMTFProcessor<GoldenPatternType>::printInfo() const {
0835   edm::LogVerbatim("OMTFReconstruction") << __PRETTY_FUNCTION__ << std::endl;
0836 
0837   ProcessorBase<GoldenPatternType>::printInfo();
0838 }
0839 
0840 template <class GoldenPatternType>
0841 void OMTFProcessor<GoldenPatternType>::saveExtrapolFactors() {
0842   //if(this->myOmtfConfig->nProcessors() == 3) //phase2
0843   extrapolMultiplier = 512;
0844 
0845   boost::property_tree::ptree tree;
0846   auto& extrFactorsTree = tree.add("ExtrapolationFactors", "");
0847   extrFactorsTree.add("<xmlattr>.multiplier", extrapolMultiplier);
0848 
0849   edm::LogVerbatim("OMTFReconstruction") << "saving extrapolFactors to ExtrapolationFactors.xml" << std::endl;
0850   for (unsigned int iRefLayer = 0; iRefLayer < extrapolFactors.size(); iRefLayer++) {
0851     for (unsigned int iLayer = 0; iLayer < extrapolFactors[iRefLayer].size(); iLayer++) {
0852       edm::LogVerbatim("OMTFReconstruction") << " iRefLayer " << iRefLayer << " iLayer " << iLayer << std::endl;
0853 
0854       auto& layerTree = extrFactorsTree.add_child("Lut", boost::property_tree::ptree());
0855       layerTree.add("<xmlattr>.RefLayer", std::to_string(iRefLayer));
0856       layerTree.add("<xmlattr>.Layer", iLayer);
0857 
0858       if (useStubQualInExtr && (iLayer == 0 || iLayer == 2 || iLayer == 4))
0859         layerTree.add("<xmlattr>.KeyType", "quality");
0860       else if (useEndcapStubsRInExtr && ((iLayer >= 6 && iLayer <= 9) || (iLayer >= 15 && iLayer <= 17)))
0861         layerTree.add("<xmlattr>.KeyType", "eta");
0862       else
0863         layerTree.add("<xmlattr>.KeyType", "none");
0864 
0865       for (auto& extrFactors : extrapolFactors[iRefLayer][iLayer]) {
0866         int norm = 1;
0867         if (!extrapolFactorsNorm[iRefLayer][iLayer].empty())
0868           norm = extrapolFactorsNorm[iRefLayer][iLayer][extrFactors.first];
0869         auto& lutVal = layerTree.add_child("LutVal", boost::property_tree::ptree());
0870         if (useEndcapStubsRInExtr && ((iLayer >= 6 && iLayer <= 9) || (iLayer >= 15 && iLayer <= 17)))
0871           lutVal.add("<xmlattr>.key", extrFactors.first);
0872         else
0873           lutVal.add("<xmlattr>.key", extrFactors.first);
0874 
0875         double value = round(extrapolMultiplier * extrFactors.second / norm);
0876         lutVal.add("<xmlattr>.value", value);
0877 
0878         edm::LogVerbatim("OMTFReconstruction")
0879             << std::setw(4) << " key = " << extrFactors.first << " extrFactors.second " << std::setw(10)
0880             << extrFactors.second << " norm " << std::setw(6) << norm << " value/norm " << std::setw(10)
0881             << extrFactors.second / norm << " value " << value << std::endl;
0882       }
0883     }
0884   }
0885 
0886   boost::property_tree::write_xml("ExtrapolationFactors.xml",
0887                                   tree,
0888                                   std::locale(),
0889                                   boost::property_tree::xml_parser::xml_writer_make_settings<std::string>(' ', 2));
0890 }
0891 
0892 template <class GoldenPatternType>
0893 void OMTFProcessor<GoldenPatternType>::loadExtrapolFactors(const std::string& filename) {
0894   boost::property_tree::ptree tree;
0895 
0896   boost::property_tree::read_xml(filename, tree);
0897 
0898   edm::LogVerbatim("OMTFReconstruction") << "loadExtrapolFactors from file " << filename << std::endl;
0899 
0900   extrapolMultiplier = tree.get<int>("ExtrapolationFactors.<xmlattr>.multiplier");
0901   edm::LogVerbatim("OMTFReconstruction") << "extrapolMultiplier " << extrapolMultiplier << std::endl;
0902 
0903   auto& lutNodes = tree.get_child("ExtrapolationFactors");
0904   for (boost::property_tree::ptree::value_type& lutNode : lutNodes) {
0905     if (lutNode.first == "Lut") {
0906       int iRefLayer = lutNode.second.get<int>("<xmlattr>.RefLayer");
0907       int iLayer = lutNode.second.get<int>("<xmlattr>.Layer");
0908       std::string keyType = lutNode.second.get<std::string>("<xmlattr>.KeyType");
0909 
0910       LogTrace("OMTFReconstruction") << "iRefLayer " << iRefLayer << " iLayer " << iLayer << " keyType " << keyType
0911                                      << std::endl;
0912 
0913       auto& valueNodes = lutNode.second;
0914       for (boost::property_tree::ptree::value_type& valueNode : valueNodes) {
0915         if (valueNode.first == "LutVal") {
0916           int key = valueNode.second.get<int>("<xmlattr>.key");
0917           float value = valueNode.second.get<float>("<xmlattr>.value");
0918           extrapolFactors.at(iRefLayer).at(iLayer)[key] = value;
0919           LogTrace("OMTFReconstruction") << "key " << key << " value " << value << std::endl;
0920         }
0921       }
0922     }
0923   }
0924 }
0925 
0926 /////////////////////////////////////////////////////////
0927 
0928 template class OMTFProcessor<GoldenPattern>;
0929 template class OMTFProcessor<GoldenPatternWithStat>;
0930 template class OMTFProcessor<GoldenPatternWithThresh>;