Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-27 02:50:19

0001 // system include files
0002 #include <memory>
0003 #include <algorithm>
0004 #include <fstream>
0005 #include <cstdio>
0006 
0007 // user include files
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/stream/EDProducer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/FileInPath.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/ParameterSet/interface/allowedValues.h"
0016 
0017 #include "DataFormats/Common/interface/View.h"
0018 #include "DataFormats/Common/interface/RefToPtr.h"
0019 #include "DataFormats/L1TParticleFlow/interface/PFCandidate.h"
0020 #include "DataFormats/L1Trigger/interface/Vertex.h"
0021 #include "DataFormats/L1Trigger/interface/VertexWord.h"
0022 
0023 #include "DataFormats/Math/interface/deltaR.h"
0024 
0025 #include "L1Trigger/Phase2L1ParticleFlow/interface/l1-converters/tkinput_ref.h"
0026 #include "L1Trigger/Phase2L1ParticleFlow/interface/l1-converters/muonGmtToL1ct_ref.h"
0027 #include "L1Trigger/Phase2L1ParticleFlow/interface/l1-converters/hgcalinput_ref.h"
0028 #include "L1Trigger/Phase2L1ParticleFlow/interface/regionizer/regionizer_base_ref.h"
0029 #include "L1Trigger/Phase2L1ParticleFlow/interface/regionizer/multififo_regionizer_ref.h"
0030 #include "L1Trigger/Phase2L1ParticleFlow/interface/regionizer/buffered_folded_multififo_regionizer_ref.h"
0031 #include "L1Trigger/Phase2L1ParticleFlow/interface/regionizer/middle_buffer_multififo_regionizer_ref.h"
0032 #include "L1Trigger/Phase2L1ParticleFlow/interface/regionizer/tdr_regionizer_ref.h"
0033 #include "L1Trigger/Phase2L1ParticleFlow/interface/pf/pfalgo2hgc_ref.h"
0034 #include "L1Trigger/Phase2L1ParticleFlow/interface/pf/pfalgo3_ref.h"
0035 #include "L1Trigger/Phase2L1ParticleFlow/interface/pf/pfalgo_dummy_ref.h"
0036 #include "L1Trigger/Phase2L1ParticleFlow/interface/puppi/linpuppi_ref.h"
0037 #include "L1Trigger/Phase2L1ParticleFlow/interface/egamma/pftkegalgo_ref.h"
0038 #include "L1Trigger/Phase2L1ParticleFlow/interface/pf/pfalgo_common_ref.h"
0039 #include "L1Trigger/Phase2L1ParticleFlow/interface/egamma/pftkegsorter_ref.h"
0040 #include "L1Trigger/Phase2L1ParticleFlow/interface/egamma/pftkegsorter_barrel_ref.h"
0041 #include "L1Trigger/Phase2L1ParticleFlow/interface/L1TCorrelatorLayer1PatternFileWriter.h"
0042 
0043 #include "DataFormats/L1TCorrelator/interface/TkElectron.h"
0044 #include "DataFormats/L1TCorrelator/interface/TkElectronFwd.h"
0045 #include "DataFormats/L1Trigger/interface/EGamma.h"
0046 #include "DataFormats/L1TCorrelator/interface/TkEm.h"
0047 #include "DataFormats/L1TCorrelator/interface/TkEmFwd.h"
0048 
0049 //--------------------------------------------------------------------------------------------------
0050 class L1TCorrelatorLayer1Producer : public edm::stream::EDProducer<> {
0051 public:
0052   explicit L1TCorrelatorLayer1Producer(const edm::ParameterSet &);
0053   ~L1TCorrelatorLayer1Producer() override;
0054 
0055   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0056 
0057 private:
0058   edm::ParameterSet config_;
0059 
0060   bool hasTracks_;
0061   edm::EDGetTokenT<l1t::PFTrackCollection> tkCands_;
0062   float trkPt_;
0063   bool emuTkVtx_;
0064   edm::EDGetTokenT<std::vector<l1t::Vertex>> extTkVtx_;
0065   edm::EDGetTokenT<std::vector<l1t::VertexWord>> tkVtxEmu_;
0066   int nVtx_;
0067 
0068   edm::EDGetTokenT<l1t::SAMuonCollection> muCands_;  // standalone muons
0069 
0070   std::vector<edm::EDGetTokenT<l1t::PFClusterCollection>> emCands_;
0071   std::vector<edm::EDGetTokenT<l1t::PFClusterCollection>> hadCands_;
0072 
0073   float emPtCut_, hadPtCut_;
0074 
0075   l1ct::Event event_;
0076   std::unique_ptr<l1ct::TrackInputEmulator> trackInput_;
0077   std::unique_ptr<l1ct::GMTMuonDecoderEmulator> muonInput_;
0078   std::unique_ptr<l1ct::HgcalClusterDecoderEmulator> hgcalInput_;
0079   std::unique_ptr<l1ct::RegionizerEmulator> regionizer_;
0080   std::unique_ptr<l1ct::PFAlgoEmulatorBase> l1pfalgo_;
0081   std::unique_ptr<l1ct::LinPuppiEmulator> l1pualgo_;
0082   std::unique_ptr<l1ct::PFTkEGAlgoEmulator> l1tkegalgo_;
0083   std::unique_ptr<l1ct::PFTkEGSorterEmulator> l1tkegsorter_;
0084 
0085   // Region dump
0086   const std::string regionDumpName_;
0087   bool writeRawHgcalCluster_;
0088   std::fstream fRegionDump_;
0089   const edm::VParameterSet patternWriterConfigs_;
0090   std::vector<std::unique_ptr<L1TCorrelatorLayer1PatternFileWriter>> patternWriters_;
0091 
0092   // region of interest debugging
0093   float debugEta_, debugPhi_, debugR_;
0094 
0095   // these are used to link items back
0096   std::unordered_map<const l1t::PFCluster *, l1t::PFClusterRef> clusterRefMap_;
0097   std::unordered_map<const l1t::PFTrack *, l1t::PFTrackRef> trackRefMap_;
0098   std::unordered_map<const l1t::SAMuon *, l1t::PFCandidate::MuonRef> muonRefMap_;
0099 
0100   // main methods
0101   void beginStream(edm::StreamID) override;
0102   void endStream() override;
0103   void produce(edm::Event &, const edm::EventSetup &) override;
0104   void addUInt(unsigned int value, std::string iLabel, edm::Event &iEvent);
0105 
0106   void initSectorsAndRegions(const edm::ParameterSet &iConfig);
0107   void initEvent(const edm::Event &e);
0108   // add object, tracking references
0109   void addTrack(const l1t::PFTrack &t, l1t::PFTrackRef ref);
0110   void addMuon(const l1t::SAMuon &t, l1t::PFCandidate::MuonRef ref);
0111   void addHadCalo(const l1t::PFCluster &t, l1t::PFClusterRef ref);
0112   void addEmCalo(const l1t::PFCluster &t, l1t::PFClusterRef ref);
0113   // add objects in already-decoded format
0114   void addDecodedTrack(l1ct::DetectorSector<l1ct::TkObjEmu> &sec, const l1t::PFTrack &t);
0115   void addDecodedMuon(l1ct::DetectorSector<l1ct::MuObjEmu> &sec, const l1t::SAMuon &t);
0116   void addDecodedHadCalo(l1ct::DetectorSector<l1ct::HadCaloObjEmu> &sec, const l1t::PFCluster &t);
0117   void addDecodedEmCalo(l1ct::DetectorSector<l1ct::EmCaloObjEmu> &sec, const l1t::PFCluster &t);
0118 
0119   void addRawHgcalCluster(l1ct::DetectorSector<ap_uint<256>> &sec, const l1t::PFCluster &c);
0120 
0121   template <class T>
0122   void rawHgcalClusterEncode(ap_uint<256> &cwrd, const l1ct::DetectorSector<T> &sec, const l1t::PFCluster &c) const {
0123     cwrd = 0;
0124     ap_ufixed<14, 12, AP_RND_CONV, AP_SAT> w_pt = c.pt();
0125     ap_uint<14> w_empt = round(c.emEt() / 0.25);
0126     constexpr float ETAPHI_LSB = M_PI / 720;
0127     ap_int<9> w_eta = round(sec.region.localEta(c.eta()) / ETAPHI_LSB);
0128     ap_int<9> w_phi = round(sec.region.localPhi(c.phi()) / ETAPHI_LSB);
0129     ap_uint<10> w_qual = c.hwQual();
0130     // NOTE: this is an arbitrary choice to keep the rounding consistent with the "addDecodedHadCalo" one
0131     ap_uint<13> w_srrtot = round(c.sigmaRR() * l1ct::Scales::SRRTOT_SCALE / l1ct::Scales::SRRTOT_LSB);
0132     ap_uint<12> w_meanz = round(c.absZBarycenter());
0133     // NOTE: the calibration can actually make hoe become negative....we add a small protection for now
0134     // We use ap_ufixed to handle saturation and rounding
0135     ap_ufixed<10, 5, AP_RND_CONV, AP_SAT> w_hoe = c.hOverE();
0136 
0137     cwrd(13, 0) = w_pt.range();
0138     cwrd(27, 14) = w_empt;
0139     cwrd(72, 64) = w_eta;
0140     cwrd(81, 73) = w_phi;
0141     cwrd(115, 106) = w_qual;
0142 
0143     // FIXME: we add the variables use by composite-ID. The definitin will have to be reviewd once the
0144     // hgc format is better defined. For now we use
0145     // hwMeanZ = word 1 bits 30-19
0146     // hwSrrTot = word 3 bits 21 - 9
0147     // hoe = word 1 bits 63-52 (currently spare in the interface)
0148     cwrd(213, 201) = w_srrtot;
0149     cwrd(94, 83) = w_meanz;
0150     cwrd(127, 116) = w_hoe.range();
0151   }
0152 
0153   // fetching outputs
0154   std::unique_ptr<l1t::PFCandidateCollection> fetchHadCalo() const;
0155   std::unique_ptr<l1t::PFCandidateCollection> fetchEmCalo() const;
0156   std::unique_ptr<l1t::PFCandidateCollection> fetchTracks() const;
0157   std::unique_ptr<l1t::PFCandidateCollection> fetchPF() const;
0158   std::unique_ptr<std::vector<l1t::PFTrack>> fetchDecodedTracks() const;
0159   void putPuppi(edm::Event &iEvent) const;
0160 
0161   void putEgStaObjects(edm::Event &iEvent, const std::string &egLablel) const;
0162   void putEgObjects(edm::Event &iEvent,
0163                     const bool writeEgSta,
0164                     const std::string &tkEmLabel,
0165                     const std::string &tkEmPerBoardLabel,
0166                     const std::string &tkEleLabel,
0167                     const std::string &tkElePerBoardLabel) const;
0168 
0169   template <typename T>
0170   void setRefs_(l1t::PFCandidate &pf, const T &p) const;
0171 
0172   void doVertexings(std::vector<float> &pvdz) const;
0173   // for multiplicities
0174   enum InputType { caloType = 0, emcaloType = 1, trackType = 2, l1muType = 3 };
0175   static constexpr const char *inputTypeName[l1muType + 1] = {"Calo", "EmCalo", "TK", "Mu"};
0176   std::unique_ptr<std::vector<unsigned>> vecSecInput(InputType i) const;
0177   std::unique_ptr<std::vector<unsigned>> vecRegInput(InputType i) const;
0178   typedef l1ct::OutputRegion::ObjType OutputType;
0179   std::unique_ptr<std::vector<unsigned>> vecOutput(OutputType i, bool usePuppi) const;
0180   std::pair<unsigned int, unsigned int> totAndMax(const std::vector<unsigned> &perRegion) const;
0181 
0182   // utilities
0183   template <typename T>
0184   static edm::ParameterDescription<edm::ParameterSetDescription> getParDesc(const std::string &name) {
0185     return edm::ParameterDescription<edm::ParameterSetDescription>(
0186         name + "Parameters", T::getParameterSetDescription(), true);
0187   }
0188 };
0189 
0190 //
0191 // constructors and destructor
0192 //
0193 L1TCorrelatorLayer1Producer::L1TCorrelatorLayer1Producer(const edm::ParameterSet &iConfig)
0194     : config_(iConfig),
0195       hasTracks_(!iConfig.getParameter<edm::InputTag>("tracks").label().empty()),
0196       tkCands_(hasTracks_ ? consumes<l1t::PFTrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))
0197                           : edm::EDGetTokenT<l1t::PFTrackCollection>()),
0198       trkPt_(iConfig.getParameter<double>("trkPtCut")),
0199       muCands_(consumes<l1t::SAMuonCollection>(iConfig.getParameter<edm::InputTag>("muons"))),
0200       emPtCut_(iConfig.getParameter<double>("emPtCut")),
0201       hadPtCut_(iConfig.getParameter<double>("hadPtCut")),
0202       regionizer_(nullptr),
0203       l1pfalgo_(nullptr),
0204       l1pualgo_(nullptr),
0205       l1tkegalgo_(nullptr),
0206       l1tkegsorter_(nullptr),
0207       regionDumpName_(iConfig.getUntrackedParameter<std::string>("dumpFileName")),
0208       writeRawHgcalCluster_(iConfig.getUntrackedParameter<bool>("writeRawHgcalCluster")),
0209       patternWriterConfigs_(iConfig.getUntrackedParameter<edm::VParameterSet>("patternWriters")),
0210       debugEta_(iConfig.getUntrackedParameter<double>("debugEta")),
0211       debugPhi_(iConfig.getUntrackedParameter<double>("debugPhi")),
0212       debugR_(iConfig.getUntrackedParameter<double>("debugR")) {
0213   produces<l1t::PFCandidateCollection>("PF");
0214   produces<l1t::PFCandidateCollection>("Puppi");
0215   produces<l1t::PFCandidateRegionalOutput>("PuppiRegional");
0216 
0217   produces<l1t::PFCandidateCollection>("EmCalo");
0218   produces<l1t::PFCandidateCollection>("Calo");
0219   produces<l1t::PFCandidateCollection>("TK");
0220 #if 0  // LATER
0221   produces<l1t::PFCandidateCollection>("TKVtx");
0222 #endif
0223   produces<std::vector<l1t::PFTrack>>("DecodedTK");
0224 
0225   for (const auto &tag : iConfig.getParameter<std::vector<edm::InputTag>>("emClusters")) {
0226     emCands_.push_back(consumes<l1t::PFClusterCollection>(tag));
0227   }
0228   for (const auto &tag : iConfig.getParameter<std::vector<edm::InputTag>>("hadClusters")) {
0229     hadCands_.push_back(consumes<l1t::PFClusterCollection>(tag));
0230   }
0231 
0232   if (hasTracks_) {
0233     const std::string &tkInAlgo = iConfig.getParameter<std::string>("trackInputConversionAlgo");
0234     if (tkInAlgo == "Emulator") {
0235       trackInput_ = std::make_unique<l1ct::TrackInputEmulator>(
0236           iConfig.getParameter<edm::ParameterSet>("trackInputConversionParameters"));
0237     } else if (tkInAlgo != "Ideal")
0238       throw cms::Exception("Configuration", "Unsupported trackInputConversionAlgo");
0239   }
0240 
0241   const std::string &muInAlgo = iConfig.getParameter<std::string>("muonInputConversionAlgo");
0242   if (muInAlgo == "Emulator") {
0243     muonInput_ = std::make_unique<l1ct::GMTMuonDecoderEmulator>(
0244         iConfig.getParameter<edm::ParameterSet>("muonInputConversionParameters"));
0245   } else if (muInAlgo != "Ideal")
0246     throw cms::Exception("Configuration", "Unsupported muonInputConversionAlgo");
0247 
0248   const std::string &hgcalInAlgo = iConfig.getParameter<std::string>("hgcalInputConversionAlgo");
0249   if (hgcalInAlgo == "Emulator") {
0250     hgcalInput_ = std::make_unique<l1ct::HgcalClusterDecoderEmulator>(
0251         iConfig.getParameter<edm::ParameterSet>("hgcalInputConversionParameters"));
0252   } else if (hgcalInAlgo != "Ideal")
0253     throw cms::Exception("Configuration", "Unsupported hgcalInputConversionAlgo");
0254 
0255   const std::string &regalgo = iConfig.getParameter<std::string>("regionizerAlgo");
0256   if (regalgo == "Ideal") {
0257     regionizer_ =
0258         std::make_unique<l1ct::RegionizerEmulator>(iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters"));
0259   } else if (regalgo == "Multififo") {
0260     regionizer_ = std::make_unique<l1ct::MultififoRegionizerEmulator>(
0261         iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters"));
0262   } else if (regalgo == "BufferedFoldedMultififo") {
0263     regionizer_ = std::make_unique<l1ct::BufferedFoldedMultififoRegionizerEmulator>(
0264         iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters"));
0265   } else if (regalgo == "MultififoBarrel") {
0266     const auto &pset = iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters");
0267     regionizer_ =
0268         std::make_unique<l1ct::MultififoRegionizerEmulator>(pset.getParameter<std::string>("barrelSetup"), pset);
0269   } else if (regalgo == "MiddleBufferMultififo") {
0270     regionizer_ = std::make_unique<l1ct::MiddleBufferMultififoRegionizerEmulator>(
0271         iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters"));
0272   } else if (regalgo == "TDR") {
0273     regionizer_ = std::make_unique<l1ct::TDRRegionizerEmulator>(
0274         iConfig.getParameter<edm::ParameterSet>("regionizerAlgoParameters"));
0275   } else
0276     throw cms::Exception("Configuration", "Unsupported regionizerAlgo");
0277 
0278   const std::string &algo = iConfig.getParameter<std::string>("pfAlgo");
0279   if (algo == "PFAlgo3") {
0280     l1pfalgo_ = std::make_unique<l1ct::PFAlgo3Emulator>(iConfig.getParameter<edm::ParameterSet>("pfAlgoParameters"));
0281   } else if (algo == "PFAlgo2HGC") {
0282     l1pfalgo_ = std::make_unique<l1ct::PFAlgo2HGCEmulator>(iConfig.getParameter<edm::ParameterSet>("pfAlgoParameters"));
0283   } else if (algo == "PFAlgoDummy") {
0284     l1pfalgo_ =
0285         std::make_unique<l1ct::PFAlgoDummyEmulator>(iConfig.getParameter<edm::ParameterSet>("pfAlgoParameters"));
0286   } else
0287     throw cms::Exception("Configuration", "Unsupported pfAlgo");
0288 
0289   const std::string &pualgo = iConfig.getParameter<std::string>("puAlgo");
0290   if (pualgo == "LinearizedPuppi") {
0291     l1pualgo_ = std::make_unique<l1ct::LinPuppiEmulator>(iConfig.getParameter<edm::ParameterSet>("puAlgoParameters"));
0292   } else
0293     throw cms::Exception("Configuration", "Unsupported puAlgo");
0294 
0295   l1tkegalgo_ = std::make_unique<l1ct::PFTkEGAlgoEmulator>(
0296       l1ct::PFTkEGAlgoEmuConfig(iConfig.getParameter<edm::ParameterSet>("tkEgAlgoParameters")));
0297 
0298   const std::string &egsortalgo = iConfig.getParameter<std::string>("tkEgSorterAlgo");
0299   if (egsortalgo == "Barrel") {
0300     l1tkegsorter_ = std::make_unique<l1ct::PFTkEGSorterBarrelEmulator>(
0301         iConfig.getParameter<edm::ParameterSet>("tkEgSorterParameters"));
0302   } else if (egsortalgo == "Endcap") {
0303     l1tkegsorter_ =
0304         std::make_unique<l1ct::PFTkEGSorterEmulator>(iConfig.getParameter<edm::ParameterSet>("tkEgSorterParameters"));
0305   } else
0306     throw cms::Exception("Configuration", "Unsupported tkEgSorterAlgo");
0307 
0308   if (l1tkegalgo_->writeEgSta())
0309     produces<BXVector<l1t::EGamma>>("L1Eg");
0310   produces<l1t::TkElectronCollection>("L1TkEle");
0311   produces<l1t::TkElectronRegionalOutput>("L1TkElePerBoard");
0312   produces<l1t::TkEmCollection>("L1TkEm");
0313   produces<l1t::TkEmRegionalOutput>("L1TkEmPerBoard");
0314 
0315   emuTkVtx_ = iConfig.getParameter<bool>("vtxCollectionEmulation");
0316   if (emuTkVtx_) {
0317     tkVtxEmu_ = consumes<std::vector<l1t::VertexWord>>(iConfig.getParameter<edm::InputTag>("vtxCollection"));
0318   } else {
0319     extTkVtx_ = consumes<std::vector<l1t::Vertex>>(iConfig.getParameter<edm::InputTag>("vtxCollection"));
0320   }
0321   nVtx_ = iConfig.getParameter<int>("nVtx");
0322 
0323   const char *iprefix[4] = {"totNReg", "maxNReg", "totNSec", "maxNSec"};
0324   for (int i = 0; i <= l1muType; ++i) {
0325     for (int ip = 0; ip < 4; ++ip) {
0326       produces<unsigned int>(std::string(iprefix[ip]) + inputTypeName[i]);
0327     }
0328     produces<std::vector<unsigned>>(std::string("vecNReg") + inputTypeName[i]);
0329     produces<std::vector<unsigned>>(std::string("vecNSec") + inputTypeName[i]);
0330   }
0331   const char *oprefix[4] = {"totNPF", "maxNPF", "totNPuppi", "maxNPuppi"};
0332   for (int i = 0; i < l1ct::OutputRegion::nPFTypes; ++i) {
0333     for (int ip = 0; ip < 4; ++ip) {
0334       produces<unsigned int>(std::string(oprefix[ip]) + l1ct::OutputRegion::objTypeName[i]);
0335     }
0336     produces<std::vector<unsigned>>(std::string("vecNPF") + l1ct::OutputRegion::objTypeName[i]);
0337     produces<std::vector<unsigned>>(std::string("vecNPuppi") + l1ct::OutputRegion::objTypeName[i]);
0338   }
0339 
0340   initSectorsAndRegions(iConfig);
0341 }
0342 
0343 L1TCorrelatorLayer1Producer::~L1TCorrelatorLayer1Producer() {}
0344 
0345 void L1TCorrelatorLayer1Producer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0346   edm::ParameterSetDescription desc;
0347   // Inputs and cuts
0348   desc.add<edm::InputTag>("tracks", edm::InputTag(""));
0349   desc.add<edm::InputTag>("muons", edm::InputTag("l1tSAMuonsGmt", "prompt"));
0350   desc.add<std::vector<edm::InputTag>>("emClusters", std::vector<edm::InputTag>());
0351   desc.add<std::vector<edm::InputTag>>("hadClusters", std::vector<edm::InputTag>());
0352   desc.add<edm::InputTag>("vtxCollection", edm::InputTag("l1tVertexFinderEmulator", "L1VerticesEmulation"));
0353   desc.add<bool>("vtxCollectionEmulation", true);
0354   desc.add<double>("emPtCut", 0.0);
0355   desc.add<double>("hadPtCut", 0.0);
0356   desc.add<double>("trkPtCut", 0.0);
0357   desc.add<int32_t>("nVtx");
0358   // Input conversion
0359   edm::EmptyGroupDescription emptyGroup;
0360   desc.ifValue(edm::ParameterDescription<std::string>("trackInputConversionAlgo", "Ideal", true),
0361                "Ideal" >> emptyGroup or "Emulator" >> getParDesc<l1ct::TrackInputEmulator>("trackInputConversion"));
0362   desc.ifValue(edm::ParameterDescription<std::string>("muonInputConversionAlgo", "Ideal", true),
0363                "Ideal" >> emptyGroup or "Emulator" >> getParDesc<l1ct::GMTMuonDecoderEmulator>("muonInputConversion"));
0364   desc.ifValue(
0365       edm::ParameterDescription<std::string>("hgcalInputConversionAlgo", "Ideal", true),
0366       "Ideal" >> emptyGroup or "Emulator" >> getParDesc<l1ct::HgcalClusterDecoderEmulator>("hgcalInputConversion"));
0367   // Regionizer
0368   auto idealRegPD = getParDesc<l1ct::RegionizerEmulator>("regionizerAlgo");
0369   auto tdrRegPD = getParDesc<l1ct::TDRRegionizerEmulator>("regionizerAlgo");
0370   auto multififoRegPD = getParDesc<l1ct::MultififoRegionizerEmulator>("regionizerAlgo");
0371   auto bfMultififoRegPD = getParDesc<l1ct::BufferedFoldedMultififoRegionizerEmulator>("regionizerAlgo");
0372   auto multififoBarrelRegPD = edm::ParameterDescription<edm::ParameterSetDescription>(
0373       "regionizerAlgoParameters", l1ct::MultififoRegionizerEmulator::getParameterSetDescriptionBarrel(), true);
0374   auto mbMultififoRegPD = getParDesc<l1ct::MiddleBufferMultififoRegionizerEmulator>("regionizerAlgo");
0375   desc.ifValue(edm::ParameterDescription<std::string>("regionizerAlgo", "Ideal", true),
0376                "Ideal" >> idealRegPD or "TDR" >> tdrRegPD or "Multififo" >> multififoRegPD or
0377                    "BufferedFoldedMultififo" >> bfMultififoRegPD or "MultififoBarrel" >> multififoBarrelRegPD or
0378                    "MiddleBufferMultififo" >> mbMultififoRegPD);
0379   // PF
0380   desc.ifValue(edm::ParameterDescription<std::string>("pfAlgo", "PFAlgo3", true),
0381                "PFAlgo3" >> getParDesc<l1ct::PFAlgo3Emulator>("pfAlgo") or
0382                    "PFAlgo2HGC" >> getParDesc<l1ct::PFAlgo2HGCEmulator>("pfAlgo") or
0383                    "PFAlgoDummy" >> getParDesc<l1ct::PFAlgoDummyEmulator>("pfAlgo"));
0384   // Puppi
0385   desc.ifValue(edm::ParameterDescription<std::string>("puAlgo", "LinearizedPuppi", true),
0386                "LinearizedPuppi" >> getParDesc<l1ct::LinPuppiEmulator>("puAlgo"));
0387   // EGamma
0388   desc.add<edm::ParameterSetDescription>("tkEgAlgoParameters", l1ct::PFTkEGAlgoEmuConfig::getParameterSetDescription());
0389   // EGamma sort
0390   desc.ifValue(edm::ParameterDescription<std::string>("tkEgSorterAlgo", "Barrel", true),
0391                "Barrel" >> getParDesc<l1ct::PFTkEGSorterBarrelEmulator>("tkEgSorter") or
0392                    "Endcap" >> getParDesc<l1ct::PFTkEGSorterEmulator>("tkEgSorter"));
0393   // geometry: calo sectors
0394   edm::ParameterSetDescription caloSectorPSD;
0395   caloSectorPSD.add<std::vector<double>>("etaBoundaries");
0396   caloSectorPSD.add<uint32_t>("phiSlices", 3);
0397   caloSectorPSD.add<double>("phiZero", 0.);
0398   desc.addVPSet("caloSectors", caloSectorPSD);
0399   // geometry: regions
0400   edm::ParameterSetDescription regionPSD;
0401   regionPSD.add<std::vector<double>>("etaBoundaries");
0402   regionPSD.add<uint32_t>("phiSlices", 9);
0403   regionPSD.add<double>("etaExtra", 0.25);
0404   regionPSD.add<double>("phiExtra", 0.25);
0405   desc.addVPSet("regions", regionPSD);
0406   // geometry: boards
0407   edm::ParameterSetDescription boardPSD;
0408   boardPSD.add<std::vector<unsigned int>>("regions");
0409   desc.addVPSet("boards", boardPSD);
0410   // dump files
0411   desc.addUntracked<std::string>("dumpFileName", "");
0412   desc.addUntracked<bool>("writeRawHgcalCluster", false);
0413   // pattern files
0414   desc.addVPSetUntracked(
0415       "patternWriters", L1TCorrelatorLayer1PatternFileWriter::getParameterSetDescription(), edm::VParameterSet());
0416   // debug
0417   desc.addUntracked<double>("debugEta", 0.);
0418   desc.addUntracked<double>("debugPhi", 0.);
0419   desc.addUntracked<double>("debugR", -1.);
0420   descriptions.add("l1tCorrelatorLayer1", desc);
0421 }
0422 
0423 void L1TCorrelatorLayer1Producer::beginStream(edm::StreamID id) {
0424   if (!regionDumpName_.empty()) {
0425     if (id == 0) {
0426       fRegionDump_.open(regionDumpName_.c_str(), std::ios::out | std::ios::binary);
0427     } else {
0428       edm::LogWarning("L1TCorrelatorLayer1Producer")
0429           << "Job running with multiple streams, but dump file will have only events on stream zero.";
0430     }
0431   }
0432   if (!patternWriterConfigs_.empty()) {
0433     if (id == 0) {
0434       for (const auto &pset : patternWriterConfigs_) {
0435         patternWriters_.emplace_back(std::make_unique<L1TCorrelatorLayer1PatternFileWriter>(pset, event_));
0436       }
0437     } else {
0438       edm::LogWarning("L1TCorrelatorLayer1Producer")
0439           << "Job running with multiple streams, but pattern files will be written only on stream zero.";
0440     }
0441   }
0442 }
0443 
0444 void L1TCorrelatorLayer1Producer::endStream() {
0445   for (auto &writer : patternWriters_) {
0446     writer->flush();
0447   }
0448 }
0449 
0450 // ------------ method called to produce the data  ------------
0451 void L1TCorrelatorLayer1Producer::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0452   // clear the regions also at the beginning, in case one event didn't complete but the job continues on
0453   initEvent(iEvent);
0454 
0455   /// ------ READ TRACKS ----
0456   if (hasTracks_) {
0457     edm::Handle<l1t::PFTrackCollection> htracks;
0458     iEvent.getByToken(tkCands_, htracks);
0459     const auto &tracks = *htracks;
0460     for (unsigned int itk = 0, ntk = tracks.size(); itk < ntk; ++itk) {
0461       const auto &tk = tracks[itk];
0462       // adding objects to PF
0463       if (debugR_ > 0 && deltaR(tk.eta(), tk.phi(), debugEta_, debugPhi_) > debugR_)
0464         continue;
0465       if (tk.pt() > trkPt_) {
0466         addTrack(tk, l1t::PFTrackRef(htracks, itk));
0467       }
0468     }
0469   }
0470 
0471   /// ------ READ MUONS ----
0472   edm::Handle<l1t::SAMuonCollection> muons;
0473   iEvent.getByToken(muCands_, muons);
0474   for (unsigned int i = 0, n = muons->size(); i < n; ++i) {
0475     const l1t::SAMuon &mu = (*muons)[i];
0476     if (debugR_ > 0 && deltaR(mu.eta(), mu.phi(), debugEta_, debugPhi_) > debugR_)
0477       continue;
0478     addMuon(mu, l1t::PFCandidate::MuonRef(muons, i));
0479   }
0480   // ------ READ CALOS -----
0481   edm::Handle<l1t::PFClusterCollection> caloHandle;
0482   for (const auto &tag : emCands_) {
0483     iEvent.getByToken(tag, caloHandle);
0484     const auto &calos = *caloHandle;
0485     for (unsigned int ic = 0, nc = calos.size(); ic < nc; ++ic) {
0486       const auto &calo = calos[ic];
0487       if (debugR_ > 0 && deltaR(calo.eta(), calo.phi(), debugEta_, debugPhi_) > debugR_)
0488         continue;
0489       if (calo.pt() > emPtCut_)
0490         addEmCalo(calo, l1t::PFClusterRef(caloHandle, ic));
0491     }
0492   }
0493   for (const auto &tag : hadCands_) {
0494     iEvent.getByToken(tag, caloHandle);
0495     const auto &calos = *caloHandle;
0496     for (unsigned int ic = 0, nc = calos.size(); ic < nc; ++ic) {
0497       const auto &calo = calos[ic];
0498       if (debugR_ > 0 && deltaR(calo.eta(), calo.phi(), debugEta_, debugPhi_) > debugR_)
0499         continue;
0500       if (calo.pt() > hadPtCut_)
0501         addHadCalo(calo, l1t::PFClusterRef(caloHandle, ic));
0502     }
0503   }
0504 
0505   regionizer_->run(event_.decoded, event_.pfinputs);
0506 
0507   // First, get a copy of the discretized and corrected inputs, and write them out
0508   iEvent.put(fetchEmCalo(), "EmCalo");
0509   iEvent.put(fetchHadCalo(), "Calo");
0510   iEvent.put(fetchTracks(), "TK");
0511   iEvent.put(fetchDecodedTracks(), "DecodedTK");
0512 
0513   // Then do the vertexing, and save it out
0514   std::vector<float> z0s;
0515   std::vector<std::pair<float, float>> ptsums;
0516   float z0 = 0;
0517   double ptsum = 0;
0518   l1t::VertexWord pvwd;
0519   // FIXME: collections seem to be already sorted
0520   if (emuTkVtx_) {
0521     edm::Handle<std::vector<l1t::VertexWord>> vtxEmuHandle;
0522     iEvent.getByToken(tkVtxEmu_, vtxEmuHandle);
0523     for (const auto &vtx : *vtxEmuHandle) {
0524       ptsums.push_back(std::pair<float, float>(vtx.pt(), vtx.z0()));
0525       if (ptsum == 0 || vtx.pt() > ptsum) {
0526         ptsum = vtx.pt();
0527         pvwd = vtx;
0528       }
0529     }
0530   } else {
0531     edm::Handle<std::vector<l1t::Vertex>> vtxHandle;
0532     iEvent.getByToken(extTkVtx_, vtxHandle);
0533     for (const auto &vtx : *vtxHandle) {
0534       ptsums.push_back(std::pair<float, float>(vtx.pt(), vtx.z0()));
0535       if (ptsum == 0 || vtx.pt() > ptsum) {
0536         ptsum = vtx.pt();
0537         z0 = vtx.z0();
0538       }
0539     }
0540     pvwd = l1t::VertexWord(1, z0, 1, ptsum, 1, 1, 1);
0541   }
0542   l1ct::PVObjEmu hwpv;
0543   hwpv.hwZ0 = l1ct::Scales::makeZ0(pvwd.z0());
0544   event_.pvs.push_back(hwpv);
0545   event_.pvs_emu.push_back(pvwd.vertexWord());
0546   //get additional vertices if requested
0547   if (nVtx_ > 1) {
0548     std::stable_sort(ptsums.begin(), ptsums.end(), [](const auto &a, const auto &b) { return a.first > b.first; });
0549     for (int i0 = 0; i0 < std::min(int(ptsums.size()), int(nVtx_)); i0++) {
0550       z0s.push_back(ptsums[i0].second);
0551     }
0552     for (unsigned int i = 1; i < z0s.size(); ++i) {
0553       l1ct::PVObjEmu hwpv;
0554       hwpv.hwZ0 = l1ct::Scales::makeZ0(z0s[i]);
0555       event_.pvs.push_back(hwpv);  //Skip emu
0556     }
0557   }
0558   // Then also save the tracks with a vertex cut
0559 #if 0
0560   iEvent.put(l1regions_.fetchTracks(/*ptmin=*/0.0, /*fromPV=*/true), "TKVtx");
0561 #endif
0562 
0563   // Then run PF in each region
0564   event_.out.resize(event_.pfinputs.size());
0565   for (unsigned int ir = 0, nr = event_.pfinputs.size(); ir < nr; ++ir) {
0566     l1pfalgo_->run(event_.pfinputs[ir], event_.out[ir]);
0567     l1pfalgo_->mergeNeutrals(event_.out[ir]);
0568     l1tkegalgo_->run(event_.pfinputs[ir], event_.out[ir]);
0569     l1tkegalgo_->runIso(event_.pfinputs[ir], event_.pvs, event_.out[ir]);
0570   }
0571 
0572   // Then run puppi (regionally)
0573   for (unsigned int ir = 0, nr = event_.pfinputs.size(); ir < nr; ++ir) {
0574     l1pualgo_->run(event_.pfinputs[ir], event_.pvs, event_.out[ir]);
0575     //l1pualgo_->runNeutralsPU(l1region, z0, -1., puGlobals);
0576   }
0577 
0578   // NOTE: This needs to happen before the EG sorting per board so that the EG objects
0579   // get a global reference to the EGSta before being mixed among differente regions
0580   std::vector<edm::Ref<BXVector<l1t::EGamma>>> egsta_refs;
0581   if (l1tkegalgo_->writeEgSta()) {
0582     putEgStaObjects(iEvent, "L1Eg");
0583   }
0584 
0585   // l1tkegsorter_->setDebug(true);
0586   for (auto &board : event_.board_out) {
0587     l1tkegsorter_->runPho(event_.pfinputs, event_.out, board.region_index, board.egphoton);
0588     l1tkegsorter_->runEle(event_.pfinputs, event_.out, board.region_index, board.egelectron);
0589   }
0590 
0591   // save PF into the event
0592   iEvent.put(fetchPF(), "PF");
0593 
0594   // and save puppi
0595   putPuppi(iEvent);
0596 
0597   // save the EG objects
0598   putEgObjects(iEvent, l1tkegalgo_->writeEgSta(), "L1TkEm", "L1TkEmPerBoard", "L1TkEle", "L1TkElePerBoard");
0599 
0600   // Then go do the multiplicities
0601   for (int i = 0; i <= l1muType; ++i) {
0602     auto vecInputs = vecSecInput(InputType(i));
0603     auto tm = totAndMax(*vecInputs);
0604     addUInt(tm.first, std::string("totNSec") + inputTypeName[i], iEvent);
0605     addUInt(tm.second, std::string("maxNSec") + inputTypeName[i], iEvent);
0606     iEvent.put(std::move(vecInputs), std::string("vecNSec") + inputTypeName[i]);
0607   }
0608   for (int i = 0; i <= l1muType; ++i) {
0609     auto vecInputs = vecRegInput(InputType(i));
0610     auto tm = totAndMax(*vecInputs);
0611     addUInt(tm.first, std::string("totNReg") + inputTypeName[i], iEvent);
0612     addUInt(tm.second, std::string("maxNReg") + inputTypeName[i], iEvent);
0613     iEvent.put(std::move(vecInputs), std::string("vecNReg") + inputTypeName[i]);
0614   }
0615   for (int i = 0; i < l1ct::OutputRegion::nPFTypes; ++i) {
0616     auto vecPF = vecOutput(OutputType(i), false);
0617     auto tmPF = totAndMax(*vecPF);
0618     addUInt(tmPF.first, std::string("totNPF") + l1ct::OutputRegion::objTypeName[i], iEvent);
0619     addUInt(tmPF.second, std::string("maxNPF") + l1ct::OutputRegion::objTypeName[i], iEvent);
0620     iEvent.put(std::move(vecPF), std::string("vecNPF") + l1ct::OutputRegion::objTypeName[i]);
0621     auto vecPuppi = vecOutput(OutputType(i), true);
0622     auto tmPuppi = totAndMax(*vecPuppi);
0623     addUInt(tmPuppi.first, std::string("totNPuppi") + l1ct::OutputRegion::objTypeName[i], iEvent);
0624     addUInt(tmPuppi.second, std::string("maxNPuppi") + l1ct::OutputRegion::objTypeName[i], iEvent);
0625     iEvent.put(std::move(vecPuppi), std::string("vecNPuppi") + l1ct::OutputRegion::objTypeName[i]);
0626   }
0627 
0628   if (fRegionDump_.is_open()) {
0629     event_.write(fRegionDump_);
0630   }
0631   for (auto &writer : patternWriters_) {
0632     writer->write(event_);
0633   }
0634 
0635   // finally clear the regions
0636   event_.clear();
0637 }
0638 
0639 void L1TCorrelatorLayer1Producer::addUInt(unsigned int value, std::string iLabel, edm::Event &iEvent) {
0640   iEvent.put(std::make_unique<unsigned>(value), iLabel);
0641 }
0642 
0643 void L1TCorrelatorLayer1Producer::initSectorsAndRegions(const edm::ParameterSet &iConfig) {
0644   // the track finder geometry is fixed
0645   unsigned int TF_phiSlices = 9;
0646   float TF_phiWidth = 2 * M_PI / TF_phiSlices;
0647   event_.decoded.track.clear();
0648   for (unsigned int ieta = 0, neta = 2; ieta < neta; ++ieta) {
0649     for (unsigned int iphi = 0; iphi < TF_phiSlices; ++iphi) {
0650       float phiCenter = reco::reducePhiRange(iphi * TF_phiWidth);
0651       event_.decoded.track.emplace_back((ieta ? 0. : -2.5), (ieta ? 2.5 : 0.0), phiCenter, TF_phiWidth);
0652       event_.raw.track.emplace_back((ieta ? 0. : -2.5), (ieta ? 2.5 : 0.0), phiCenter, TF_phiWidth);
0653     }
0654   }
0655 
0656   event_.decoded.emcalo.clear();
0657   event_.decoded.hadcalo.clear();
0658   event_.raw.hgcalcluster.clear();
0659 
0660   for (const edm::ParameterSet &preg : iConfig.getParameter<edm::VParameterSet>("caloSectors")) {
0661     std::vector<double> etaBoundaries = preg.getParameter<std::vector<double>>("etaBoundaries");
0662     if (!std::is_sorted(etaBoundaries.begin(), etaBoundaries.end()))
0663       throw cms::Exception("Configuration", "caloSectors.etaBoundaries not sorted\n");
0664     unsigned int phiSlices = preg.getParameter<uint32_t>("phiSlices");
0665     float phiWidth = 2 * M_PI / phiSlices;
0666     if (phiWidth > 2 * l1ct::Scales::maxAbsPhi())
0667       throw cms::Exception("Configuration", "caloSectors phi range too large for phi_t data type");
0668     double phiZero = preg.getParameter<double>("phiZero");
0669     for (unsigned int ieta = 0, neta = etaBoundaries.size() - 1; ieta < neta; ++ieta) {
0670       float etaWidth = etaBoundaries[ieta + 1] - etaBoundaries[ieta];
0671       if (etaWidth > 2 * l1ct::Scales::maxAbsEta())
0672         throw cms::Exception("Configuration", "caloSectors eta range too large for eta_t data type");
0673       for (unsigned int iphi = 0; iphi < phiSlices; ++iphi) {
0674         float phiCenter = reco::reducePhiRange(iphi * phiWidth + phiZero);
0675         event_.decoded.hadcalo.emplace_back(etaBoundaries[ieta], etaBoundaries[ieta + 1], phiCenter, phiWidth);
0676         event_.decoded.emcalo.emplace_back(etaBoundaries[ieta], etaBoundaries[ieta + 1], phiCenter, phiWidth);
0677         event_.raw.hgcalcluster.emplace_back(etaBoundaries[ieta], etaBoundaries[ieta + 1], phiCenter, phiWidth);
0678       }
0679     }
0680   }
0681 
0682   event_.decoded.muon.region = l1ct::PFRegionEmu(0., 0.);  // centered at (0,0)
0683   event_.raw.muon.region = l1ct::PFRegionEmu(0., 0.);      // centered at (0,0)
0684 
0685   event_.pfinputs.clear();
0686   for (const edm::ParameterSet &preg : iConfig.getParameter<edm::VParameterSet>("regions")) {
0687     std::vector<double> etaBoundaries = preg.getParameter<std::vector<double>>("etaBoundaries");
0688     if (!std::is_sorted(etaBoundaries.begin(), etaBoundaries.end()))
0689       throw cms::Exception("Configuration", "regions.etaBoundaries not sorted\n");
0690     unsigned int phiSlices = preg.getParameter<uint32_t>("phiSlices");
0691     float etaExtra = preg.getParameter<double>("etaExtra");
0692     float phiExtra = preg.getParameter<double>("phiExtra");
0693     float phiWidth = 2 * M_PI / phiSlices;
0694     for (unsigned int ieta = 0, neta = etaBoundaries.size() - 1; ieta < neta; ++ieta) {
0695       for (unsigned int iphi = 0; iphi < phiSlices; ++iphi) {
0696         float phiCenter = reco::reducePhiRange(iphi * phiWidth);  //align with L1 TrackFinder phi sector indexing
0697         event_.pfinputs.emplace_back(
0698             etaBoundaries[ieta], etaBoundaries[ieta + 1], phiCenter, phiWidth, etaExtra, phiExtra);
0699       }
0700     }
0701   }
0702 
0703   event_.board_out.clear();
0704   const edm::VParameterSet &board_params = iConfig.getParameter<edm::VParameterSet>("boards");
0705   event_.board_out.resize(board_params.size());
0706   for (unsigned int bidx = 0; bidx < board_params.size(); bidx++) {
0707     event_.board_out[bidx].region_index = board_params[bidx].getParameter<std::vector<unsigned int>>("regions");
0708     float etaBoard = 0.;
0709     float phiBoard = 0.;
0710     for (auto ridx : event_.board_out[bidx].region_index) {
0711       etaBoard += event_.pfinputs[ridx].region.floatEtaCenter();
0712       phiBoard += event_.pfinputs[ridx].region.floatPhiCenter();
0713     }
0714     event_.board_out[bidx].eta = etaBoard / event_.board_out[bidx].region_index.size();
0715     event_.board_out[bidx].phi = phiBoard / event_.board_out[bidx].region_index.size();
0716   }
0717 }
0718 
0719 void L1TCorrelatorLayer1Producer::initEvent(const edm::Event &iEvent) {
0720   event_.clear();
0721   event_.run = iEvent.id().run();
0722   event_.lumi = iEvent.id().luminosityBlock();
0723   event_.event = iEvent.id().event();
0724   clusterRefMap_.clear();
0725   trackRefMap_.clear();
0726   muonRefMap_.clear();
0727 }
0728 
0729 void L1TCorrelatorLayer1Producer::addTrack(const l1t::PFTrack &t, l1t::PFTrackRef ref) {
0730   auto &rawsectors = event_.raw.track;
0731   auto &sectors = event_.decoded.track;
0732   assert(sectors.size() == 18 && rawsectors.size() == 18);
0733   int isec = t.track()->phiSector() + (t.eta() >= 0 ? 9 : 0);
0734   rawsectors[isec].obj.push_back(t.trackWord().getTrackWord());
0735   addDecodedTrack(sectors[isec], t);
0736   trackRefMap_[&t] = ref;
0737 }
0738 void L1TCorrelatorLayer1Producer::addMuon(const l1t::SAMuon &mu, l1t::PFCandidate::MuonRef ref) {
0739   event_.raw.muon.obj.emplace_back(mu.word());
0740   addDecodedMuon(event_.decoded.muon, mu);
0741   muonRefMap_[&mu] = ref;
0742 }
0743 void L1TCorrelatorLayer1Producer::addHadCalo(const l1t::PFCluster &c, l1t::PFClusterRef ref) {
0744   int sidx = 0;
0745   for (auto &sec : event_.decoded.hadcalo) {
0746     if (sec.region.contains(c.eta(), c.phi())) {
0747       addDecodedHadCalo(sec, c);
0748       if (writeRawHgcalCluster_)
0749         addRawHgcalCluster(event_.raw.hgcalcluster[sidx], c);
0750     }
0751     sidx++;
0752   }
0753   clusterRefMap_[&c] = ref;
0754 }
0755 void L1TCorrelatorLayer1Producer::addEmCalo(const l1t::PFCluster &c, l1t::PFClusterRef ref) {
0756   for (auto &sec : event_.decoded.emcalo) {
0757     if (sec.region.contains(c.eta(), c.phi())) {
0758       addDecodedEmCalo(sec, c);
0759     }
0760   }
0761   clusterRefMap_[&c] = ref;
0762 }
0763 
0764 void L1TCorrelatorLayer1Producer::addDecodedTrack(l1ct::DetectorSector<l1ct::TkObjEmu> &sec, const l1t::PFTrack &t) {
0765   std::pair<l1ct::TkObjEmu, bool> tkAndSel;
0766   if (trackInput_) {
0767     tkAndSel = trackInput_->decodeTrack(t.trackWord().getTrackWord(), sec.region);
0768   } else {
0769     tkAndSel.first.hwPt = l1ct::Scales::makePtFromFloat(t.pt());
0770     tkAndSel.first.hwEta =
0771         l1ct::Scales::makeGlbEta(t.caloEta()) -
0772         sec.region.hwEtaCenter;  // important to enforce that the region boundary is on a discrete value
0773     tkAndSel.first.hwPhi = l1ct::Scales::makePhi(sec.region.localPhi(t.caloPhi()));
0774     tkAndSel.first.hwCharge = t.charge() > 0;
0775     tkAndSel.first.hwQuality = t.quality();
0776     tkAndSel.first.hwDEta = l1ct::Scales::makeEta(t.eta() - t.caloEta());
0777     tkAndSel.first.hwDPhi = l1ct::Scales::makePhi(std::abs(reco::deltaPhi(t.phi(), t.caloPhi())));
0778     tkAndSel.first.hwZ0 = l1ct::Scales::makeZ0(t.vertex().Z());
0779     tkAndSel.first.hwDxy = 0;
0780     tkAndSel.second = t.quality() > 0;
0781   }
0782   // CMSSW-only extra info
0783   tkAndSel.first.hwChi2 = round(t.chi2() * 10);
0784   tkAndSel.first.simPt = t.pt();
0785   tkAndSel.first.simCaloEta = t.caloEta();
0786   tkAndSel.first.simCaloPhi = t.caloPhi();
0787   tkAndSel.first.simVtxEta = t.eta();
0788   tkAndSel.first.simVtxPhi = t.phi();
0789   tkAndSel.first.simZ0 = t.vertex().Z();
0790   tkAndSel.first.simD0 = t.vertex().Rho();
0791   tkAndSel.first.src = &t;
0792 
0793   // If the track fails, we set its pT to zero, so that the decoded tracks are still aligned with the raw tracks
0794   // Downstream, the regionizer will just ignore zero-momentum tracks
0795   if (!tkAndSel.second)
0796     tkAndSel.first.hwPt = 0;
0797   sec.obj.push_back(tkAndSel.first);
0798 }
0799 
0800 void L1TCorrelatorLayer1Producer::addDecodedMuon(l1ct::DetectorSector<l1ct::MuObjEmu> &sec, const l1t::SAMuon &t) {
0801   l1ct::MuObjEmu mu;
0802   if (muonInput_) {
0803     mu = muonInput_->decode(t.word());
0804   } else {
0805     mu.hwPt = l1ct::Scales::makePtFromFloat(t.pt());
0806     mu.hwEta = l1ct::Scales::makeGlbEta(t.eta());  // IMPORTANT: input is in global coordinates!
0807     mu.hwPhi = l1ct::Scales::makeGlbPhi(t.phi());
0808     mu.hwCharge = !t.hwCharge();
0809     mu.hwQuality = t.hwQual() / 2;
0810     mu.hwDEta = 0;
0811     mu.hwDPhi = 0;
0812     mu.hwZ0 = l1ct::Scales::makeZ0(t.vertex().Z());
0813     mu.hwDxy = 0;  // Dxy not defined yet
0814   }
0815   mu.src = &t;
0816   sec.obj.push_back(mu);
0817 }
0818 
0819 void L1TCorrelatorLayer1Producer::addDecodedHadCalo(l1ct::DetectorSector<l1ct::HadCaloObjEmu> &sec,
0820                                                     const l1t::PFCluster &c) {
0821   l1ct::HadCaloObjEmu calo;
0822   ap_uint<256> word = 0;
0823   rawHgcalClusterEncode(word, sec, c);
0824   if (hgcalInput_) {
0825     calo = hgcalInput_->decode(word);
0826   } else {
0827     calo.hwPt = l1ct::Scales::makePtFromFloat(c.pt());
0828     calo.hwEta = l1ct::Scales::makeGlbEta(c.eta()) -
0829                  sec.region.hwEtaCenter;  // important to enforce that the region boundary is on a discrete value
0830     calo.hwPhi = l1ct::Scales::makePhi(sec.region.localPhi(c.phi()));
0831     calo.hwEmPt = l1ct::Scales::makePtFromFloat(c.emEt());
0832     calo.hwEmID = c.hwEmID();
0833     calo.hwSrrTot = l1ct::Scales::makeSrrTot(c.sigmaRR());
0834     calo.hwMeanZ = c.absZBarycenter() < 320. ? l1ct::meanz_t(0) : l1ct::Scales::makeMeanZ(c.absZBarycenter());
0835     calo.hwHoe = l1ct::Scales::makeHoe(c.hOverE());
0836   }
0837   calo.src = &c;
0838   sec.obj.push_back(calo);
0839 }
0840 
0841 void L1TCorrelatorLayer1Producer::addRawHgcalCluster(l1ct::DetectorSector<ap_uint<256>> &sec, const l1t::PFCluster &c) {
0842   ap_uint<256> cwrd = 0;
0843   rawHgcalClusterEncode(cwrd, sec, c);
0844   sec.obj.push_back(cwrd);
0845 }
0846 
0847 void L1TCorrelatorLayer1Producer::addDecodedEmCalo(l1ct::DetectorSector<l1ct::EmCaloObjEmu> &sec,
0848                                                    const l1t::PFCluster &c) {
0849   l1ct::EmCaloObjEmu calo;
0850   // set the endcap-sepcific variables to default value:
0851   calo.clear();
0852   calo.hwPt = l1ct::Scales::makePtFromFloat(c.pt());
0853   calo.hwEta = l1ct::Scales::makeGlbEta(c.eta()) -
0854                sec.region.hwEtaCenter;  // important to enforce that the region boundary is on a discrete value
0855   calo.hwPhi = l1ct::Scales::makePhi(sec.region.localPhi(c.phi()));
0856   calo.hwPtErr = l1ct::Scales::makePtFromFloat(c.ptError());
0857   calo.hwEmID = c.hwEmID();
0858   calo.src = &c;
0859   sec.obj.push_back(calo);
0860 }
0861 
0862 template <typename T>
0863 void L1TCorrelatorLayer1Producer::setRefs_(l1t::PFCandidate &pf, const T &p) const {
0864   if (p.srcCluster) {
0865     auto match = clusterRefMap_.find(p.srcCluster);
0866     if (match == clusterRefMap_.end()) {
0867       throw cms::Exception("CorruptData") << "Invalid cluster pointer in PF candidate id " << p.intId() << " pt "
0868                                           << p.floatPt() << " eta " << p.floatEta() << " phi " << p.floatPhi();
0869     }
0870     pf.setPFCluster(match->second);
0871   }
0872   if (p.srcTrack) {
0873     auto match = trackRefMap_.find(p.srcTrack);
0874     if (match == trackRefMap_.end()) {
0875       throw cms::Exception("CorruptData") << "Invalid track pointer in PF candidate id " << p.intId() << " pt "
0876                                           << p.floatPt() << " eta " << p.floatEta() << " phi " << p.floatPhi();
0877     }
0878     pf.setPFTrack(match->second);
0879   }
0880   if (p.srcMu) {
0881     auto match = muonRefMap_.find(p.srcMu);
0882     if (match == muonRefMap_.end()) {
0883       throw cms::Exception("CorruptData") << "Invalid muon pointer in PF candidate id " << p.intId() << " pt "
0884                                           << p.floatPt() << " eta " << p.floatEta() << " phi " << p.floatPhi();
0885     }
0886     pf.setMuon(match->second);
0887   }
0888 }
0889 
0890 template <>
0891 void L1TCorrelatorLayer1Producer::setRefs_<l1ct::PFNeutralObjEmu>(l1t::PFCandidate &pf,
0892                                                                   const l1ct::PFNeutralObjEmu &p) const {
0893   if (p.srcCluster) {
0894     auto match = clusterRefMap_.find(p.srcCluster);
0895     if (match == clusterRefMap_.end()) {
0896       throw cms::Exception("CorruptData") << "Invalid cluster pointer in PF candidate id " << p.intId() << " pt "
0897                                           << p.floatPt() << " eta " << p.floatEta() << " phi " << p.floatPhi();
0898     }
0899     pf.setPFCluster(match->second);
0900   }
0901 }
0902 
0903 template <>
0904 void L1TCorrelatorLayer1Producer::setRefs_<l1ct::HadCaloObjEmu>(l1t::PFCandidate &pf,
0905                                                                 const l1ct::HadCaloObjEmu &p) const {
0906   if (p.src) {
0907     auto match = clusterRefMap_.find(p.src);
0908     if (match == clusterRefMap_.end()) {
0909       throw cms::Exception("CorruptData") << "Invalid cluster pointer in hadcalo candidate  pt " << p.floatPt()
0910                                           << " eta " << p.floatEta() << " phi " << p.floatPhi();
0911     }
0912     pf.setPFCluster(match->second);
0913   }
0914 }
0915 
0916 template <>
0917 void L1TCorrelatorLayer1Producer::setRefs_<l1ct::EmCaloObjEmu>(l1t::PFCandidate &pf,
0918                                                                const l1ct::EmCaloObjEmu &p) const {
0919   if (p.src) {
0920     auto match = clusterRefMap_.find(p.src);
0921     if (match == clusterRefMap_.end()) {
0922       throw cms::Exception("CorruptData") << "Invalid cluster pointer in emcalo candidate  pt " << p.floatPt()
0923                                           << " eta " << p.floatEta() << " phi " << p.floatPhi();
0924     }
0925     pf.setPFCluster(match->second);
0926   }
0927 }
0928 
0929 template <>
0930 void L1TCorrelatorLayer1Producer::setRefs_<l1ct::TkObjEmu>(l1t::PFCandidate &pf, const l1ct::TkObjEmu &p) const {
0931   if (p.src) {
0932     auto match = trackRefMap_.find(p.src);
0933     if (match == trackRefMap_.end()) {
0934       throw cms::Exception("CorruptData") << "Invalid track pointer in track candidate  pt " << p.floatPt() << " eta "
0935                                           << p.floatEta() << " phi " << p.floatPhi();
0936     }
0937     pf.setPFTrack(match->second);
0938   }
0939 }
0940 
0941 std::unique_ptr<l1t::PFCandidateCollection> L1TCorrelatorLayer1Producer::fetchHadCalo() const {
0942   auto ret = std::make_unique<l1t::PFCandidateCollection>();
0943   for (const auto &r : event_.pfinputs) {
0944     const auto &reg = r.region;
0945     for (const auto &p : r.hadcalo) {
0946       if (p.hwPt == 0 || !reg.isFiducial(p))
0947         continue;
0948       reco::Particle::PolarLorentzVector p4(p.floatPt(), reg.floatGlbEtaOf(p), reg.floatGlbPhiOf(p), 0.13f);
0949       l1t::PFCandidate::ParticleType type = p.hwIsEM() ? l1t::PFCandidate::Photon : l1t::PFCandidate::NeutralHadron;
0950       ret->emplace_back(type, 0, p4, 1, p.intPt(), p.intEta(), p.intPhi());
0951       ret->back().setHwEmID(p.hwEmID);
0952       setRefs_(ret->back(), p);
0953     }
0954   }
0955   return ret;
0956 }
0957 std::unique_ptr<l1t::PFCandidateCollection> L1TCorrelatorLayer1Producer::fetchEmCalo() const {
0958   auto ret = std::make_unique<l1t::PFCandidateCollection>();
0959   for (const auto &r : event_.pfinputs) {
0960     const auto &reg = r.region;
0961     for (const auto &p : r.emcalo) {
0962       if (p.hwPt == 0 || !reg.isFiducial(p))
0963         continue;
0964       reco::Particle::PolarLorentzVector p4(p.floatPt(), reg.floatGlbEtaOf(p), reg.floatGlbPhiOf(p), 0.13f);
0965       ret->emplace_back(l1t::PFCandidate::Photon, 0, p4, 1, p.intPt(), p.intEta(), p.intPhi());
0966       ret->back().setHwEmID(p.hwEmID);
0967       setRefs_(ret->back(), p);
0968     }
0969   }
0970   return ret;
0971 }
0972 std::unique_ptr<l1t::PFCandidateCollection> L1TCorrelatorLayer1Producer::fetchTracks() const {
0973   auto ret = std::make_unique<l1t::PFCandidateCollection>();
0974   for (const auto &r : event_.pfinputs) {
0975     const auto &reg = r.region;
0976     for (const auto &p : r.track) {
0977       if (p.hwPt == 0 || !reg.isFiducial(p))
0978         continue;
0979       reco::Particle::PolarLorentzVector p4(
0980           p.floatPt(), reg.floatGlbEta(p.hwVtxEta()), reg.floatGlbPhi(p.hwVtxPhi()), 0.13f);
0981       ret->emplace_back(l1t::PFCandidate::ChargedHadron, p.intCharge(), p4, 1, p.intPt(), p.intEta(), p.intPhi());
0982       setRefs_(ret->back(), p);
0983     }
0984   }
0985   return ret;
0986 }
0987 
0988 std::unique_ptr<std::vector<l1t::PFTrack>> L1TCorrelatorLayer1Producer::fetchDecodedTracks() const {
0989   auto ret = std::make_unique<std::vector<l1t::PFTrack>>();
0990   for (const auto &r : event_.decoded.track) {
0991     const auto &reg = r.region;
0992     for (const auto &p : r.obj) {
0993       if (p.hwPt == 0 || !reg.isFiducial(p))
0994         continue;
0995       reco::Particle::PolarLorentzVector p4(
0996           p.floatPt(), reg.floatGlbEta(p.hwVtxEta()), reg.floatGlbPhi(p.hwVtxPhi()), 0);
0997 
0998       reco::Particle::Point vtx(0, 0, p.floatZ0());
0999 
1000       ret->emplace_back(l1t::PFTrack(p.intCharge(),
1001                                      reco::Particle::LorentzVector(p4),
1002                                      vtx,
1003                                      p.src->track(),
1004                                      0,
1005                                      reg.floatGlbEta(p.hwEta),
1006                                      reg.floatGlbPhi(p.hwPhi),
1007                                      -1,
1008                                      -1,
1009                                      p.hwQuality.to_int(),
1010                                      false,
1011                                      p.intPt(),
1012                                      p.intEta(),
1013                                      p.intPhi()));
1014     }
1015   }
1016   return ret;
1017 }
1018 
1019 std::unique_ptr<l1t::PFCandidateCollection> L1TCorrelatorLayer1Producer::fetchPF() const {
1020   auto ret = std::make_unique<l1t::PFCandidateCollection>();
1021   for (unsigned int ir = 0, nr = event_.pfinputs.size(); ir < nr; ++ir) {
1022     const auto &reg = event_.pfinputs[ir].region;
1023     for (const auto &p : event_.out[ir].pfcharged) {
1024       if (p.hwPt == 0 || !reg.isFiducial(p))
1025         continue;
1026       reco::Particle::PolarLorentzVector p4(
1027           p.floatPt(), reg.floatGlbEta(p.hwVtxEta()), reg.floatGlbPhi(p.hwVtxPhi()), 0.13f);
1028       l1t::PFCandidate::ParticleType type = l1t::PFCandidate::ChargedHadron;
1029       if (p.hwId.isMuon())
1030         type = l1t::PFCandidate::Muon;
1031       else if (p.hwId.isElectron())
1032         type = l1t::PFCandidate::Electron;
1033       ret->emplace_back(type, p.intCharge(), p4, 1, p.intPt(), p.intEta(), p.intPhi());
1034       ret->back().setZ0(p.floatZ0());
1035       ret->back().setDxy(p.floatDxy());
1036       ret->back().setHwZ0(p.hwZ0);
1037       ret->back().setHwDxy(p.hwDxy);
1038       ret->back().setHwTkQuality(p.hwTkQuality);
1039       ret->back().setCaloEta(reg.floatGlbEtaOf(p));
1040       ret->back().setCaloPhi(reg.floatGlbPhiOf(p));
1041 
1042       setRefs_(ret->back(), p);
1043     }
1044     for (const auto &p : event_.out[ir].pfneutral) {
1045       if (p.hwPt == 0 || !reg.isFiducial(p))
1046         continue;
1047       reco::Particle::PolarLorentzVector p4(p.floatPt(), reg.floatGlbEtaOf(p), reg.floatGlbPhiOf(p), 0.13f);
1048       l1t::PFCandidate::ParticleType type =
1049           p.hwId.isPhoton() ? l1t::PFCandidate::Photon : l1t::PFCandidate::NeutralHadron;
1050       ret->emplace_back(type, 0, p4, 1, p.intPt(), p.intEta(), p.intPhi());
1051       ret->back().setHwEmID(p.hwEmID);
1052       ret->back().setCaloEta(reg.floatGlbEtaOf(p));
1053       ret->back().setCaloPhi(reg.floatGlbPhiOf(p));
1054       setRefs_(ret->back(), p);
1055     }
1056   }
1057   return ret;
1058 }
1059 
1060 void L1TCorrelatorLayer1Producer::putPuppi(edm::Event &iEvent) const {
1061   auto refprod = iEvent.getRefBeforePut<l1t::PFCandidateCollection>("Puppi");
1062   auto coll = std::make_unique<l1t::PFCandidateCollection>();
1063   auto reg = std::make_unique<l1t::PFCandidateRegionalOutput>(refprod);
1064   std::vector<int> nobj;
1065   for (unsigned int ir = 0, nr = event_.pfinputs.size(); ir < nr; ++ir) {
1066     nobj.clear();
1067     for (const auto &p : event_.out[ir].puppi) {
1068       if (p.hwPt == 0)
1069         continue;
1070       // note: Puppi candidates are already in global coordinates & fiducial-only!
1071       l1t::PFCandidate::ParticleType type;
1072       float mass = 0.13f;
1073       if (p.hwId.charged()) {
1074         if (p.hwId.isMuon()) {
1075           type = l1t::PFCandidate::Muon;
1076           mass = 0.105;
1077         } else if (p.hwId.isElectron()) {
1078           type = l1t::PFCandidate::Electron;
1079           mass = 0.005;
1080         } else
1081           type = l1t::PFCandidate::ChargedHadron;
1082       } else {
1083         type = p.hwId.isPhoton() ? l1t::PFCandidate::Photon : l1t::PFCandidate::NeutralHadron;
1084         mass = p.hwId.isPhoton() ? 0.0 : 0.5;
1085       }
1086       reco::Particle::PolarLorentzVector p4(p.floatPt(), p.floatEta(), p.floatPhi(), mass);
1087       coll->emplace_back(type, p.intCharge(), p4, p.floatPuppiW(), p.intPt(), p.intEta(), p.intPhi());
1088       if (p.hwId.charged()) {
1089         coll->back().setZ0(p.floatZ0());
1090         coll->back().setDxy(p.floatDxy());
1091         coll->back().setHwZ0(p.hwZ0());
1092         coll->back().setHwDxy(p.hwDxy());
1093         coll->back().setHwTkQuality(p.hwTkQuality());
1094       } else {
1095         coll->back().setHwPuppiWeight(p.hwPuppiW());
1096         coll->back().setHwEmID(p.hwEmID());
1097       }
1098       coll->back().setEncodedPuppi64(p.pack().to_uint64());
1099       setRefs_(coll->back(), p);
1100       nobj.push_back(coll->size() - 1);
1101     }
1102     reg->addRegion(nobj, event_.pfinputs[ir].region.floatEtaCenter(), event_.pfinputs[ir].region.floatPhiCenter());
1103   }
1104   iEvent.put(std::move(coll), "Puppi");
1105   iEvent.put(std::move(reg), "PuppiRegional");
1106 }
1107 
1108 void L1TCorrelatorLayer1Producer::putEgStaObjects(edm::Event &iEvent, const std::string &egLablel) const {
1109   auto egs = std::make_unique<BXVector<l1t::EGamma>>();
1110   // FIXME: in case more BXes are introduced shuld probably use egs->key(egs->end(bx));
1111 
1112   for (unsigned int ir = 0, nr = event_.pfinputs.size(); ir < nr; ++ir) {
1113     const auto &reg = event_.pfinputs[ir].region;
1114 
1115     // EG standalone objects
1116     for (unsigned int ieg = 0, neg = event_.out[ir].egsta.size(); ieg < neg; ++ieg) {
1117       const auto &p = event_.out[ir].egsta[ieg];
1118       if (p.hwPt == 0 || !reg.isFiducial(p))
1119         continue;
1120       l1t::EGamma eg(
1121           reco::Candidate::PolarLorentzVector(p.floatPt(), reg.floatGlbEta(p.hwEta), reg.floatGlbPhi(p.hwPhi), 0.));
1122       eg.setHwQual(p.hwQual);
1123       egs->push_back(0, eg);
1124     }
1125   }
1126 
1127   iEvent.put(std::move(egs), egLablel);
1128 }
1129 
1130 void L1TCorrelatorLayer1Producer::putEgObjects(edm::Event &iEvent,
1131                                                const bool writeEgSta,
1132                                                const std::string &tkEmLabel,
1133                                                const std::string &tkEmPerBoardLabel,
1134                                                const std::string &tkEleLabel,
1135                                                const std::string &tkElePerBoardLabel) const {
1136   auto tkems = std::make_unique<l1t::TkEmCollection>();
1137   auto tkemRefProd = iEvent.getRefBeforePut<l1t::TkEmCollection>(tkEmLabel);
1138   auto tkemPerBoard = std::make_unique<l1t::TkEmRegionalOutput>(tkemRefProd);
1139   auto tkeles = std::make_unique<l1t::TkElectronCollection>();
1140   auto tkeleRefProd = iEvent.getRefBeforePut<l1t::TkElectronCollection>(tkEleLabel);
1141   auto tkelePerBoard = std::make_unique<l1t::TkElectronRegionalOutput>(tkeleRefProd);
1142 
1143   // TkEG objects are written out after the per-board sorting.
1144   // The mapping to each board is saved into the regionalmap for further (stage-2 consumption)
1145   std::vector<int> nele_obj;
1146   std::vector<int> npho_obj;
1147 
1148   for (const auto &board : event_.board_out) {
1149     npho_obj.clear();
1150     for (const auto &egiso : board.egphoton) {
1151       if (egiso.hwPt == 0)
1152         continue;
1153 
1154       reco::Candidate::PolarLorentzVector mom(egiso.floatPt(), egiso.floatEta(), egiso.floatPhi(), 0.);
1155 
1156       l1t::TkEm tkem(reco::Candidate::LorentzVector(mom),
1157                      egiso.srcCluster->constituentsAndFractions()[0].first,
1158                      egiso.floatRelIso(l1ct::EGIsoObjEmu::IsoType::TkIso),
1159                      egiso.floatRelIso(l1ct::EGIsoObjEmu::IsoType::TkIsoPV));
1160       tkem.setHwQual(egiso.hwQual);
1161       tkem.setPFIsol(egiso.floatRelIso(l1ct::EGIsoObjEmu::IsoType::PfIso));
1162       tkem.setPFIsolPV(egiso.floatRelIso(l1ct::EGIsoObjEmu::IsoType::PfIsoPV));
1163       tkem.setEgBinaryWord(egiso.pack(), l1t::TkEm::HWEncoding::CT);
1164       tkems->push_back(tkem);
1165       npho_obj.push_back(tkems->size() - 1);
1166     }
1167     tkemPerBoard->addRegion(npho_obj, board.eta, board.phi);
1168 
1169     nele_obj.clear();
1170     for (const auto &egele : board.egelectron) {
1171       if (egele.hwPt == 0)
1172         continue;
1173 
1174       reco::Candidate::PolarLorentzVector mom(egele.floatPt(), egele.floatVtxEta(), egele.floatVtxPhi(), 0.);
1175 
1176       l1t::TkElectron tkele(reco::Candidate::LorentzVector(mom),
1177                             egele.srcCluster->constituentsAndFractions()[0].first,
1178                             edm::refToPtr(egele.srcTrack->track()),
1179                             egele.floatRelIso(l1ct::EGIsoEleObjEmu::IsoType::TkIso));
1180       tkele.setHwQual(egele.hwQual);
1181       tkele.setPFIsol(egele.floatRelIso(l1ct::EGIsoEleObjEmu::IsoType::PfIso));
1182       tkele.setEgBinaryWord(egele.pack(), l1t::TkElectron::HWEncoding::CT);
1183       tkele.setIdScore(egele.floatIDScore());
1184       tkele.setCharge(egele.intCharge());
1185       tkele.setTrkzVtx(egele.floatZ0());
1186       tkeles->push_back(tkele);
1187       nele_obj.push_back(tkeles->size() - 1);
1188     }
1189     tkelePerBoard->addRegion(nele_obj, board.eta, board.phi);
1190   }
1191 
1192   iEvent.put(std::move(tkems), tkEmLabel);
1193   iEvent.put(std::move(tkemPerBoard), tkEmPerBoardLabel);
1194   iEvent.put(std::move(tkeles), tkEleLabel);
1195   iEvent.put(std::move(tkelePerBoard), tkElePerBoardLabel);
1196 }
1197 
1198 std::unique_ptr<std::vector<unsigned>> L1TCorrelatorLayer1Producer::vecSecInput(InputType t) const {
1199   auto v = std::make_unique<std::vector<unsigned>>();
1200   {
1201     switch (t) {
1202       case caloType:
1203         for (const auto &s : event_.decoded.hadcalo)
1204           v->push_back(s.size());
1205         break;
1206       case emcaloType:
1207         for (const auto &s : event_.decoded.emcalo)
1208           v->push_back(s.size());
1209         break;
1210       case trackType:
1211         for (const auto &s : event_.decoded.track)
1212           v->push_back(s.size());
1213         break;
1214       case l1muType:
1215         v->push_back(event_.decoded.muon.size());
1216         break;
1217     }
1218   }
1219   return v;
1220 }
1221 
1222 std::unique_ptr<std::vector<unsigned>> L1TCorrelatorLayer1Producer::vecRegInput(InputType t) const {
1223   auto v = std::make_unique<std::vector<unsigned>>();
1224   for (const auto &reg : event_.pfinputs) {
1225     switch (t) {
1226       case caloType:
1227         v->push_back(reg.hadcalo.size());
1228         break;
1229       case emcaloType:
1230         v->push_back(reg.emcalo.size());
1231         break;
1232       case trackType:
1233         v->push_back(reg.track.size());
1234         break;
1235       case l1muType:
1236         v->push_back(reg.muon.size());
1237         break;
1238     }
1239   }
1240   return v;
1241 }
1242 
1243 std::unique_ptr<std::vector<unsigned>> L1TCorrelatorLayer1Producer::vecOutput(OutputType i, bool usePuppi) const {
1244   auto v = std::make_unique<std::vector<unsigned>>();
1245   for (const auto &reg : event_.out) {
1246     v->push_back(reg.nObj(i, usePuppi));
1247   }
1248   return v;
1249 }
1250 std::pair<unsigned int, unsigned int> L1TCorrelatorLayer1Producer::totAndMax(
1251     const std::vector<unsigned> &perRegion) const {
1252   unsigned int ntot = 0, nmax = 0;
1253   for (unsigned ni : perRegion) {
1254     ntot += ni;
1255     nmax = std::max(nmax, ni);
1256   }
1257   return std::make_pair(ntot, nmax);
1258 }
1259 //define this as a plug-in
1260 #include "FWCore/Framework/interface/MakerMacros.h"
1261 DEFINE_FWK_MODULE(L1TCorrelatorLayer1Producer);