CMS 3D CMS Logo

TOoLLiPProducer.cc
Go to the documentation of this file.
7 
13 
15 
16 #include <cmath>
17 #include <vector>
18 
19 #include <string>
20 #include "ap_fixed.h"
21 #include "hls4ml/emulator.h"
22 
23 using namespace l1t;
24 
26 public:
27  explicit TOoLLiPProducer(const edm::ParameterSet&);
28  ~TOoLLiPProducer() override = default;
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
31 
32 private:
33  std::unique_ptr<JetId> fJetId_;
34  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
35 
37  bool const fUseRawPt_;
38  double const fMinPt_;
39  double const fMaxEta_;
40  unsigned int const fMaxJets_;
41  int const fNParticles_;
43 
44  hls4mlEmulator::ModelLoader loader;
45  std::shared_ptr<hls4mlEmulator::Model> model;
46 };
47 
49  : jets_(consumes<edm::View<l1t::PFJet>>(cfg.getParameter<edm::InputTag>("jets"))),
50  fUseRawPt_(cfg.getParameter<bool>("useRawPt")),
51  fMinPt_(cfg.getParameter<double>("minPt")),
52  fMaxEta_(cfg.getParameter<double>("maxEta")),
53  fMaxJets_(cfg.getParameter<int>("maxJets")),
54  fNParticles_(cfg.getParameter<int>("nParticles")),
55  fVtxEmu_(consumes<std::vector<l1t::VertexWord>>(cfg.getParameter<edm::InputTag>("vtx"))),
56  loader(hls4mlEmulator::ModelLoader(cfg.getParameter<string>("TOoLLiPVersion"))) {
57  model = loader.load_model();
58  fJetId_ = std::make_unique<JetId>(
59  cfg.getParameter<std::string>("NNInput"), cfg.getParameter<std::string>("NNOutput"), model, fNParticles_);
60  produces<edm::ValueMap<float>>("L1PFLLPJets");
61 }
62 
65  iEvent.getByToken(jets_, jets);
66  float vz = 0.;
67  double ptsum = 0;
69  iEvent.getByToken(fVtxEmu_, vtxEmuHandle);
70  for (const auto& vtx : *vtxEmuHandle) {
71  if (ptsum == 0 || vtx.pt() > ptsum) {
72  ptsum = vtx.pt();
73  vz = vtx.z0();
74  }
75  }
76 
77  std::vector<float> LLPScores;
78  for (const auto& srcjet : *jets) {
79  if (((fUseRawPt_ ? srcjet.rawPt() : srcjet.pt()) < fMinPt_) || std::abs(srcjet.eta()) > fMaxEta_ ||
80  LLPScores.size() >= fMaxJets_) {
81  LLPScores.push_back(-1.);
82  continue;
83  }
84  ap_fixed<16, 6> LLPScore = fJetId_->computeFixed(srcjet, vz, fUseRawPt_);
85  LLPScores.push_back(LLPScore);
86  }
87 
88  auto outT = std::make_unique<edm::ValueMap<float>>();
89  edm::ValueMap<float>::Filler fillerT(*outT);
90  fillerT.insert(jets, LLPScores.begin(), LLPScores.end());
91  fillerT.fill();
92 
93  iEvent.put(std::move(outT), "L1PFLLPJets");
94 }
95 
98  desc.add<edm::InputTag>("jets", edm::InputTag("scPFL1Puppi"));
99  desc.add<bool>("useRawPt", true);
100  desc.add<std::string>("TOoLLiPVersion", std::string("TOoLLiP_v1"));
101  desc.add<std::string>("NNInput", "input:0");
102  desc.add<std::string>("NNOutput", "sequential/dense_2/Sigmoid");
103  desc.add<int>("maxJets", 10);
104  desc.add<int>("nParticles", 10);
105  desc.add<double>("minPt", 20);
106  desc.add<double>("maxEta", 2.4);
107  desc.add<edm::InputTag>("vtx", edm::InputTag("L1VertexFinderEmulator", "L1VerticesEmulation"));
108  descriptions.add("TOoLLiPProducer", desc);
109 }
110 
112 DEFINE_FWK_MODULE(TOoLLiPProducer);
TOoLLiPProducer(const edm::ParameterSet &)
hls4mlEmulator::ModelLoader loader
double const fMinPt_
std::shared_ptr< hls4mlEmulator::Model > model
int const fNParticles_
delete x;
Definition: CaloConfig.h:22
unsigned int const fMaxJets_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< l1t::PFJet > > const jets_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double const fMaxEta_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
std::unique_ptr< JetId > fJetId_
bool const fUseRawPt_
edm::EDGetTokenT< std::vector< l1t::VertexWord > > const fVtxEmu_
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override