Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    L1Trigger/L1TTrackMatch
0004 // Class:      L1TrackSelectionProducer
0005 //
0006 /**\class L1TrackSelectionProducer L1TrackSelectionProducer.cc L1Trigger/L1TTrackMatch/plugins/L1TrackSelectionProducer.cc
0007 
0008  Description: Selects a set of L1Tracks based on a set of predefined criteria.
0009 
0010  Implementation:
0011      Inputs:
0012          std::vector<TTTrack> - Each floating point TTTrack inside this collection inherits from
0013                                 a bit-accurate TTTrack_TrackWord, used for emulation purposes.
0014      Outputs:
0015          std::vector<TTTrack> - A collection of TTTracks selected from cuts on the TTTrack properties
0016          std::vector<TTTrack> - A collection of TTTracks selected from cuts on the TTTrack_TrackWord properties
0017 */
0018 //
0019 // Original Author:  Alexx Perloff
0020 //         Created:  Thu, 16 Dec 2021 19:02:50 GMT
0021 //
0022 // Updates: Claire Savard (claire.savard@colorado.edu), Nov. 2023
0023 //
0024 
0025 // system include files
0026 #include <algorithm>
0027 #include <memory>
0028 #include <string>
0029 #include <vector>
0030 
0031 // Xilinx HLS includes
0032 #include <ap_fixed.h>
0033 #include <ap_int.h>
0034 
0035 // user include files
0036 #include "DataFormats/Common/interface/Handle.h"
0037 #include "DataFormats/Common/interface/Ref.h"
0038 #include "DataFormats/Common/interface/RefVector.h"
0039 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0040 #include "DataFormats/L1Trigger/interface/Vertex.h"
0041 #include "DataFormats/L1Trigger/interface/VertexWord.h"
0042 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0043 #include "CommonTools/Utils/interface/AndSelector.h"
0044 #include "CommonTools/Utils/interface/EtaRangeSelector.h"
0045 #include "CommonTools/Utils/interface/MinSelector.h"
0046 #include "CommonTools/Utils/interface/MinFunctionSelector.h"
0047 #include "CommonTools/Utils/interface/MinNumberSelector.h"
0048 #include "CommonTools/Utils/interface/PtMinSelector.h"
0049 #include "CommonTools/Utils/interface/Selection.h"
0050 #include "FWCore/Framework/interface/Frameworkfwd.h"
0051 #include "FWCore/Framework/interface/global/EDProducer.h"
0052 #include "FWCore/Framework/interface/Event.h"
0053 #include "FWCore/Framework/interface/MakerMacros.h"
0054 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0055 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0056 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0057 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0058 #include "FWCore/Utilities/interface/EDMException.h"
0059 #include "FWCore/Utilities/interface/StreamID.h"
0060 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0061 
0062 //
0063 // class declaration
0064 //
0065 
0066 class L1TrackSelectionProducer : public edm::global::EDProducer<> {
0067 public:
0068   explicit L1TrackSelectionProducer(const edm::ParameterSet&);
0069   ~L1TrackSelectionProducer() override;
0070 
0071   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0072 
0073 private:
0074   // ----------constants, enums and typedefs ---------
0075   // Relevant constants for the converted track word
0076   enum TrackBitWidths {
0077     kPtSize = TTTrack_TrackWord::TrackBitWidths::kRinvSize - 1,  // Width of pt
0078     kPtMagSize = 9,                                              // Width of pt magnitude (unsigned)
0079     kEtaSize = TTTrack_TrackWord::TrackBitWidths::kTanlSize,     // Width of eta
0080     kEtaMagSize = 3,                                             // Width of eta magnitude (signed)
0081   };
0082 
0083   typedef TTTrack<Ref_Phase2TrackerDigi_> L1Track;
0084   typedef std::vector<L1Track> TTTrackCollection;
0085   typedef edm::Handle<TTTrackCollection> TTTrackCollectionHandle;
0086   typedef edm::Ref<TTTrackCollection> TTTrackRef;
0087   typedef edm::RefVector<TTTrackCollection> TTTrackRefCollection;
0088   typedef std::unique_ptr<TTTrackRefCollection> TTTrackRefCollectionUPtr;
0089 
0090   // ----------member functions ----------------------
0091   void printDebugInfo(const TTTrackCollectionHandle& l1TracksHandle,
0092                       const TTTrackRefCollectionUPtr& vTTTrackOutput,
0093                       const TTTrackRefCollectionUPtr& vTTTrackEmulationOutput) const;
0094   void printTrackInfo(edm::LogInfo& log, const L1Track& track, bool printEmulation = false) const;
0095   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0096 
0097   // ----------selectors -----------------------------
0098   // Based on recommendations from https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideGenericSelectors
0099   struct TTTrackPtMinSelector {
0100     TTTrackPtMinSelector(double ptMin) : ptMin_(ptMin) {}
0101     TTTrackPtMinSelector(const edm::ParameterSet& cfg) : ptMin_(cfg.template getParameter<double>("ptMin")) {}
0102     bool operator()(const L1Track& t) const { return t.momentum().perp() >= ptMin_; }
0103 
0104   private:
0105     double ptMin_;
0106   };
0107   struct TTTrackWordPtMinSelector {
0108     TTTrackWordPtMinSelector(double ptMin) : ptMin_(ptMin) {}
0109     TTTrackWordPtMinSelector(const edm::ParameterSet& cfg) : ptMin_(cfg.template getParameter<double>("ptMin")) {}
0110     bool operator()(const L1Track& t) const {
0111       ap_uint<TrackBitWidths::kPtSize> ptEmulationBits = t.getTrackWord()(
0112           TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1, TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
0113       ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize> ptEmulation;
0114       ptEmulation.V = ptEmulationBits.range();
0115       return ptEmulation.to_double() >= ptMin_;
0116     }
0117 
0118   private:
0119     double ptMin_;
0120   };
0121   struct TTTrackAbsEtaMaxSelector {
0122     TTTrackAbsEtaMaxSelector(double absEtaMax) : absEtaMax_(absEtaMax) {}
0123     TTTrackAbsEtaMaxSelector(const edm::ParameterSet& cfg)
0124         : absEtaMax_(cfg.template getParameter<double>("absEtaMax")) {}
0125     bool operator()(const L1Track& t) const { return std::abs(t.momentum().eta()) <= absEtaMax_; }
0126 
0127   private:
0128     double absEtaMax_;
0129   };
0130   struct TTTrackWordAbsEtaMaxSelector {
0131     TTTrackWordAbsEtaMaxSelector(double absEtaMax) : absEtaMax_(absEtaMax) {}
0132     TTTrackWordAbsEtaMaxSelector(const edm::ParameterSet& cfg)
0133         : absEtaMax_(cfg.template getParameter<double>("absEtaMax")) {}
0134     bool operator()(const L1Track& t) const {
0135       TTTrack_TrackWord::tanl_t etaEmulationBits = t.getTanlWord();
0136       ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize> etaEmulation;
0137       etaEmulation.V = etaEmulationBits.range();
0138       return std::abs(etaEmulation.to_double()) <= absEtaMax_;
0139     }
0140 
0141   private:
0142     double absEtaMax_;
0143   };
0144   struct TTTrackAbsZ0MaxSelector {
0145     TTTrackAbsZ0MaxSelector(double absZ0Max) : absZ0Max_(absZ0Max) {}
0146     TTTrackAbsZ0MaxSelector(const edm::ParameterSet& cfg) : absZ0Max_(cfg.template getParameter<double>("absZ0Max")) {}
0147     bool operator()(const L1Track& t) const { return std::abs(t.z0()) <= absZ0Max_; }
0148 
0149   private:
0150     double absZ0Max_;
0151   };
0152   struct TTTrackWordAbsZ0MaxSelector {
0153     TTTrackWordAbsZ0MaxSelector(double absZ0Max) : absZ0Max_(absZ0Max) {}
0154     TTTrackWordAbsZ0MaxSelector(const edm::ParameterSet& cfg)
0155         : absZ0Max_(cfg.template getParameter<double>("absZ0Max")) {}
0156     bool operator()(const L1Track& t) const {
0157       double floatZ0 = t.undigitizeSignedValue(
0158           t.getZ0Bits(), TTTrack_TrackWord::TrackBitWidths::kZ0Size, TTTrack_TrackWord::stepZ0, 0.0);
0159       return std::abs(floatZ0) <= absZ0Max_;
0160     }
0161 
0162   private:
0163     double absZ0Max_;
0164   };
0165   struct TTTrackNStubsMinSelector {
0166     TTTrackNStubsMinSelector(double nStubsMin) : nStubsMin_(nStubsMin) {}
0167     TTTrackNStubsMinSelector(const edm::ParameterSet& cfg)
0168         : nStubsMin_(cfg.template getParameter<double>("nStubsMin")) {}
0169     bool operator()(const L1Track& t) const { return t.getStubRefs().size() >= nStubsMin_; }
0170 
0171   private:
0172     double nStubsMin_;
0173   };
0174   struct TTTrackWordNStubsMinSelector {
0175     TTTrackWordNStubsMinSelector(double nStubsMin) : nStubsMin_(nStubsMin) {}
0176     TTTrackWordNStubsMinSelector(const edm::ParameterSet& cfg)
0177         : nStubsMin_(cfg.template getParameter<double>("nStubsMin")) {}
0178     bool operator()(const L1Track& t) const { return t.getNStubs() >= nStubsMin_; }
0179 
0180   private:
0181     double nStubsMin_;
0182   };
0183   struct TTTrackNPSStubsMinSelector {
0184     TTTrackNPSStubsMinSelector(double nStubsMin, const TrackerTopology& tTopo)
0185         : nPSStubsMin_(nStubsMin), tTopo_(tTopo) {}
0186     TTTrackNPSStubsMinSelector(const edm::ParameterSet& cfg, const TrackerTopology& tTopo)
0187         : nPSStubsMin_(cfg.template getParameter<double>("nPSStubsMin")), tTopo_(tTopo) {}
0188     bool operator()(const L1Track& t) const {
0189       int nPSStubs = 0;
0190       for (const auto& stub : t.getStubRefs()) {
0191         DetId detId(stub->getDetId());
0192         if (detId.det() == DetId::Detector::Tracker) {
0193           if ((detId.subdetId() == StripSubdetector::TOB && tTopo_.tobLayer(detId) <= 3) ||
0194               (detId.subdetId() == StripSubdetector::TID && tTopo_.tidRing(detId) <= 9))
0195             nPSStubs++;
0196         }
0197       }
0198       return nPSStubs >= nPSStubsMin_;
0199     }
0200 
0201   private:
0202     double nPSStubsMin_;
0203     const TrackerTopology& tTopo_;
0204   };
0205   struct TTTrackPromptMVAMinSelector {
0206     TTTrackPromptMVAMinSelector(double promptMVAMin) : promptMVAMin_(promptMVAMin) {}
0207     TTTrackPromptMVAMinSelector(const edm::ParameterSet& cfg)
0208         : promptMVAMin_(cfg.template getParameter<double>("promptMVAMin")) {}
0209     bool operator()(const L1Track& t) const { return t.trkMVA1() >= promptMVAMin_; }
0210 
0211   private:
0212     double promptMVAMin_;
0213   };
0214   struct TTTrackWordPromptMVAMinSelector {
0215     TTTrackWordPromptMVAMinSelector(double promptMVAMin) : promptMVAMin_(promptMVAMin) {}
0216     TTTrackWordPromptMVAMinSelector(const edm::ParameterSet& cfg)
0217         : promptMVAMin_(cfg.template getParameter<double>("promptMVAMin")) {}
0218     bool operator()(const L1Track& t) const {
0219       return t.trkMVA1() >= promptMVAMin_;
0220     }  //change when mva bins in word are set
0221 
0222   private:
0223     double promptMVAMin_;
0224   };
0225   struct TTTrackBendChi2MaxSelector {
0226     TTTrackBendChi2MaxSelector(double bendChi2Max) : bendChi2Max_(bendChi2Max) {}
0227     TTTrackBendChi2MaxSelector(const edm::ParameterSet& cfg)
0228         : bendChi2Max_(cfg.template getParameter<double>("reducedBendChi2Max")) {}
0229     bool operator()(const L1Track& t) const { return t.stubPtConsistency() < bendChi2Max_; }
0230 
0231   private:
0232     double bendChi2Max_;
0233   };
0234   struct TTTrackWordBendChi2MaxSelector {
0235     TTTrackWordBendChi2MaxSelector(double bendChi2Max) : bendChi2Max_(bendChi2Max) {}
0236     TTTrackWordBendChi2MaxSelector(const edm::ParameterSet& cfg)
0237         : bendChi2Max_(cfg.template getParameter<double>("reducedBendChi2Max")) {}
0238     bool operator()(const L1Track& t) const { return t.getBendChi2() < bendChi2Max_; }
0239 
0240   private:
0241     double bendChi2Max_;
0242   };
0243   struct TTTrackChi2RZMaxSelector {
0244     TTTrackChi2RZMaxSelector(double reducedChi2RZMax) : reducedChi2RZMax_(reducedChi2RZMax) {}
0245     TTTrackChi2RZMaxSelector(const edm::ParameterSet& cfg)
0246         : reducedChi2RZMax_(cfg.template getParameter<double>("reducedChi2RZMax")) {}
0247     bool operator()(const L1Track& t) const { return t.chi2ZRed() < reducedChi2RZMax_; }
0248 
0249   private:
0250     double reducedChi2RZMax_;
0251   };
0252   struct TTTrackWordChi2RZMaxSelector {
0253     TTTrackWordChi2RZMaxSelector(double reducedChi2RZMax) : reducedChi2RZMax_(reducedChi2RZMax) {}
0254     TTTrackWordChi2RZMaxSelector(const edm::ParameterSet& cfg)
0255         : reducedChi2RZMax_(cfg.template getParameter<double>("reducedChi2RZMax")) {}
0256     bool operator()(const L1Track& t) const { return t.getChi2RZ() < reducedChi2RZMax_; }
0257 
0258   private:
0259     double reducedChi2RZMax_;
0260   };
0261   struct TTTrackChi2RPhiMaxSelector {
0262     TTTrackChi2RPhiMaxSelector(double reducedChi2RPhiMax) : reducedChi2RPhiMax_(reducedChi2RPhiMax) {}
0263     TTTrackChi2RPhiMaxSelector(const edm::ParameterSet& cfg)
0264         : reducedChi2RPhiMax_(cfg.template getParameter<double>("reducedChi2RPhiMax")) {}
0265     bool operator()(const L1Track& t) const { return t.chi2XYRed() < reducedChi2RPhiMax_; }
0266 
0267   private:
0268     double reducedChi2RPhiMax_;
0269   };
0270   struct TTTrackWordChi2RPhiMaxSelector {
0271     TTTrackWordChi2RPhiMaxSelector(double reducedChi2RPhiMax) : reducedChi2RPhiMax_(reducedChi2RPhiMax) {}
0272     TTTrackWordChi2RPhiMaxSelector(const edm::ParameterSet& cfg)
0273         : reducedChi2RPhiMax_(cfg.template getParameter<double>("reducedChi2RPhiMax")) {}
0274     bool operator()(const L1Track& t) const { return t.getChi2RPhi() < reducedChi2RPhiMax_; }
0275 
0276   private:
0277     double reducedChi2RPhiMax_;
0278   };
0279   struct TTTrackChi2RZMaxNstubSelector {
0280     TTTrackChi2RZMaxNstubSelector(double reducedChi2RZMaxNstub4, double reducedChi2RZMaxNstub5)
0281         : reducedChi2RZMaxNstub4_(reducedChi2RZMaxNstub4), reducedChi2RZMaxNstub5_(reducedChi2RZMaxNstub5) {}
0282     TTTrackChi2RZMaxNstubSelector(const edm::ParameterSet& cfg)
0283         : reducedChi2RZMaxNstub4_(cfg.template getParameter<double>("reducedChi2RZMaxNstub4")),
0284           reducedChi2RZMaxNstub5_(cfg.template getParameter<double>("reducedChi2RZMaxNstub5")) {}
0285     bool operator()(const L1Track& t) const {
0286       return (((t.chi2ZRed() < reducedChi2RZMaxNstub4_) && (t.getStubRefs().size() == 4)) ||
0287               ((t.chi2ZRed() < reducedChi2RZMaxNstub5_) && (t.getStubRefs().size() > 4)));
0288     }
0289 
0290   private:
0291     double reducedChi2RZMaxNstub4_;
0292     double reducedChi2RZMaxNstub5_;
0293   };
0294   struct TTTrackWordChi2RZMaxNstubSelector {
0295     TTTrackWordChi2RZMaxNstubSelector(double reducedChi2RZMaxNstub4, double reducedChi2RZMaxNstub5)
0296         : reducedChi2RZMaxNstub4_(reducedChi2RZMaxNstub4), reducedChi2RZMaxNstub5_(reducedChi2RZMaxNstub5) {}
0297     TTTrackWordChi2RZMaxNstubSelector(const edm::ParameterSet& cfg)
0298         : reducedChi2RZMaxNstub4_(cfg.template getParameter<double>("reducedChi2RZMaxNstub4")),
0299           reducedChi2RZMaxNstub5_(cfg.template getParameter<double>("reducedChi2RZMaxNstub5")) {}
0300     bool operator()(const L1Track& t) const {
0301       return (((t.getChi2RZ() < reducedChi2RZMaxNstub4_) && (t.getNStubs() == 4)) ||
0302               ((t.getChi2RZ() < reducedChi2RZMaxNstub5_) && (t.getNStubs() > 4)));
0303     }
0304 
0305   private:
0306     double reducedChi2RZMaxNstub4_;
0307     double reducedChi2RZMaxNstub5_;
0308   };
0309   struct TTTrackChi2RPhiMaxNstubSelector {
0310     TTTrackChi2RPhiMaxNstubSelector(double reducedChi2RPhiMaxNstub4, double reducedChi2RPhiMaxNstub5)
0311         : reducedChi2RPhiMaxNstub4_(reducedChi2RPhiMaxNstub4), reducedChi2RPhiMaxNstub5_(reducedChi2RPhiMaxNstub5) {}
0312     TTTrackChi2RPhiMaxNstubSelector(const edm::ParameterSet& cfg)
0313         : reducedChi2RPhiMaxNstub4_(cfg.template getParameter<double>("reducedChi2RPhiMaxNstub4")),
0314           reducedChi2RPhiMaxNstub5_(cfg.template getParameter<double>("reducedChi2RPhiMaxNstub5")) {}
0315     bool operator()(const L1Track& t) const {
0316       return (((t.chi2XYRed() < reducedChi2RPhiMaxNstub4_) && (t.getStubRefs().size() == 4)) ||
0317               ((t.chi2XYRed() < reducedChi2RPhiMaxNstub5_) && (t.getStubRefs().size() > 4)));
0318     }
0319 
0320   private:
0321     double reducedChi2RPhiMaxNstub4_;
0322     double reducedChi2RPhiMaxNstub5_;
0323   };
0324   struct TTTrackWordChi2RPhiMaxNstubSelector {  // using simulated chi2 since not implemented in track word, updates needed
0325     TTTrackWordChi2RPhiMaxNstubSelector(double reducedChi2RPhiMaxNstub4, double reducedChi2RPhiMaxNstub5)
0326         : reducedChi2RPhiMaxNstub4_(reducedChi2RPhiMaxNstub4), reducedChi2RPhiMaxNstub5_(reducedChi2RPhiMaxNstub5) {}
0327     TTTrackWordChi2RPhiMaxNstubSelector(const edm::ParameterSet& cfg)
0328         : reducedChi2RPhiMaxNstub4_(cfg.template getParameter<double>("reducedChi2RPhiMaxNstub4")),
0329           reducedChi2RPhiMaxNstub5_(cfg.template getParameter<double>("reducedChi2RPhiMaxNstub5")) {}
0330     bool operator()(const L1Track& t) const {
0331       return (((t.getChi2RPhi() < reducedChi2RPhiMaxNstub4_) && (t.getNStubs() == 4)) ||
0332               ((t.getChi2RPhi() < reducedChi2RPhiMaxNstub5_) && (t.getNStubs() > 4)));
0333     }
0334 
0335   private:
0336     double reducedChi2RPhiMaxNstub4_;
0337     double reducedChi2RPhiMaxNstub5_;
0338   };
0339   struct TTTrackBendChi2MaxNstubSelector {
0340     TTTrackBendChi2MaxNstubSelector(double reducedBendChi2MaxNstub4, double reducedBendChi2MaxNstub5)
0341         : reducedBendChi2MaxNstub4_(reducedBendChi2MaxNstub4), reducedBendChi2MaxNstub5_(reducedBendChi2MaxNstub5) {}
0342     TTTrackBendChi2MaxNstubSelector(const edm::ParameterSet& cfg)
0343         : reducedBendChi2MaxNstub4_(cfg.template getParameter<double>("reducedBendChi2MaxNstub4")),
0344           reducedBendChi2MaxNstub5_(cfg.template getParameter<double>("reducedBendChi2MaxNstub5")) {}
0345     bool operator()(const L1Track& t) const {
0346       return (((t.stubPtConsistency() < reducedBendChi2MaxNstub4_) && (t.getStubRefs().size() == 4)) ||
0347               ((t.stubPtConsistency() < reducedBendChi2MaxNstub5_) && (t.getStubRefs().size() > 4)));
0348     }
0349 
0350   private:
0351     double reducedBendChi2MaxNstub4_;
0352     double reducedBendChi2MaxNstub5_;
0353   };
0354   struct TTTrackWordBendChi2MaxNstubSelector {
0355     TTTrackWordBendChi2MaxNstubSelector(double reducedBendChi2MaxNstub4, double reducedBendChi2MaxNstub5)
0356         : reducedBendChi2MaxNstub4_(reducedBendChi2MaxNstub4), reducedBendChi2MaxNstub5_(reducedBendChi2MaxNstub5) {}
0357     TTTrackWordBendChi2MaxNstubSelector(const edm::ParameterSet& cfg)
0358         : reducedBendChi2MaxNstub4_(cfg.template getParameter<double>("reducedBendChi2MaxNstub4")),
0359           reducedBendChi2MaxNstub5_(cfg.template getParameter<double>("reducedBendChi2MaxNstub5")) {}
0360     bool operator()(const L1Track& t) const {
0361       return (((t.getBendChi2() < reducedBendChi2MaxNstub4_) && (t.getNStubs() == 4)) ||
0362               ((t.getBendChi2() < reducedBendChi2MaxNstub5_) && (t.getNStubs() > 4)));
0363     }
0364 
0365   private:
0366     double reducedBendChi2MaxNstub4_;
0367     double reducedBendChi2MaxNstub5_;
0368   };
0369 
0370   typedef AndSelector<TTTrackPtMinSelector, TTTrackAbsEtaMaxSelector, TTTrackAbsZ0MaxSelector, TTTrackNStubsMinSelector>
0371       TTTrackPtMinEtaMaxZ0MaxNStubsMinSelector;
0372   typedef AndSelector<TTTrackWordPtMinSelector,
0373                       TTTrackWordAbsEtaMaxSelector,
0374                       TTTrackWordAbsZ0MaxSelector,
0375                       TTTrackWordNStubsMinSelector>
0376       TTTrackWordPtMinEtaMaxZ0MaxNStubsMinSelector;
0377   typedef AndSelector<TTTrackBendChi2MaxSelector, TTTrackChi2RZMaxSelector, TTTrackChi2RPhiMaxSelector>
0378       TTTrackBendChi2Chi2RZChi2RPhiMaxSelector;
0379   typedef AndSelector<TTTrackWordBendChi2MaxSelector, TTTrackWordChi2RZMaxSelector, TTTrackWordChi2RPhiMaxSelector>
0380       TTTrackWordBendChi2Chi2RZChi2RPhiMaxSelector;
0381   typedef AndSelector<TTTrackChi2RZMaxNstubSelector, TTTrackChi2RPhiMaxNstubSelector, TTTrackBendChi2MaxNstubSelector>
0382       TTTrackChi2MaxNstubSelector;
0383   typedef AndSelector<TTTrackWordChi2RZMaxNstubSelector,
0384                       TTTrackWordChi2RPhiMaxNstubSelector,
0385                       TTTrackWordBendChi2MaxNstubSelector>
0386       TTTrackWordChi2MaxNstubSelector;
0387 
0388   // ----------member data ---------------------------
0389   const edm::EDGetTokenT<TTTrackCollection> l1TracksToken_;
0390   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
0391   const std::string outputCollectionName_;
0392   const edm::ParameterSet cutSet_;
0393   const double ptMin_, absEtaMax_, absZ0Max_, promptMVAMin_, bendChi2Max_, reducedChi2RZMax_, reducedChi2RPhiMax_;
0394   const double reducedChi2RZMaxNstub4_, reducedChi2RZMaxNstub5_, reducedChi2RPhiMaxNstub4_, reducedChi2RPhiMaxNstub5_,
0395       reducedBendChi2MaxNstub4_, reducedBendChi2MaxNstub5_;
0396   const int nStubsMin_, nPSStubsMin_;
0397   bool processSimulatedTracks_, processEmulatedTracks_;
0398   int debug_;
0399 };
0400 
0401 //
0402 // constructors and destructor
0403 //
0404 L1TrackSelectionProducer::L1TrackSelectionProducer(const edm::ParameterSet& iConfig)
0405     : l1TracksToken_(consumes<TTTrackCollection>(iConfig.getParameter<edm::InputTag>("l1TracksInputTag"))),
0406       tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))),
0407       outputCollectionName_(iConfig.getParameter<std::string>("outputCollectionName")),
0408       cutSet_(iConfig.getParameter<edm::ParameterSet>("cutSet")),
0409 
0410       ptMin_(cutSet_.getParameter<double>("ptMin")),
0411       absEtaMax_(cutSet_.getParameter<double>("absEtaMax")),
0412       absZ0Max_(cutSet_.getParameter<double>("absZ0Max")),
0413       promptMVAMin_(cutSet_.getParameter<double>("promptMVAMin")),
0414       bendChi2Max_(cutSet_.getParameter<double>("reducedBendChi2Max")),
0415       reducedChi2RZMax_(cutSet_.getParameter<double>("reducedChi2RZMax")),
0416       reducedChi2RPhiMax_(cutSet_.getParameter<double>("reducedChi2RPhiMax")),
0417       reducedChi2RZMaxNstub4_(cutSet_.getParameter<double>("reducedChi2RZMaxNstub4")),
0418       reducedChi2RZMaxNstub5_(cutSet_.getParameter<double>("reducedChi2RZMaxNstub5")),
0419       reducedChi2RPhiMaxNstub4_(cutSet_.getParameter<double>("reducedChi2RPhiMaxNstub4")),
0420       reducedChi2RPhiMaxNstub5_(cutSet_.getParameter<double>("reducedChi2RPhiMaxNstub5")),
0421       reducedBendChi2MaxNstub4_(cutSet_.getParameter<double>("reducedBendChi2MaxNstub4")),
0422       reducedBendChi2MaxNstub5_(cutSet_.getParameter<double>("reducedBendChi2MaxNstub5")),
0423       nStubsMin_(cutSet_.getParameter<int>("nStubsMin")),
0424       nPSStubsMin_(cutSet_.getParameter<int>("nPSStubsMin")),
0425       processSimulatedTracks_(iConfig.getParameter<bool>("processSimulatedTracks")),
0426       processEmulatedTracks_(iConfig.getParameter<bool>("processEmulatedTracks")),
0427       debug_(iConfig.getParameter<int>("debug")) {
0428   // Confirm the the configuration makes sense
0429   if (!processSimulatedTracks_ && !processEmulatedTracks_) {
0430     throw cms::Exception("You must process at least one of the track collections (simulated or emulated).");
0431   }
0432 
0433   if (processSimulatedTracks_) {
0434     produces<TTTrackRefCollection>(outputCollectionName_);
0435   }
0436   if (processEmulatedTracks_) {
0437     produces<TTTrackRefCollection>(outputCollectionName_ + "Emulation");
0438   }
0439 }
0440 
0441 L1TrackSelectionProducer::~L1TrackSelectionProducer() {}
0442 
0443 //
0444 // member functions
0445 //
0446 
0447 void L1TrackSelectionProducer::printDebugInfo(const TTTrackCollectionHandle& l1TracksHandle,
0448                                               const TTTrackRefCollectionUPtr& vTTTrackOutput,
0449                                               const TTTrackRefCollectionUPtr& vTTTrackEmulationOutput) const {
0450   edm::LogInfo log("L1TrackSelectionProducer");
0451   log << "The original track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are ... \n";
0452   for (const auto& track : *l1TracksHandle) {
0453     printTrackInfo(log, track, debug_ >= 4);
0454   }
0455   log << "\t---\n\tNumber of tracks in this selection = " << l1TracksHandle->size() << "\n\n";
0456   if (processSimulatedTracks_) {
0457     log << "The selected track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are ... \n";
0458     for (const auto& track : *vTTTrackOutput) {
0459       printTrackInfo(log, *track, debug_ >= 4);
0460     }
0461     log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackOutput->size() << "\n\n";
0462   }
0463   if (processEmulatedTracks_) {
0464     log << "The emulation selected track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are "
0465            "... \n";
0466     for (const auto& track : *vTTTrackEmulationOutput) {
0467       printTrackInfo(log, *track, debug_ >= 4);
0468     }
0469     log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackEmulationOutput->size() << "\n\n";
0470   }
0471   if (processSimulatedTracks_ && processEmulatedTracks_) {
0472     TTTrackRefCollection inSimButNotEmu;
0473     TTTrackRefCollection inEmuButNotSim;
0474     std::set_difference(vTTTrackOutput->begin(),
0475                         vTTTrackOutput->end(),
0476                         vTTTrackEmulationOutput->begin(),
0477                         vTTTrackEmulationOutput->end(),
0478                         std::back_inserter(inSimButNotEmu));
0479     std::set_difference(vTTTrackEmulationOutput->begin(),
0480                         vTTTrackEmulationOutput->end(),
0481                         vTTTrackOutput->begin(),
0482                         vTTTrackOutput->end(),
0483                         std::back_inserter(inEmuButNotSim));
0484     log << "The set of tracks selected via cuts on the simulated values which are not in the set of tracks selected "
0485            "by cutting on the emulated values ... \n";
0486     for (const auto& track : inSimButNotEmu) {
0487       printTrackInfo(log, *track, debug_ >= 3);
0488     }
0489     log << "\t---\n\tNumber of tracks in this selection = " << inSimButNotEmu.size() << "\n\n"
0490         << "The set of tracks selected via cuts on the emulated values which are not in the set of tracks selected "
0491            "by cutting on the simulated values ... \n";
0492     for (const auto& track : inEmuButNotSim) {
0493       printTrackInfo(log, *track, debug_ >= 3);
0494     }
0495     log << "\t---\n\tNumber of tracks in this selection = " << inEmuButNotSim.size() << "\n\n";
0496   }
0497 }
0498 
0499 void L1TrackSelectionProducer::printTrackInfo(edm::LogInfo& log, const L1Track& track, bool printEmulation) const {
0500   log << "\t(" << track.momentum().perp() << ", " << track.momentum().eta() << ", " << track.momentum().phi() << ", "
0501       << track.getStubRefs().size() << ", " << track.stubPtConsistency() << ", " << track.chi2ZRed() << ", "
0502       << track.chi2XYRed() << ", " << track.z0() << ")\n";
0503 
0504   if (printEmulation) {
0505     ap_uint<TrackBitWidths::kPtSize> ptEmulationBits = track.getTrackWord()(
0506         TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1, TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
0507     ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize> ptEmulation;
0508     ptEmulation.V = ptEmulationBits.range();
0509     TTTrack_TrackWord::tanl_t etaEmulationBits = track.getTanlWord();
0510     ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize> etaEmulation;
0511     etaEmulation.V = etaEmulationBits.range();
0512     double floatTkZ0 = track.undigitizeSignedValue(
0513         track.getZ0Bits(), TTTrack_TrackWord::TrackBitWidths::kZ0Size, TTTrack_TrackWord::stepZ0, 0.0);
0514     double floatTkPhi = track.undigitizeSignedValue(
0515         track.getPhiBits(), TTTrack_TrackWord::TrackBitWidths::kPhiSize, TTTrack_TrackWord::stepPhi0, 0.0);
0516     log << "\t\t(" << ptEmulation.to_double() << ", " << etaEmulation.to_double() << ", " << floatTkPhi << ", "
0517         << track.getNStubs() << ", " << track.getBendChi2() << ", " << track.getChi2RZ() << ", " << track.getChi2RPhi()
0518         << ", " << floatTkZ0 << ")\n";
0519   }
0520 }
0521 
0522 // ------------ method called to produce the data  ------------
0523 void L1TrackSelectionProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0524   auto vTTTrackOutput = std::make_unique<TTTrackRefCollection>();
0525   auto vTTTrackEmulationOutput = std::make_unique<TTTrackRefCollection>();
0526 
0527   // Tracker Topology
0528   const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
0529 
0530   TTTrackCollectionHandle l1TracksHandle;
0531 
0532   iEvent.getByToken(l1TracksToken_, l1TracksHandle);
0533   size_t nOutputApproximate = l1TracksHandle->size();
0534   if (processSimulatedTracks_) {
0535     vTTTrackOutput->reserve(nOutputApproximate);
0536   }
0537   if (processEmulatedTracks_) {
0538     vTTTrackEmulationOutput->reserve(nOutputApproximate);
0539   }
0540 
0541   TTTrackPtMinEtaMaxZ0MaxNStubsMinSelector kinSel(ptMin_, absEtaMax_, absZ0Max_, nStubsMin_);
0542   TTTrackWordPtMinEtaMaxZ0MaxNStubsMinSelector kinSelEmu(ptMin_, absEtaMax_, absZ0Max_, nStubsMin_);
0543   TTTrackBendChi2Chi2RZChi2RPhiMaxSelector chi2Sel(bendChi2Max_, reducedChi2RZMax_, reducedChi2RPhiMax_);
0544   TTTrackWordBendChi2Chi2RZChi2RPhiMaxSelector chi2SelEmu(bendChi2Max_, reducedChi2RZMax_, reducedChi2RPhiMax_);
0545   TTTrackNPSStubsMinSelector nPSStubsSel(nPSStubsMin_, tTopo);
0546   TTTrackPromptMVAMinSelector mvaSel(promptMVAMin_);
0547   TTTrackWordPromptMVAMinSelector mvaSelEmu(promptMVAMin_);
0548   TTTrackChi2MaxNstubSelector chi2NstubSel({reducedChi2RZMaxNstub4_, reducedChi2RZMaxNstub5_},
0549                                            {reducedChi2RPhiMaxNstub4_, reducedChi2RPhiMaxNstub5_},
0550                                            {reducedBendChi2MaxNstub4_, reducedBendChi2MaxNstub5_});
0551   TTTrackWordChi2MaxNstubSelector chi2NstubSelEmu({reducedChi2RZMaxNstub4_, reducedChi2RZMaxNstub5_},
0552                                                   {reducedChi2RPhiMaxNstub4_, reducedChi2RPhiMaxNstub5_},
0553                                                   {reducedBendChi2MaxNstub4_, reducedBendChi2MaxNstub5_});
0554 
0555   for (size_t i = 0; i < nOutputApproximate; i++) {
0556     const auto& track = l1TracksHandle->at(i);
0557 
0558     // Select tracks based on the floating point TTTrack
0559     if (processSimulatedTracks_ && kinSel(track) && nPSStubsSel(track) && chi2Sel(track) && mvaSel(track) &&
0560         chi2NstubSel(track)) {
0561       vTTTrackOutput->push_back(TTTrackRef(l1TracksHandle, i));
0562     }
0563 
0564     // Select tracks based on the bitwise accurate TTTrack_TrackWord
0565     if (processEmulatedTracks_ && kinSelEmu(track) && chi2SelEmu(track) && mvaSelEmu(track) && chi2NstubSelEmu(track)) {
0566       vTTTrackEmulationOutput->push_back(TTTrackRef(l1TracksHandle, i));
0567     }
0568   }
0569 
0570   if (debug_ >= 2) {
0571     printDebugInfo(l1TracksHandle, vTTTrackOutput, vTTTrackEmulationOutput);
0572   }
0573 
0574   // Put the outputs into the event
0575   if (processSimulatedTracks_) {
0576     iEvent.put(std::move(vTTTrackOutput), outputCollectionName_);
0577   }
0578   if (processEmulatedTracks_) {
0579     iEvent.put(std::move(vTTTrackEmulationOutput), outputCollectionName_ + "Emulation");
0580   }
0581 }
0582 
0583 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0584 void L1TrackSelectionProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0585   //L1TrackSelectionProducer
0586   edm::ParameterSetDescription desc;
0587   desc.add<edm::InputTag>("l1TracksInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
0588   desc.add<std::string>("outputCollectionName", "Level1TTTracksSelected");
0589   {
0590     edm::ParameterSetDescription descCutSet;
0591     descCutSet.add<double>("ptMin", 2.0)->setComment("pt must be greater than this value, [GeV]");
0592     descCutSet.add<double>("absEtaMax", 2.4)->setComment("absolute value of eta must be less than this value");
0593     descCutSet.add<double>("absZ0Max", 15.0)->setComment("z0 must be less than this value, [cm]");
0594     descCutSet.add<int>("nStubsMin", 4)->setComment("number of stubs must be greater than or equal to this value");
0595     descCutSet.add<int>("nPSStubsMin", 0)
0596         ->setComment("number of stubs in the PS Modules must be greater than or equal to this value");
0597 
0598     descCutSet.add<double>("promptMVAMin", -1.0)->setComment("MVA must be greater than this value");
0599     descCutSet.add<double>("reducedBendChi2Max", 2.25)->setComment("bend chi2 must be less than this value");
0600     descCutSet.add<double>("reducedChi2RZMax", 5.0)->setComment("chi2rz/dof must be less than this value");
0601     descCutSet.add<double>("reducedChi2RPhiMax", 20.0)->setComment("chi2rphi/dof must be less than this value");
0602     descCutSet.add<double>("reducedChi2RZMaxNstub4", 999.9)
0603         ->setComment("chi2rz/dof must be less than this value in nstub==4");
0604     descCutSet.add<double>("reducedChi2RZMaxNstub5", 999.9)
0605         ->setComment("chi2rz/dof must be less than this value in nstub>4");
0606     descCutSet.add<double>("reducedChi2RPhiMaxNstub4", 999.9)
0607         ->setComment("chi2rphi/dof must be less than this value in nstub==4");
0608     descCutSet.add<double>("reducedChi2RPhiMaxNstub5", 999.9)
0609         ->setComment("chi2rphi/dof must be less than this value in nstub>4");
0610     descCutSet.add<double>("reducedBendChi2MaxNstub4", 999.9)
0611         ->setComment("bend chi2 must be less than this value in nstub==4");
0612     descCutSet.add<double>("reducedBendChi2MaxNstub5", 999.9)
0613         ->setComment("bend chi2 must be less than this value in nstub>4");
0614 
0615     desc.add<edm::ParameterSetDescription>("cutSet", descCutSet);
0616   }
0617   desc.add<bool>("processSimulatedTracks", true)
0618       ->setComment("return selected tracks after cutting on the floating point values");
0619   desc.add<bool>("processEmulatedTracks", true)
0620       ->setComment("return selected tracks after cutting on the bitwise emulated values");
0621   desc.add<int>("debug", 0)->setComment("Verbosity levels: 0, 1, 2, 3");
0622   descriptions.addWithDefaultLabel(desc);
0623 }
0624 
0625 //define this as a plug-in
0626 DEFINE_FWK_MODULE(L1TrackSelectionProducer);