CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTScoutingPFProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTrigger/JetMET
4 // Class: HLTScoutingPFProducer
5 //
11 //
12 // Original Author: Dustin James Anderson
13 // Created: Fri, 12 Jun 2015 15:49:20 GMT
14 //
15 //
16 
17 // system include files
18 #include <memory>
19 
20 // user include files
26 
35 
39 
41 
43 
45 
47 public:
48  explicit HLTScoutingPFProducer(const edm::ParameterSet &);
49  ~HLTScoutingPFProducer() override;
50 
51  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
52 
53 private:
54  void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const final;
55 
62 
63  const double pfJetPtCut_;
64  const double pfJetEtaCut_;
65  const double pfCandidatePtCut_;
66  const double pfCandidateEtaCut_;
67  const int mantissaPrecision_;
68 
69  const bool doJetTags_;
70  const bool doCandidates_;
71  const bool doMet_;
72  const bool doTrackRelVars_;
73  const bool doCandIndsForJets_;
74 };
75 
76 //
77 // constructors and destructor
78 //
80  : pfJetCollection_(consumes(iConfig.getParameter<edm::InputTag>("pfJetCollection"))),
81  pfJetTagCollection_(consumes(iConfig.getParameter<edm::InputTag>("pfJetTagCollection"))),
82  pfCandidateCollection_(consumes(iConfig.getParameter<edm::InputTag>("pfCandidateCollection"))),
83  vertexCollection_(consumes(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
84  metCollection_(consumes(iConfig.getParameter<edm::InputTag>("metCollection"))),
85  rho_(consumes(iConfig.getParameter<edm::InputTag>("rho"))),
86  pfJetPtCut_(iConfig.getParameter<double>("pfJetPtCut")),
87  pfJetEtaCut_(iConfig.getParameter<double>("pfJetEtaCut")),
88  pfCandidatePtCut_(iConfig.getParameter<double>("pfCandidatePtCut")),
89  pfCandidateEtaCut_(iConfig.getParameter<double>("pfCandidateEtaCut")),
90  mantissaPrecision_(iConfig.getParameter<int>("mantissaPrecision")),
91  doJetTags_(iConfig.getParameter<bool>("doJetTags")),
92  doCandidates_(iConfig.getParameter<bool>("doCandidates")),
93  doMet_(iConfig.getParameter<bool>("doMet")),
94  doTrackRelVars_(iConfig.getParameter<bool>("doTrackRelVars")),
95  doCandIndsForJets_(iConfig.getParameter<bool>("doCandIndsForJets")) {
96  //register products
97  produces<Run3ScoutingPFJetCollection>();
98  produces<Run3ScoutingParticleCollection>();
99  produces<double>("rho");
100  produces<double>("pfMetPt");
101  produces<double>("pfMetPhi");
102 }
103 
105 
106 // ------------ method called to produce the data ------------
108  using namespace edm;
109 
110  //get vertices
112  auto outVertices = std::make_unique<Run3ScoutingVertexCollection>();
114  for (auto const &vtx : *vertexCollection) {
115  outVertices->emplace_back(MiniFloatConverter::reduceMantissaToNbitsRounding(vtx.x(), mantissaPrecision_),
121  vtx.tracksSize(),
123  vtx.ndof(),
124  vtx.isValid());
125  }
126  }
127 
128  //get rho
130  auto outRho = std::make_unique<double>(-999);
131  if (iEvent.getByToken(rho_, rho)) {
132  *outRho = *rho;
133  }
134 
135  //get MET
137  auto outMetPt = std::make_unique<double>(-999);
138  auto outMetPhi = std::make_unique<double>(-999);
139  if (doMet_ && iEvent.getByToken(metCollection_, metCollection)) {
140  outMetPt = std::make_unique<double>(metCollection->front().pt());
141  outMetPhi = std::make_unique<double>(metCollection->front().phi());
142  }
143 
144  //get PF candidates
146  auto outPFCandidates = std::make_unique<Run3ScoutingParticleCollection>();
148  for (auto const &cand : *pfCandidateCollection) {
149  if (cand.pt() > pfCandidatePtCut_ && std::abs(cand.eta()) < pfCandidateEtaCut_) {
150  int vertex_index = -1;
151  int index_counter = 0;
152  double dr2 = 0.0001;
153  for (auto const &vtx : *outVertices) {
154  double tmp_dr2 = pow(vtx.x() - cand.vx(), 2) + pow(vtx.y() - cand.vy(), 2) + pow(vtx.z() - cand.vz(), 2);
155  if (tmp_dr2 < dr2) {
156  dr2 = tmp_dr2;
157  vertex_index = index_counter;
158  }
159  if (dr2 == 0.0)
160  break;
161  ++index_counter;
162  }
163  float normchi2{0}, dz{0}, dxy{0}, dzError{0}, dxyError{0}, trk_pt{0}, trk_eta{0}, trk_phi{0};
164  uint8_t lostInnerHits{0}, quality{0};
165  if (doTrackRelVars_) {
166  const auto *trk = cand.bestTrack();
167  if (trk != nullptr) {
169  lostInnerHits = btagbtvdeep::lost_inner_hits_from_pfcand(cand);
174 
175  if (not vertexCollection->empty()) {
176  const reco::Vertex &pv = (*vertexCollection)[0];
177 
178  dz = trk->dz(pv.position());
181 
182  dxy = trk->dxy(pv.position());
185  }
186  } else {
188  }
189  }
190  outPFCandidates->emplace_back(
195  cand.pdgId(),
196  vertex_index,
197  normchi2,
198  dz,
199  dxy,
200  dzError,
201  dxyError,
202  lostInnerHits,
203  quality,
204  trk_pt,
205  trk_eta,
206  trk_phi);
207  }
208  }
209  }
210 
211  //get PF jets
213  auto outPFJets = std::make_unique<Run3ScoutingPFJetCollection>();
215  //get PF jet tags
217  bool haveJetTags = false;
218  if (doJetTags_ && iEvent.getByToken(pfJetTagCollection_, pfJetTagCollection)) {
219  haveJetTags = true;
220  }
221 
222  for (auto const &jet : *pfJetCollection) {
223  if (jet.pt() < pfJetPtCut_ || std::abs(jet.eta()) > pfJetEtaCut_)
224  continue;
225  //find the jet tag corresponding to the jet
226  float tagValue = -20;
227  float minDR2 = 0.01;
228  if (haveJetTags) {
229  for (auto const &tag : *pfJetTagCollection) {
230  float dR2 = reco::deltaR2(jet, *(tag.first));
231  if (dR2 < minDR2) {
232  minDR2 = dR2;
233  tagValue = tag.second;
234  }
235  }
236  }
237  //get the PF constituents of the jet
238  std::vector<int> candIndices;
240  for (auto const &cand : jet.getPFConstituents()) {
241  if (not(cand.isNonnull() and cand.isAvailable())) {
242  throw cms::Exception("HLTScoutingPFProducer")
243  << "invalid reference to reco::PFCandidate from reco::PFJet::getPFConstituents()";
244  }
245  if (cand->pt() > pfCandidatePtCut_ && std::abs(cand->eta()) < pfCandidateEtaCut_) {
246  //search for the candidate in the collection
247  float minDR2 = 0.0001;
248  int matchIndex = -1;
249  int outIndex = 0;
250  for (auto &outCand : *outPFCandidates) {
251  auto const dR2 = reco::deltaR2(*cand, outCand);
252  if (dR2 < minDR2) {
253  minDR2 = dR2;
254  matchIndex = outIndex;
255  }
256  if (minDR2 == 0) {
257  break;
258  }
259  outIndex++;
260  }
261  candIndices.push_back(matchIndex);
262  }
263  }
264  }
265  outPFJets->emplace_back(
278  jet.chargedHadronMultiplicity(),
279  jet.neutralHadronMultiplicity(),
280  jet.photonMultiplicity(),
281  jet.electronMultiplicity(),
282  jet.muonMultiplicity(),
283  jet.HFHadronMultiplicity(),
284  jet.HFEMMultiplicity(),
288  candIndices);
289  }
290  }
291 
292  //put output
293  iEvent.put(std::move(outPFCandidates));
294  iEvent.put(std::move(outPFJets));
295  iEvent.put(std::move(outRho), "rho");
296  iEvent.put(std::move(outMetPt), "pfMetPt");
297  iEvent.put(std::move(outMetPhi), "pfMetPhi");
298 }
299 
300 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
303  desc.add<edm::InputTag>("pfJetCollection", edm::InputTag("hltAK4PFJets"));
304  desc.add<edm::InputTag>("pfJetTagCollection", edm::InputTag("hltDeepCombinedSecondaryVertexBJetTagsPF"));
305  desc.add<edm::InputTag>("pfCandidateCollection", edm::InputTag("hltParticleFlow"));
306  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
307  desc.add<edm::InputTag>("metCollection", edm::InputTag("hltPFMETProducer"));
308  desc.add<edm::InputTag>("rho", edm::InputTag("hltFixedGridRhoFastjetAll"));
309  desc.add<double>("pfJetPtCut", 20.0);
310  desc.add<double>("pfJetEtaCut", 3.0);
311  desc.add<double>("pfCandidatePtCut", 0.6);
312  desc.add<double>("pfCandidateEtaCut", 5.0);
313  desc.add<int>("mantissaPrecision", 10)->setComment("default float16, change to 23 for float32");
314  desc.add<bool>("doJetTags", true);
315  desc.add<bool>("doCandidates", true);
316  desc.add<bool>("doMet", true);
317  desc.add<bool>("doTrackRelVars", true);
318  desc.add<bool>("doCandIndsForJets", false);
319  descriptions.addWithDefaultLabel(desc);
320 }
321 
322 // declare this class as a framework plugin
float quality_from_pfcand(const reco::PFCandidate &pfcand)
Definition: deep_helpers.cc:96
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const edm::EDGetTokenT< reco::PFMETCollection > metCollection_
~HLTScoutingPFProducer() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
uint32_t const *__restrict__ Quality * quality
tuple pfJetTagCollection
const Point & position() const
position
Definition: Vertex.h:127
tuple vertexCollection
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollection_
tuple pfCandidateCollection
const edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
tuple pfJetCollection
int iEvent
Definition: GenABIO.cc:224
void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const final
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::EDGetTokenT< reco::PFJetCollection > pfJetCollection_
HLTScoutingPFProducer(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< reco::JetTagCollection > pfJetTagCollection_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
const edm::EDGetTokenT< double > rho_
float lost_inner_hits_from_pfcand(const reco::PFCandidate &pfcand)
static float reduceMantissaToNbitsRounding(const float &f)
Definition: libminifloat.h:79
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29