CMS 3D CMS Logo

DeepMETSonicProducer.cc
Go to the documentation of this file.
10 
11 using namespace deepmet_helper;
12 
14 public:
15  explicit DeepMETSonicProducer(const edm::ParameterSet&);
16  void acquire(edm::Event const& iEvent, edm::EventSetup const& iSetup, Input& iInput) override;
17  void produce(edm::Event& iEvent, edm::EventSetup const& iSetup, Output const& iOutput) override;
18  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
19 
20 private:
22  const float norm_;
23  const bool ignore_leptons_;
24  const unsigned int max_n_pf_;
25  const float scale_;
26  float px_leptons_;
27  float py_leptons_;
28 };
29 
31  : TritonEDProducer<>(cfg),
32  pf_token_(consumes<std::vector<pat::PackedCandidate>>(cfg.getParameter<edm::InputTag>("pf_src"))),
33  norm_(cfg.getParameter<double>("norm_factor")),
34  ignore_leptons_(cfg.getParameter<bool>("ignore_leptons")),
35  max_n_pf_(cfg.getParameter<unsigned int>("max_n_pf")),
36  scale_(1.0 / norm_) {
37  produces<pat::METCollection>();
38 }
39 
41  // one event per batch
42  client_->setBatchSize(1);
43  px_leptons_ = 0.;
44  py_leptons_ = 0.;
45 
46  auto const& pfs = iEvent.get(pf_token_);
47 
48  auto& input = iInput.at("input");
49  auto pfdata = input.allocate<float>();
50  auto& vpfdata = (*pfdata)[0];
51 
52  auto& input_cat0 = iInput.at("input_cat0");
53  auto pfchg = input_cat0.allocate<float>();
54  auto& vpfchg = (*pfchg)[0];
55 
56  auto& input_cat1 = iInput.at("input_cat1");
57  auto pfpdgId = input_cat1.allocate<float>();
58  auto& vpfpdgId = (*pfpdgId)[0];
59 
60  auto& input_cat2 = iInput.at("input_cat2");
61  auto pffromPV = input_cat2.allocate<float>();
62  auto& vpffromPV = (*pffromPV)[0];
63 
64  size_t i_pf = 0;
65  for (const auto& pf : pfs) {
66  if (ignore_leptons_) {
67  int pdg_id = std::abs(pf.pdgId());
68  if (pdg_id == 11 || pdg_id == 13) {
69  px_leptons_ += pf.px();
70  py_leptons_ += pf.py();
71  continue;
72  }
73  }
74 
75  // PF keys [b'PF_dxy', b'PF_dz', b'PF_eta', b'PF_mass', b'PF_pt', b'PF_puppiWeight', b'PF_px', b'PF_py']
76  vpfdata.push_back(pf.dxy());
77  vpfdata.push_back(pf.dz());
78  vpfdata.push_back(pf.eta());
79  vpfdata.push_back(pf.mass());
80  vpfdata.push_back(scale_and_rm_outlier(pf.pt(), scale_));
81  vpfdata.push_back(pf.puppiWeight());
82  vpfdata.push_back(scale_and_rm_outlier(pf.px(), scale_));
83  vpfdata.push_back(scale_and_rm_outlier(pf.py(), scale_));
84 
85  vpfchg.push_back(charge_embedding.at(pf.charge()));
86 
87  vpfpdgId.push_back(pdg_id_embedding.at(pf.pdgId()));
88 
89  vpffromPV.push_back(pf.fromPV());
90 
91  ++i_pf;
92  if (i_pf == max_n_pf_) {
93  edm::LogWarning("acquire")
94  << "<DeepMETSonicProducer::acquire>:" << std::endl
95  << " The number of particles is equal to or exceeds the maximum considerable for DeepMET" << std::endl;
96  break;
97  }
98  }
99 
100  // fill the remaining with zeros
101  // resize the vector to 4500 for zero-padding
102  vpfdata.resize(8 * max_n_pf_);
103  vpfchg.resize(max_n_pf_);
104  vpfpdgId.resize(max_n_pf_);
105  vpffromPV.resize(max_n_pf_);
106 
107  input.toServer(pfdata);
108  input_cat0.toServer(pfchg);
109  input_cat1.toServer(pfpdgId);
110  input_cat2.toServer(pffromPV);
111 }
112 
114  const auto& output1 = iOutput.begin()->second;
115  const auto& outputs = output1.fromServer<float>();
116 
117  // outputs are px and py
118  float px = outputs[0][0] * norm_;
119  float py = outputs[0][1] * norm_;
120 
121  // subtract the lepton pt contribution
122  px -= px_leptons_;
123  py -= py_leptons_;
124 
125  LogDebug("produce") << "<DeepMETSonicProducer::produce>:" << std::endl
126  << " MET from DeepMET Sonic Producer is MET_x " << px << " and MET_y " << py << std::endl;
127 
128  auto pf_mets = std::make_unique<pat::METCollection>();
129  const reco::Candidate::LorentzVector p4(px, py, 0., std::hypot(px, py));
130  pf_mets->emplace_back(reco::MET(p4, {}));
131  iEvent.put(std::move(pf_mets));
132 }
133 
137  desc.add<edm::InputTag>("pf_src", edm::InputTag("packedPFCandidates"));
138  desc.add<bool>("ignore_leptons", false);
139  desc.add<double>("norm_factor", 50.);
140  desc.add<unsigned int>("max_n_pf", 4500);
141  descriptions.add("deepMETSonicProducer", desc);
142 }
143 
const unsigned int max_n_pf_
void acquire(edm::Event const &iEvent, edm::EventSetup const &iSetup, Input &iInput) override
void produce(edm::Event &iEvent, edm::EventSetup const &iSetup, Output const &iOutput) override
static const std::unordered_map< int, int32_t > charge_embedding
Definition: DeepMETHelp.h:10
static std::string const input
Definition: EdmProvDump.cc:50
Definition: HeavyIon.h:7
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
Definition: MET.h:41
const edm::EDGetTokenT< std::vector< pat::PackedCandidate > > pf_token_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static const std::unordered_map< int, int32_t > pdg_id_embedding
Definition: DeepMETHelp.h:11
DeepMETSonicProducer(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
HLT enums.
float scale_and_rm_outlier(float val, float scale)
Definition: DeepMETHelper.cc:4
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)