CMS 3D CMS Logo

L1MHtPFProducer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <numeric>
3 
10 
13 
16 
17 // bitwise emulation headers
19 
21 public:
22  explicit L1MhtPfProducer(const edm::ParameterSet&);
23  ~L1MhtPfProducer() override;
24 
25 private:
26  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
28  float minJetPt;
29  float maxJetEta;
30 
31  std::vector<l1ct::Jet> convertEDMToHW(std::vector<l1t::PFJet> edmJets) const;
32  std::vector<l1t::EtSum> convertHWToEDM(l1ct::Sum hwSums) const;
33 };
34 
36  : jetsToken(consumes<std::vector<l1t::PFJet>>(cfg.getParameter<edm::InputTag>("jets"))),
37  minJetPt(cfg.getParameter<double>("minJetPt")),
38  maxJetEta(cfg.getParameter<double>("maxJetEta")) {
39  produces<std::vector<l1t::EtSum>>();
40 }
41 
43  // Get the jets from the event
44  l1t::PFJetCollection edmJets = iEvent.get(jetsToken);
45 
46  // Apply pT and eta selections
47  l1t::PFJetCollection edmJetsFiltered;
48  std::copy_if(edmJets.begin(), edmJets.end(), std::back_inserter(edmJetsFiltered), [&](auto jet) {
49  return jet.pt() > minJetPt && std::abs(jet.eta()) < maxJetEta;
50  });
51 
52  // Run the emulation
53  std::vector<l1ct::Jet> hwJets = convertEDMToHW(edmJetsFiltered); // convert to the emulator format
54  l1ct::Sum hwSums = htmht(hwJets); // call the emulator
55  std::vector<l1t::EtSum> edmSums = convertHWToEDM(hwSums); // convert back to edm format
56 
57  // Put the sums in the event
58  std::unique_ptr<std::vector<l1t::EtSum>> mhtCollection(new std::vector<l1t::EtSum>(0));
59  mhtCollection->push_back(edmSums.at(0)); // HT
60  mhtCollection->push_back(edmSums.at(1)); // MHT
61 
62  iEvent.put(std::move(mhtCollection));
63 }
64 
65 std::vector<l1ct::Jet> L1MhtPfProducer::convertEDMToHW(std::vector<l1t::PFJet> edmJets) const {
66  std::vector<l1ct::Jet> hwJets;
67  std::for_each(edmJets.begin(), edmJets.end(), [&](l1t::PFJet jet) {
68  l1ct::Jet hwJet = l1ct::Jet::unpack(jet.encodedJet());
69  hwJets.push_back(hwJet);
70  });
71  return hwJets;
72 }
73 
74 std::vector<l1t::EtSum> L1MhtPfProducer::convertHWToEDM(l1ct::Sum hwSums) const {
75  std::vector<l1t::EtSum> edmSums;
76 
78  htVector.SetPt(hwSums.hwSumPt.to_double());
79  htVector.SetPhi(0);
80  htVector.SetEta(0);
81 
83  mhtVector.SetPt(hwSums.hwPt.to_double());
84  mhtVector.SetPhi(l1ct::Scales::floatPhi(hwSums.hwPhi));
85  mhtVector.SetEta(0);
86 
89 
90  edmSums.push_back(ht);
91  edmSums.push_back(mht);
92  return edmSums;
93 }
94 
96 
Definition: jets.h:12
pt_t hwPt
Definition: sums.h:11
l1ct::Sum htmht(std::vector< l1ct::Jet > jets)
delete x;
Definition: CaloConfig.h:22
glbphi_t hwPhi
Definition: sums.h:12
std::vector< l1t::EtSum > convertHWToEDM(l1ct::Sum hwSums) const
std::vector< l1t::PFJet > PFJetCollection
Definition: PFJet.h:49
L1MhtPfProducer(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
~L1MhtPfProducer() override
edm::EDGetTokenT< std::vector< l1t::PFJet > > jetsToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static Jet unpack(const std::array< uint64_t, 2 > &src)
Definition: jets.h:68
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
HLT enums.
Definition: sums.h:10
float floatPhi(phi_t phi)
Definition: datatypes.h:158
pt_t hwSumPt
Definition: sums.h:13
std::vector< l1ct::Jet > convertEDMToHW(std::vector< l1t::PFJet > edmJets) const
def move(src, dest)
Definition: eostools.py:511
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38