Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-12 04:16:40

0001 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0002 #include "DataFormats/Common/interface/Handle.h"
0003 #include "FWCore/Utilities/interface/EDGetToken.h"
0004 
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "DataFormats/Common/interface/View.h"
0007 
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/global/EDProducer.h"
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "DataFormats/L1Trigger/interface/P2GTCandidate.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "L1Trigger/Phase2L1GT/interface/L1GTScales.h"
0016 
0017 #include "DataFormats/L1Trigger/interface/TkJetWord.h"
0018 #include "DataFormats/L1Trigger/interface/VertexWord.h"
0019 
0020 #include "DataFormats/L1TMuonPhase2/interface/SAMuon.h"
0021 #include "DataFormats/L1TMuonPhase2/interface/TrackerMuon.h"
0022 
0023 #include "DataFormats/L1TParticleFlow/interface/PFJet.h"
0024 #include "DataFormats/L1TCorrelator/interface/TkEmFwd.h"
0025 #include "DataFormats/L1TCorrelator/interface/TkEm.h"
0026 #include "DataFormats/L1TCorrelator/interface/TkElectronFwd.h"
0027 #include "DataFormats/L1TCorrelator/interface/TkElectron.h"
0028 #include "DataFormats/L1TParticleFlow/interface/PFTau.h"
0029 #include "DataFormats/L1TParticleFlow/interface/gt_datatypes.h"
0030 
0031 #include "DataFormats/L1Trigger/interface/EtSum.h"
0032 
0033 #include "L1Trigger/L1TTrackMatch/interface/L1TkHTMissEmulatorProducer.h"
0034 
0035 #include <vector>
0036 #include <array>
0037 #include <string>
0038 #include <type_traits>
0039 
0040 namespace l1t {
0041 
0042   class L1GTProducer : public edm::global::EDProducer<> {
0043   public:
0044     explicit L1GTProducer(const edm::ParameterSet &);
0045     ~L1GTProducer() override = default;
0046 
0047     static void fillDescriptions(edm::ConfigurationDescriptions &);
0048 
0049   private:
0050     void produceGTTPromptJets(edm::Event &event) const;
0051     void produceGTTDisplacedJets(edm::Event &event) const;
0052     void produceGTTPromptHtSum(edm::Event &event) const;
0053     void produceGTTDisplacedHtSum(edm::Event &event) const;
0054     void produceGTTEtSum(edm::Event &event) const;
0055     void produceGTTPrimaryVert(edm::Event &event) const;
0056 
0057     void produceGMTSaPromptMuons(edm::Event &event) const;
0058     void produceGMTSaDisplacedMuons(edm::Event &event) const;
0059     void produceGMTTkMuons(edm::Event &event) const;
0060 
0061     void produceCL2JetsSC4(edm::Event &event) const;
0062     void produceCL2JetsSC8(edm::Event &event) const;
0063     void produceCL2Photons(edm::Event &event) const;
0064     void produceCL2Electrons(edm::Event &event) const;
0065     void produceCL2Taus(edm::Event &event) const;
0066     void produceCL2EtSum(edm::Event &event) const;
0067     void produceCl2HtSum(edm::Event &event) const;
0068 
0069     void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
0070 
0071     const L1GTScales scales_;
0072 
0073     const edm::EDGetTokenT<TkJetWordCollection> gttPromptJetToken_;
0074     const edm::EDGetTokenT<TkJetWordCollection> gttDisplacedJetToken_;
0075     const edm::EDGetTokenT<std::vector<l1t::EtSum>> gttPromptHtSumToken_;
0076     const edm::EDGetTokenT<std::vector<l1t::EtSum>> gttDisplacedHtSumToken_;
0077     const edm::EDGetTokenT<std::vector<l1t::EtSum>> gttEtSumToken_;
0078     const edm::EDGetTokenT<VertexWordCollection> gttPrimaryVertexToken_;
0079 
0080     const edm::EDGetTokenT<SAMuonCollection> gmtSaPromptMuonToken_;
0081     const edm::EDGetTokenT<SAMuonCollection> gmtSaDisplacedMuonToken_;
0082     const edm::EDGetTokenT<TrackerMuonCollection> gmtTkMuonToken_;
0083 
0084     const edm::EDGetTokenT<PFJetCollection> cl2JetSC4Token_;
0085     const edm::EDGetTokenT<PFJetCollection> cl2JetSC8Token_;
0086     const edm::EDGetTokenT<TkEmCollection> cl2PhotonToken_;
0087     const edm::EDGetTokenT<TkElectronCollection> cl2ElectronToken_;
0088     const edm::EDGetTokenT<PFTauCollection> cl2TauToken_;
0089     const edm::EDGetTokenT<std::vector<l1t::EtSum>> cl2EtSumToken_;
0090     const edm::EDGetTokenT<std::vector<l1t::EtSum>> cl2HtSumToken_;
0091   };
0092 
0093   L1GTProducer::L1GTProducer(const edm::ParameterSet &config)
0094       : scales_(config.getParameter<edm::ParameterSet>("scales")),
0095         gttPromptJetToken_(consumes<TkJetWordCollection>(config.getParameter<edm::InputTag>("GTTPromptJets"))),
0096         gttDisplacedJetToken_(consumes<TkJetWordCollection>(config.getParameter<edm::InputTag>("GTTDisplacedJets"))),
0097         gttPromptHtSumToken_(consumes<std::vector<l1t::EtSum>>(config.getParameter<edm::InputTag>("GTTPromptHtSum"))),
0098         gttDisplacedHtSumToken_(
0099             consumes<std::vector<l1t::EtSum>>(config.getParameter<edm::InputTag>("GTTDisplacedHtSum"))),
0100         gttEtSumToken_(consumes<std::vector<l1t::EtSum>>(config.getParameter<edm::InputTag>("GTTEtSum"))),
0101         gttPrimaryVertexToken_(consumes<VertexWordCollection>(config.getParameter<edm::InputTag>("GTTPrimaryVert"))),
0102         gmtSaPromptMuonToken_(consumes<SAMuonCollection>(config.getParameter<edm::InputTag>("GMTSaPromptMuons"))),
0103         gmtSaDisplacedMuonToken_(consumes<SAMuonCollection>(config.getParameter<edm::InputTag>("GMTSaDisplacedMuons"))),
0104         gmtTkMuonToken_(consumes<TrackerMuonCollection>(config.getParameter<edm::InputTag>("GMTTkMuons"))),
0105         cl2JetSC4Token_(consumes<PFJetCollection>(config.getParameter<edm::InputTag>("CL2JetsSC4"))),
0106         cl2JetSC8Token_(consumes<PFJetCollection>(config.getParameter<edm::InputTag>("CL2JetsSC8"))),
0107         cl2PhotonToken_(consumes<TkEmCollection>(config.getParameter<edm::InputTag>("CL2Photons"))),
0108         cl2ElectronToken_(consumes<TkElectronCollection>(config.getParameter<edm::InputTag>("CL2Electrons"))),
0109         cl2TauToken_(consumes<PFTauCollection>(config.getParameter<edm::InputTag>("CL2Taus"))),
0110         cl2EtSumToken_(consumes<std::vector<l1t::EtSum>>(config.getParameter<edm::InputTag>("CL2EtSum"))),
0111         cl2HtSumToken_(consumes<std::vector<l1t::EtSum>>(config.getParameter<edm::InputTag>("CL2HtSum"))) {
0112     produces<P2GTCandidateCollection>("GTTPromptJets");
0113     produces<P2GTCandidateCollection>("GTTDisplacedJets");
0114     produces<P2GTCandidateCollection>("GTTPromptHtSum");
0115     produces<P2GTCandidateCollection>("GTTDisplacedHtSum");
0116     produces<P2GTCandidateCollection>("GTTEtSum");
0117     produces<P2GTCandidateCollection>("GTTPrimaryVert");
0118 
0119     produces<P2GTCandidateCollection>("GMTSaPromptMuons");
0120     produces<P2GTCandidateCollection>("GMTSaDisplacedMuons");
0121     produces<P2GTCandidateCollection>("GMTTkMuons");
0122 
0123     produces<P2GTCandidateCollection>("CL2JetsSC4");
0124     produces<P2GTCandidateCollection>("CL2JetsSC8");
0125     produces<P2GTCandidateCollection>("CL2Photons");
0126     produces<P2GTCandidateCollection>("CL2Electrons");
0127     produces<P2GTCandidateCollection>("CL2Taus");
0128     produces<P2GTCandidateCollection>("CL2EtSum");
0129     produces<P2GTCandidateCollection>("CL2HtSum");
0130   }
0131 
0132   void L1GTProducer::fillDescriptions(edm::ConfigurationDescriptions &description) {
0133     edm::ParameterSetDescription desc;
0134 
0135     edm::ParameterSetDescription scalesDesc;
0136     L1GTScales::fillPSetDescription(scalesDesc);
0137     desc.add<edm::ParameterSetDescription>("scales", scalesDesc);
0138 
0139     desc.add<edm::InputTag>("GTTPromptJets");
0140     desc.add<edm::InputTag>("GTTDisplacedJets");
0141     desc.add<edm::InputTag>("GTTPromptHtSum");
0142     desc.add<edm::InputTag>("GTTDisplacedHtSum");
0143     desc.add<edm::InputTag>("GTTEtSum");
0144     desc.add<edm::InputTag>("GTTPrimaryVert");
0145 
0146     desc.add<edm::InputTag>("GMTSaPromptMuons");
0147     desc.add<edm::InputTag>("GMTSaDisplacedMuons");
0148     desc.add<edm::InputTag>("GMTTkMuons");
0149 
0150     desc.add<edm::InputTag>("CL2JetsSC4");
0151     desc.add<edm::InputTag>("CL2JetsSC8");
0152     desc.add<edm::InputTag>("CL2Photons");
0153     desc.add<edm::InputTag>("CL2Electrons");
0154     desc.add<edm::InputTag>("CL2Taus");
0155     desc.add<edm::InputTag>("CL2EtSum");
0156     desc.add<edm::InputTag>("CL2HtSum");
0157 
0158     description.addWithDefaultLabel(desc);
0159   }
0160 
0161   void L1GTProducer::produceGTTPrimaryVert(edm::Event &event) const {
0162     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0163     const VertexWordCollection &collection = event.get(gttPrimaryVertexToken_);
0164     for (std::size_t i = 0; i < collection.size() && i < 10; i++) {
0165       const VertexWord &obj = collection[i];
0166       int hwZ0 = obj.z0Word().V.to_int() * 5;
0167       P2GTCandidate gtObj(
0168           0, reco::ParticleState::PolarLorentzVector(), reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0169       gtObj.hwZ0_ = hwZ0;
0170       gtObj.hwQualityScore_ = obj.qualityWord().V.to_int();
0171       gtObj.hwSum_pT_pv_ = obj.multiplicityWord().V.to_int();
0172       gtObj.hwNumber_of_tracks_in_pv_ = obj.multiplicityWord().V.to_int();
0173       gtObj.hwNumber_of_tracks_not_in_pv_ = obj.inverseMultiplicityWord().V.to_int();
0174       gtObj.objectType_ = P2GTCandidate::GTTPrimaryVert;
0175 
0176       outputCollection->push_back(gtObj);
0177     }
0178     event.put(std::move(outputCollection), "GTTPrimaryVert");
0179   }
0180 
0181   void L1GTProducer::produceGTTPromptJets(edm::Event &event) const {
0182     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0183     const TkJetWordCollection &collection = event.get(gttPromptJetToken_);
0184     for (std::size_t i = 0; i < collection.size() && i < 12; i++) {
0185       const TkJetWord &obj = collection[i];
0186       int hwZ0 = obj.z0Word().V.to_int() << 7;
0187       P2GTCandidate gtObj(0,
0188                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.ptWord().V.to_int()),
0189                                                                   scales_.to_eta(obj.glbEtaWord().V.to_int()),
0190                                                                   scales_.to_phi(obj.glbPhiWord().V.to_int()),
0191                                                                   0),
0192                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0193       gtObj.hwPT_ = obj.ptWord().V.to_int();
0194       gtObj.hwPhi_ = obj.glbPhiWord().V.to_int();
0195       gtObj.hwEta_ = obj.glbEtaWord().V.to_int();
0196       gtObj.hwZ0_ = hwZ0;
0197       gtObj.hwNumber_of_tracks_ = obj.ntWord().V.to_int();
0198       gtObj.objectType_ = P2GTCandidate::GTTPromptJets;
0199 
0200       outputCollection->push_back(gtObj);
0201     }
0202     event.put(std::move(outputCollection), "GTTPromptJets");
0203   }
0204 
0205   void L1GTProducer::produceGTTDisplacedJets(edm::Event &event) const {
0206     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0207     const TkJetWordCollection &collection = event.get(gttDisplacedJetToken_);
0208     for (std::size_t i = 0; i < collection.size() && i < 12; i++) {
0209       const TkJetWord &obj = collection[i];
0210       int hwZ0 = obj.z0Word().V.to_int() << 7;
0211       P2GTCandidate gtObj(0,
0212                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.ptWord().V.to_int()),
0213                                                                   scales_.to_eta(obj.glbEtaWord().V.to_int()),
0214                                                                   scales_.to_phi(obj.glbPhiWord().V.to_int()),
0215                                                                   0),
0216                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0217       gtObj.hwPT_ = obj.ptWord().V.to_int();
0218       gtObj.hwPhi_ = obj.glbPhiWord().V.to_int();
0219       gtObj.hwEta_ = obj.glbEtaWord().V.to_int();
0220       gtObj.hwZ0_ = hwZ0;
0221       gtObj.hwNumber_of_tracks_ = obj.ntWord().V.to_int();
0222       gtObj.objectType_ = P2GTCandidate::GTTDisplacedJets;
0223 
0224       outputCollection->push_back(gtObj);
0225     }
0226     event.put(std::move(outputCollection), "GTTDisplacedJets");
0227   }
0228 
0229   void L1GTProducer::produceGTTPromptHtSum(edm::Event &event) const {
0230     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0231     const std::vector<l1t::EtSum> &collection = event.get(gttPromptHtSumToken_);
0232     if (collection.size() > 0) {
0233       const l1t::EtSum &obj = collection[0];
0234       l1tmhtemu::EtMiss htMiss;
0235       htMiss.Et = obj.p4().energy();
0236       P2GTCandidate gtObj(0,
0237                           reco::ParticleState::PolarLorentzVector(
0238                               scales_.to_pT(htMiss.Et.V.to_int()), 0, scales_.to_phi(obj.hwPhi()), 0));
0239 
0240       gtObj.hwPT_ = htMiss.Et.V.to_int();
0241       gtObj.hwPhi_ = obj.hwPhi();
0242       gtObj.hwScalarSumPT_ = obj.hwPt();
0243       gtObj.objectType_ = P2GTCandidate::GTTPromptHtSum;
0244 
0245       outputCollection->push_back(gtObj);
0246     }
0247 
0248     event.put(std::move(outputCollection), "GTTPromptHtSum");
0249   }
0250 
0251   void L1GTProducer::produceGTTDisplacedHtSum(edm::Event &event) const {
0252     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0253     const std::vector<l1t::EtSum> &collection = event.get(gttDisplacedHtSumToken_);
0254     if (collection.size() > 0) {
0255       const l1t::EtSum &obj = collection[0];
0256       l1tmhtemu::EtMiss htMiss;
0257       htMiss.Et = obj.p4().energy();
0258       P2GTCandidate gtObj(0,
0259                           reco::ParticleState::PolarLorentzVector(
0260                               scales_.to_pT(htMiss.Et.V.to_int()), 0, scales_.to_phi(obj.hwPhi()), 0));
0261 
0262       gtObj.hwPT_ = htMiss.Et.V.to_int();
0263       gtObj.hwPhi_ = obj.hwPhi();
0264       gtObj.hwScalarSumPT_ = obj.hwPt();
0265       gtObj.objectType_ = P2GTCandidate::GTTDisplacedHtSum;
0266 
0267       outputCollection->push_back(gtObj);
0268     }
0269 
0270     event.put(std::move(outputCollection), "GTTDisplacedHtSum");
0271   }
0272 
0273   void L1GTProducer::produceGTTEtSum(edm::Event &event) const {
0274     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0275     const std::vector<l1t::EtSum> &collection = event.get(gttEtSumToken_);
0276     if (collection.size() > 0) {
0277       const l1t::EtSum &obj = collection[0];
0278       P2GTCandidate gtObj(
0279           0, reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.hwPt()), 0, scales_.to_phi(obj.hwPhi()), 0));
0280       gtObj.hwPT_ = obj.hwPt();
0281       gtObj.hwPhi_ = obj.hwPhi();
0282       gtObj.objectType_ = P2GTCandidate::GTTEtSum;
0283 
0284       outputCollection->push_back(gtObj);
0285     }
0286 
0287     event.put(std::move(outputCollection), "GTTEtSum");
0288   }
0289 
0290   void L1GTProducer::produceGMTSaPromptMuons(edm::Event &event) const {
0291     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0292     const SAMuonCollection &collection = event.get(gmtSaPromptMuonToken_);
0293     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0294       const SAMuon &obj = collection[i];
0295       int hwZ0 = obj.apZ0().to_int() << 12;
0296       P2GTCandidate gtObj(scales_.to_chg(obj.apCharge().to_int()),
0297                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.apPt().to_int()),
0298                                                                   scales_.to_eta(obj.apEta().to_int()),
0299                                                                   scales_.to_phi(obj.apPhi().to_int()),
0300                                                                   0),
0301                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0302       gtObj.hwPT_ = obj.apPt().to_int();
0303       gtObj.hwPhi_ = obj.apPhi().to_int();
0304       gtObj.hwEta_ = obj.apEta().to_int();
0305       gtObj.hwZ0_ = hwZ0;
0306       gtObj.hwQualityFlags_ = obj.apQualFlags().to_int();
0307       gtObj.hwCharge_ = obj.apCharge().to_int();
0308       gtObj.hwD0_ = obj.apD0().to_int();
0309       gtObj.objectType_ = P2GTCandidate::GMTSaPromptMuons;
0310 
0311       outputCollection->push_back(gtObj);
0312     }
0313     event.put(std::move(outputCollection), "GMTSaPromptMuons");
0314   }
0315 
0316   void L1GTProducer::produceGMTSaDisplacedMuons(edm::Event &event) const {
0317     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0318     const SAMuonCollection &collection = event.get(gmtSaDisplacedMuonToken_);
0319     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0320       const SAMuon &obj = collection[i];
0321       int hwZ0 = obj.apZ0().to_int() << 12;
0322       P2GTCandidate gtObj(scales_.to_chg(obj.apCharge().to_int()),
0323                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.apPt().to_int()),
0324                                                                   scales_.to_eta(obj.apEta().to_int()),
0325                                                                   scales_.to_phi(obj.apPhi().to_int()),
0326                                                                   0),
0327                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0328       gtObj.hwPT_ = obj.apPt().to_int();
0329       gtObj.hwPhi_ = obj.apPhi().to_int();
0330       gtObj.hwEta_ = obj.apEta().to_int();
0331       gtObj.hwZ0_ = hwZ0;
0332       gtObj.hwQualityFlags_ = obj.apQualFlags().to_int();
0333       gtObj.hwCharge_ = obj.apCharge().to_int();
0334       gtObj.hwD0_ = obj.apD0().to_int();
0335       gtObj.objectType_ = P2GTCandidate::GMTSaDisplacedMuons;
0336 
0337       outputCollection->push_back(gtObj);
0338     }
0339     event.put(std::move(outputCollection), "GMTSaDisplacedMuons");
0340   }
0341 
0342   void L1GTProducer::produceGMTTkMuons(edm::Event &event) const {
0343     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0344     const TrackerMuonCollection &collection = event.get(gmtTkMuonToken_);
0345     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0346       const TrackerMuon &obj = collection[i];
0347       int hwZ0 = obj.apZ0().to_int() << 7;
0348       P2GTCandidate gtObj(scales_.to_chg(obj.apCharge().to_int()),
0349                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(obj.apPt().to_int()),
0350                                                                   scales_.to_eta(obj.apEta().to_int()),
0351                                                                   scales_.to_phi(obj.apPhi().to_int()),
0352                                                                   0),
0353                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0354       gtObj.hwPT_ = obj.apPt().to_int();
0355       gtObj.hwPhi_ = obj.apPhi().to_int();
0356       gtObj.hwEta_ = obj.apEta().to_int();
0357       gtObj.hwZ0_ = hwZ0;
0358       gtObj.hwQualityFlags_ = obj.apQualFlags().to_int();
0359       gtObj.hwIsolationPT_ = obj.apIso().to_int();
0360       gtObj.hwCharge_ = obj.apCharge().to_int();
0361       gtObj.hwD0_ = obj.apD0().to_int();
0362       gtObj.hwBeta_ = obj.apBeta().to_int();
0363       gtObj.objectType_ = P2GTCandidate::GMTTkMuons;
0364 
0365       outputCollection->push_back(gtObj);
0366     }
0367     event.put(std::move(outputCollection), "GMTTkMuons");
0368   }
0369 
0370   void L1GTProducer::produceCL2JetsSC4(edm::Event &event) const {
0371     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0372     const PFJetCollection &collection = event.get(cl2JetSC4Token_);
0373     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0374       l1gt::Jet gtJet = l1gt::Jet::unpack(collection[i].getHWJetGT());
0375       int hwZ0 = gtJet.z0.V.to_int() << 7;
0376       P2GTCandidate gtObj(0,
0377                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(gtJet.v3.pt.V.to_int()),
0378                                                                   scales_.to_eta(gtJet.v3.eta.V.to_int()),
0379                                                                   scales_.to_phi(gtJet.v3.phi.V.to_int()),
0380                                                                   0),
0381                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0382       gtObj.hwPT_ = gtJet.v3.pt.V.to_int();
0383       gtObj.hwPhi_ = gtJet.v3.phi.V.to_int();
0384       gtObj.hwEta_ = gtJet.v3.eta.V.to_int();
0385       gtObj.hwZ0_ = hwZ0;
0386       gtObj.objectType_ = P2GTCandidate::CL2JetsSC4;
0387 
0388       outputCollection->push_back(gtObj);
0389     }
0390     event.put(std::move(outputCollection), "CL2JetsSC4");
0391   }
0392 
0393   void L1GTProducer::produceCL2JetsSC8(edm::Event &event) const {
0394     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0395     const PFJetCollection &collection = event.get(cl2JetSC8Token_);
0396     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0397       l1gt::Jet gtJet = l1gt::Jet::unpack(collection[i].getHWJetGT());
0398       int hwZ0 = gtJet.z0.V.to_int() << 7;
0399       P2GTCandidate gtObj(0,
0400                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(gtJet.v3.pt.V.to_int()),
0401                                                                   scales_.to_eta(gtJet.v3.eta.V.to_int()),
0402                                                                   scales_.to_phi(gtJet.v3.phi.V.to_int()),
0403                                                                   0),
0404                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0405       gtObj.hwPT_ = gtJet.v3.pt.V.to_int();
0406       gtObj.hwPhi_ = gtJet.v3.phi.V.to_int();
0407       gtObj.hwEta_ = gtJet.v3.eta.V.to_int();
0408       gtObj.hwZ0_ = hwZ0;
0409       gtObj.objectType_ = P2GTCandidate::CL2JetsSC8;
0410 
0411       outputCollection->push_back(gtObj);
0412     }
0413     event.put(std::move(outputCollection), "CL2JetsSC8");
0414   }
0415 
0416   void L1GTProducer::produceCL2Photons(edm::Event &event) const {
0417     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0418     const TkEmCollection &collection = event.get(cl2PhotonToken_);
0419     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0420       l1gt::Photon gtPhoton = collection[i].hwObj();
0421       P2GTCandidate gtObj(0,
0422                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(gtPhoton.v3.pt.V.to_int()),
0423                                                                   scales_.to_eta(gtPhoton.v3.eta.V.to_int()),
0424                                                                   scales_.to_phi(gtPhoton.v3.phi.V.to_int()),
0425                                                                   0));
0426       gtObj.hwPT_ = gtPhoton.v3.pt.V.to_int();
0427       gtObj.hwPhi_ = gtPhoton.v3.phi.V.to_int();
0428       gtObj.hwEta_ = gtPhoton.v3.eta.V.to_int();
0429       gtObj.hwIsolationPT_ = gtPhoton.isolationPT.V.to_int();
0430       gtObj.hwQualityFlags_ = gtPhoton.qualityFlags.V.to_int();
0431       gtObj.objectType_ = P2GTCandidate::CL2Photons;
0432 
0433       outputCollection->push_back(gtObj);
0434     }
0435     event.put(std::move(outputCollection), "CL2Photons");
0436   }
0437 
0438   void L1GTProducer::produceCL2Electrons(edm::Event &event) const {
0439     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0440     const TkElectronCollection &collection = event.get(cl2ElectronToken_);
0441     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0442       l1gt::Electron gtElectron = collection[i].hwObj();
0443       int hwZ0 = gtElectron.z0.V.to_int() << 7;
0444       P2GTCandidate gtObj(scales_.to_chg(gtElectron.charge.V.to_int()),
0445                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(gtElectron.v3.pt.V.to_int()),
0446                                                                   scales_.to_eta(gtElectron.v3.eta.V.to_int()),
0447                                                                   scales_.to_phi(gtElectron.v3.phi.V.to_int()),
0448                                                                   0),
0449                           reco::ParticleState::Point(0, 0, scales_.to_z0(hwZ0)));
0450       gtObj.hwPT_ = gtElectron.v3.pt.V.to_int();
0451       gtObj.hwPhi_ = gtElectron.v3.phi.V.to_int();
0452       gtObj.hwEta_ = gtElectron.v3.eta.V.to_int();
0453       gtObj.hwZ0_ = hwZ0;
0454       gtObj.hwIsolationPT_ = gtElectron.isolationPT.V.to_int();
0455       gtObj.hwQualityFlags_ = gtElectron.qualityFlags.V.to_int();
0456       gtObj.hwCharge_ = gtElectron.charge.V.to_int();
0457       gtObj.objectType_ = P2GTCandidate::CL2Electrons;
0458 
0459       outputCollection->push_back(gtObj);
0460     }
0461     event.put(std::move(outputCollection), "CL2Electrons");
0462   }
0463 
0464   void L1GTProducer::produceCL2Taus(edm::Event &event) const {
0465     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0466     const PFTauCollection &collection = event.get(cl2TauToken_);
0467     for (size_t i = 0; i < collection.size() && i < 12; i++) {
0468       l1gt::Tau gtTau = collection[i].getHWTauGT();
0469       P2GTCandidate gtObj(scales_.to_chg(gtTau.charge.V.to_int()),
0470                           reco::ParticleState::PolarLorentzVector(scales_.to_pT(gtTau.v3.pt.V.to_int()),
0471                                                                   scales_.to_eta(gtTau.v3.eta.V.to_int()),
0472                                                                   scales_.to_phi(gtTau.v3.phi.V.to_int()),
0473                                                                   0));
0474       gtObj.hwPT_ = gtTau.v3.pt.V.to_int();
0475       gtObj.hwPhi_ = gtTau.v3.phi.V.to_int();
0476       gtObj.hwEta_ = gtTau.v3.eta.V.to_int();
0477       gtObj.hwSeed_pT_ = gtTau.seed_pt.V.to_int();
0478       gtObj.hwSeed_z0_ = gtTau.seed_z0.V.to_int();
0479       gtObj.hwCharge_ = gtTau.charge.V.to_int();
0480       gtObj.hwType_ = gtTau.type.V.to_int();
0481       gtObj.hwQualityScore_ = gtTau.quality.V.to_int();
0482       gtObj.objectType_ = P2GTCandidate::CL2Taus;
0483 
0484       outputCollection->push_back(gtObj);
0485     }
0486     event.put(std::move(outputCollection), "CL2Taus");
0487   }
0488 
0489   void L1GTProducer::produceCL2EtSum(edm::Event &event) const {
0490     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0491     const std::vector<EtSum> &collection = event.get(cl2EtSumToken_);
0492     const EtSum &met = collection[0];
0493 
0494     l1gt::Sum sum{true /* valid */, met.pt(), met.phi() / l1gt::Scales::ETAPHI_LSB, 0 /* scalar sum */};
0495 
0496     P2GTCandidate gtObj(0,
0497                         reco::ParticleState::PolarLorentzVector(
0498                             scales_.to_pT(sum.vector_pt.V.to_int()), 0, scales_.to_phi(sum.vector_phi.V.to_int()), 0));
0499     gtObj.hwPT_ = sum.vector_pt.V.to_int();
0500     gtObj.hwPhi_ = sum.vector_phi.V.to_int();
0501     gtObj.hwScalarSumPT_ = sum.scalar_pt.V.to_int();
0502     gtObj.objectType_ = P2GTCandidate::CL2EtSum;
0503 
0504     outputCollection->push_back(gtObj);
0505     event.put(std::move(outputCollection), "CL2EtSum");
0506   }
0507 
0508   void L1GTProducer::produceCl2HtSum(edm::Event &event) const {
0509     std::unique_ptr<P2GTCandidateCollection> outputCollection = std::make_unique<P2GTCandidateCollection>();
0510     const std::vector<EtSum> &collection = event.get(cl2HtSumToken_);
0511     const EtSum &ht = collection[0];
0512     const EtSum &mht = collection[1];
0513 
0514     P2GTCandidate gtObj(
0515         0, reco::ParticleState::PolarLorentzVector(scales_.to_pT(mht.hwPt()), 0, scales_.to_phi(mht.hwPhi()), 0));
0516     gtObj.hwPT_ = mht.hwPt();
0517     gtObj.hwPhi_ = mht.hwPhi();
0518     gtObj.hwScalarSumPT_ = ht.hwPt();
0519     gtObj.objectType_ = P2GTCandidate::CL2HtSum;
0520 
0521     outputCollection->push_back(gtObj);
0522     event.put(std::move(outputCollection), "CL2HtSum");
0523   }
0524 
0525   void L1GTProducer::produce(edm::StreamID, edm::Event &event, const edm::EventSetup &setup) const {
0526     produceGTTPromptJets(event);
0527     produceGTTDisplacedJets(event);
0528     produceGTTPromptHtSum(event);
0529     produceGTTDisplacedHtSum(event);
0530     produceGTTEtSum(event);
0531     produceGTTPrimaryVert(event);
0532 
0533     produceGMTSaPromptMuons(event);
0534     produceGMTSaDisplacedMuons(event);
0535     produceGMTTkMuons(event);
0536 
0537     produceCL2JetsSC4(event);
0538     produceCL2JetsSC8(event);
0539     produceCL2Photons(event);
0540     produceCL2Electrons(event);
0541     produceCL2Taus(event);
0542     produceCL2EtSum(event);
0543     produceCl2HtSum(event);
0544   }
0545 }  // namespace l1t
0546 
0547 using namespace l1t;
0548 
0549 DEFINE_FWK_MODULE(L1GTProducer);