File indexing completed on 2024-04-06 12:21:28
0001 #include <vector>
0002 #include <numeric>
0003
0004
0005
0006 #include "FWCore/Framework/interface/global/EDProducer.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010
0011 #include "DataFormats/L1TParticleFlow/interface/PFCandidate.h"
0012 #include "DataFormats/L1TParticleFlow/interface/PFJet.h"
0013
0014 #include "DataFormats/L1Trigger/interface/EtSum.h"
0015 #include "DataFormats/Math/interface/LorentzVector.h"
0016
0017
0018 #include "L1Trigger/Phase2L1ParticleFlow/interface/jetmet/L1PFHtEmulator.h"
0019
0020 class L1MhtPfProducer : public edm::global::EDProducer<> {
0021 public:
0022 explicit L1MhtPfProducer(const edm::ParameterSet&);
0023 ~L1MhtPfProducer() override;
0024
0025 private:
0026 void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0027 edm::EDGetTokenT<std::vector<l1t::PFJet>> jetsToken;
0028 float minJetPt;
0029 float maxJetEta;
0030
0031 std::vector<l1ct::Jet> convertEDMToHW(std::vector<l1t::PFJet> edmJets) const;
0032 std::vector<l1t::EtSum> convertHWToEDM(l1ct::Sum hwSums) const;
0033 };
0034
0035 L1MhtPfProducer::L1MhtPfProducer(const edm::ParameterSet& cfg)
0036 : jetsToken(consumes<std::vector<l1t::PFJet>>(cfg.getParameter<edm::InputTag>("jets"))),
0037 minJetPt(cfg.getParameter<double>("minJetPt")),
0038 maxJetEta(cfg.getParameter<double>("maxJetEta")) {
0039 produces<std::vector<l1t::EtSum>>();
0040 }
0041
0042 void L1MhtPfProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0043
0044 l1t::PFJetCollection edmJets = iEvent.get(jetsToken);
0045
0046 std::vector<l1ct::Jet> hwJets = convertEDMToHW(edmJets);
0047
0048 std::vector<l1ct::Jet> hwJetsFiltered;
0049 std::copy_if(hwJets.begin(), hwJets.end(), std::back_inserter(hwJetsFiltered), [&](auto jet) {
0050 return jet.hwPt > l1ct::Scales::makePtFromFloat(minJetPt) &&
0051 std::abs(jet.hwEta) < l1ct::Scales::makeGlbEta(maxJetEta);
0052 });
0053
0054 l1ct::Sum hwSums = htmht(hwJetsFiltered);
0055 std::vector<l1t::EtSum> edmSums = convertHWToEDM(hwSums);
0056
0057
0058 std::unique_ptr<std::vector<l1t::EtSum>> mhtCollection(new std::vector<l1t::EtSum>(0));
0059 mhtCollection->push_back(edmSums.at(0));
0060 mhtCollection->push_back(edmSums.at(1));
0061
0062 iEvent.put(std::move(mhtCollection));
0063 }
0064
0065 std::vector<l1ct::Jet> L1MhtPfProducer::convertEDMToHW(std::vector<l1t::PFJet> edmJets) const {
0066 std::vector<l1ct::Jet> hwJets;
0067 std::for_each(edmJets.begin(), edmJets.end(), [&](l1t::PFJet jet) {
0068 l1ct::Jet hwJet = l1ct::Jet::unpack(jet.getHWJetCT());
0069 hwJets.push_back(hwJet);
0070 });
0071 return hwJets;
0072 }
0073
0074 std::vector<l1t::EtSum> L1MhtPfProducer::convertHWToEDM(l1ct::Sum hwSums) const {
0075 std::vector<l1t::EtSum> edmSums;
0076
0077 reco::Candidate::PolarLorentzVector htVector;
0078 l1gt::Sum gtSum = hwSums.toGT();
0079 htVector.SetPt(l1gt::Scales::floatPt(gtSum.scalar_pt));
0080 htVector.SetPhi(0);
0081 htVector.SetEta(0);
0082
0083 reco::Candidate::PolarLorentzVector mhtVector;
0084 mhtVector.SetPt(l1gt::Scales::floatPt(gtSum.vector_pt));
0085 mhtVector.SetPhi(l1gt::Scales::floatPhi(gtSum.vector_phi));
0086 mhtVector.SetEta(0);
0087
0088 l1t::EtSum ht(htVector, l1t::EtSum::EtSumType::kTotalHt, gtSum.scalar_pt.bits_to_uint64());
0089 l1t::EtSum mht(mhtVector, l1t::EtSum::EtSumType::kMissingHt, gtSum.vector_pt.bits_to_uint64(), 0, gtSum.vector_phi);
0090
0091 edmSums.push_back(ht);
0092 edmSums.push_back(mht);
0093 return edmSums;
0094 }
0095
0096 L1MhtPfProducer::~L1MhtPfProducer() {}
0097
0098 #include "FWCore/Framework/interface/MakerMacros.h"
0099 DEFINE_FWK_MODULE(L1MhtPfProducer);