CMS 3D CMS Logo

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 public:
46  explicit HLTScoutingPFProducer(const edm::ParameterSet &);
47  ~HLTScoutingPFProducer() override;
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
50 
51 private:
52  void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const final;
53 
60 
61  const double pfJetPtCut;
62  const double pfJetEtaCut;
63  const double pfCandidatePtCut;
64  const double pfCandidateEtaCut;
65  const int mantissaPrecision;
66 
67  const bool doJetTags;
68  const bool doCandidates;
69  const bool doMet;
70 };
71 
72 //
73 // constructors and destructor
74 //
76  : pfJetCollection_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("pfJetCollection"))),
77  pfJetTagCollection_(consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("pfJetTagCollection"))),
78  pfCandidateCollection_(
79  consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandidateCollection"))),
80  vertexCollection_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
81  metCollection_(consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("metCollection"))),
82  rho_(consumes<double>(iConfig.getParameter<edm::InputTag>("rho"))),
83  pfJetPtCut(iConfig.getParameter<double>("pfJetPtCut")),
84  pfJetEtaCut(iConfig.getParameter<double>("pfJetEtaCut")),
85  pfCandidatePtCut(iConfig.getParameter<double>("pfCandidatePtCut")),
86  pfCandidateEtaCut(iConfig.getParameter<double>("pfCandidateEtaCut")),
87  mantissaPrecision(iConfig.getParameter<int>("mantissaPrecision")),
88  doJetTags(iConfig.getParameter<bool>("doJetTags")),
89  doCandidates(iConfig.getParameter<bool>("doCandidates")),
90  doMet(iConfig.getParameter<bool>("doMet")) {
91  //register products
92  produces<ScoutingPFJetCollection>();
93  produces<ScoutingParticleCollection>();
94  produces<double>("rho");
95  produces<double>("pfMetPt");
96  produces<double>("pfMetPhi");
97 }
98 
100 
101 // ------------ method called to produce the data ------------
103  using namespace edm;
104 
105  //get vertices
107  std::unique_ptr<ScoutingVertexCollection> outVertices(new ScoutingVertexCollection());
108  if (iEvent.getByToken(vertexCollection_, vertexCollection)) {
109  for (auto &vtx : *vertexCollection) {
110  outVertices->emplace_back(vtx.x(),
111  vtx.y(),
112  vtx.z(),
113  vtx.zError(),
114  vtx.xError(),
115  vtx.yError(),
116  vtx.tracksSize(),
117  vtx.chi2(),
118  vtx.ndof(),
119  vtx.isValid());
120  }
121  }
122 
123  //get rho
125  std::unique_ptr<double> outRho(new double(-999));
126  if (iEvent.getByToken(rho_, rho)) {
127  outRho.reset(new double(*rho));
128  }
129 
130  //get MET
132  std::unique_ptr<double> outMetPt(new double(-999));
133  std::unique_ptr<double> outMetPhi(new double(-999));
134  if (doMet && iEvent.getByToken(metCollection_, metCollection)) {
135  outMetPt.reset(new double(metCollection->front().pt()));
136  outMetPhi.reset(new double(metCollection->front().phi()));
137  }
138 
139  //get PF candidates
141  std::unique_ptr<ScoutingParticleCollection> outPFCandidates(new ScoutingParticleCollection());
143  for (auto &cand : *pfCandidateCollection) {
144  if (cand.pt() > pfCandidatePtCut && std::abs(cand.eta()) < pfCandidateEtaCut) {
145  int vertex_index = -1;
146  int index_counter = 0;
147  double dr2 = 0.0001;
148  for (auto &vtx : *outVertices) {
149  double tmp_dr2 = pow(vtx.x() - cand.vx(), 2) + pow(vtx.y() - cand.vy(), 2) + pow(vtx.z() - cand.vz(), 2);
150  if (tmp_dr2 < dr2) {
151  dr2 = tmp_dr2;
152  vertex_index = index_counter;
153  }
154  if (dr2 == 0.0)
155  break;
156  ++index_counter;
157  }
158 
159  outPFCandidates->emplace_back(MiniFloatConverter::reduceMantissaToNbitsRounding(cand.pt(), mantissaPrecision),
163  cand.pdgId(),
164  vertex_index);
165  }
166  }
167  }
168 
169  //get PF jets
171  std::unique_ptr<ScoutingPFJetCollection> outPFJets(new ScoutingPFJetCollection());
172  if (iEvent.getByToken(pfJetCollection_, pfJetCollection)) {
173  //get PF jet tags
175  bool haveJetTags = false;
177  haveJetTags = true;
178  }
179 
180  for (auto &jet : *pfJetCollection) {
181  if (jet.pt() < pfJetPtCut || std::abs(jet.eta()) > pfJetEtaCut)
182  continue;
183  //find the jet tag corresponding to the jet
184  float tagValue = -20;
185  float minDR2 = 0.01;
186  if (haveJetTags) {
187  for (auto &tag : *pfJetTagCollection) {
188  float dR2 = reco::deltaR2(jet, *(tag.first));
189  if (dR2 < minDR2) {
190  minDR2 = dR2;
191  tagValue = tag.second;
192  }
193  }
194  }
195  //get the PF constituents of the jet
196  std::vector<int> candIndices;
197  if (doCandidates) {
198  for (auto &cand : jet.getPFConstituents()) {
199  if (cand->pt() > pfCandidatePtCut && std::abs(cand->eta()) < pfCandidateEtaCut) {
200  //search for the candidate in the collection
201  float minDR2 = 0.0001;
202  int matchIndex = -1;
203  int outIndex = 0;
204  for (auto &outCand : *outPFCandidates) {
205  float dR2 = pow(cand->eta() - outCand.eta(), 2) + pow(cand->phi() - outCand.phi(), 2);
206  if (dR2 < minDR2) {
207  minDR2 = dR2;
208  matchIndex = outIndex;
209  }
210  if (minDR2 == 0) {
211  break;
212  }
213  outIndex++;
214  }
215  candIndices.push_back(matchIndex);
216  }
217  }
218  }
219  outPFJets->emplace_back(jet.pt(),
220  jet.eta(),
221  jet.phi(),
222  jet.mass(),
223  jet.jetArea(),
224  jet.chargedHadronEnergy(),
225  jet.neutralHadronEnergy(),
226  jet.photonEnergy(),
227  jet.electronEnergy(),
228  jet.muonEnergy(),
229  jet.HFHadronEnergy(),
230  jet.HFEMEnergy(),
231  jet.chargedHadronMultiplicity(),
232  jet.neutralHadronMultiplicity(),
233  jet.photonMultiplicity(),
234  jet.electronMultiplicity(),
235  jet.muonMultiplicity(),
236  jet.HFHadronMultiplicity(),
237  jet.HFEMMultiplicity(),
238  jet.hoEnergy(),
239  tagValue,
240  0.0,
241  candIndices);
242  }
243  }
244 
245  //put output
246  iEvent.put(std::move(outPFCandidates));
247  iEvent.put(std::move(outPFJets));
248  iEvent.put(std::move(outRho), "rho");
249  iEvent.put(std::move(outMetPt), "pfMetPt");
250  iEvent.put(std::move(outMetPhi), "pfMetPhi");
251 }
252 
253 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
256  desc.add<edm::InputTag>("pfJetCollection", edm::InputTag("hltAK4PFJets"));
257  desc.add<edm::InputTag>("pfJetTagCollection", edm::InputTag("hltDeepCombinedSecondaryVertexBJetTagsPF"));
258  desc.add<edm::InputTag>("pfCandidateCollection", edm::InputTag("hltParticleFlow"));
259  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
260  desc.add<edm::InputTag>("metCollection", edm::InputTag("hltPFMETProducer"));
261  desc.add<edm::InputTag>("rho", edm::InputTag("hltFixedGridRhoFastjetAll"));
262  desc.add<double>("pfJetPtCut", 20.0);
263  desc.add<double>("pfJetEtaCut", 3.0);
264  desc.add<double>("pfCandidatePtCut", 0.6);
265  desc.add<double>("pfCandidateEtaCut", 5.0);
266  desc.add<int>("mantissaPrecision", 23);
267  desc.add<bool>("doJetTags", true);
268  desc.add<bool>("doCandidates", true);
269  desc.add<bool>("doMet", true);
270  descriptions.add("hltScoutingPFProducer", desc);
271 }
272 
273 // declare this class as a framework plugin
edm::StreamID
Definition: StreamID.h:30
JetTag.h
HLTScoutingPFProducer::pfCandidatePtCut
const double pfCandidatePtCut
Definition: HLTScoutingPFProducer.cc:63
electrons_cff.bool
bool
Definition: electrons_cff.py:372
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
HLTScoutingPFProducer
Definition: HLTScoutingPFProducer.cc:44
HLTScoutingPFProducer::pfJetEtaCut
const double pfJetEtaCut
Definition: HLTScoutingPFProducer.cc:62
HLTScoutingPFProducer::doCandidates
const bool doCandidates
Definition: HLTScoutingPFProducer.cc:68
HLTScoutingPFProducer::~HLTScoutingPFProducer
~HLTScoutingPFProducer() override
PFCandidate.h
susyDQM_cfi.metCollection
metCollection
Definition: susyDQM_cfi.py:11
edm::EDGetTokenT< reco::PFJetCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
PFMETCollection
Collection of PF MET.
HLT_2018_cff.doCandidates
doCandidates
Definition: HLT_2018_cff.py:88444
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
PFJet.h
HLTScoutingPFProducer::mantissaPrecision
const int mantissaPrecision
Definition: HLTScoutingPFProducer.cc:65
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
HLTScoutingPFProducer::pfJetCollection_
const edm::EDGetTokenT< reco::PFJetCollection > pfJetCollection_
Definition: HLTScoutingPFProducer.cc:54
HLT_2018_cff.pfCandidatePtCut
pfCandidatePtCut
Definition: HLT_2018_cff.py:88456
HLT_2018_cff.pfJetPtCut
pfJetPtCut
Definition: HLT_2018_cff.py:88451
HLTScoutingPFProducer::produce
void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const final
Definition: HLTScoutingPFProducer.cc:102
PFMETCollection.h
HLT_2018_cff.mantissaPrecision
mantissaPrecision
Definition: HLT_2018_cff.py:88452
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
ScoutingPFJet.h
edm::Handle< reco::VertexCollection >
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
HLTScoutingPFProducer::pfJetPtCut
const double pfJetPtCut
Definition: HLTScoutingPFProducer.cc:61
deltaR.h
HLTScoutingPFProducer::pfCandidateEtaCut
const double pfCandidateEtaCut
Definition: HLTScoutingPFProducer.cc:64
MakerMacros.h
HLTScoutingPFProducer::pfJetTagCollection_
const edm::EDGetTokenT< reco::JetTagCollection > pfJetTagCollection_
Definition: HLTScoutingPFProducer.cc:55
HLTScoutingPFProducer::metCollection_
const edm::EDGetTokenT< reco::PFMETCollection > metCollection_
Definition: HLTScoutingPFProducer.cc:58
HLTScoutingPFProducer::HLTScoutingPFProducer
HLTScoutingPFProducer(const edm::ParameterSet &)
Definition: HLTScoutingPFProducer.cc:75
HLT_2018_cff.doJetTags
doJetTags
Definition: HLT_2018_cff.py:88454
libminifloat.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
GlobalPosition_Frontier_DevDB_cff.tag
tag
Definition: GlobalPosition_Frontier_DevDB_cff.py:11
ScoutingVertex.h
ExoticaDQM_cfi.pfJetCollection
pfJetCollection
Definition: ExoticaDQM_cfi.py:19
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
DDAxes::rho
badGlobalMuonTaggersAOD_cff.vtx
vtx
Definition: badGlobalMuonTaggersAOD_cff.py:5
Vertex.h
HLT_2018_cff.pfJetEtaCut
pfJetEtaCut
Definition: HLT_2018_cff.py:88449
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
MiniFloatConverter::reduceMantissaToNbitsRounding
static float reduceMantissaToNbitsRounding(const float &f)
Definition: libminifloat.h:102
edm::ParameterSet
Definition: ParameterSet.h:36
Event.h
ScoutingVertexCollection
std::vector< ScoutingVertex > ScoutingVertexCollection
Definition: ScoutingVertex.h:69
HLT_2018_cff.pfCandidateEtaCut
pfCandidateEtaCut
Definition: HLT_2018_cff.py:88445
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
PFMET.h
cand
Definition: decayParser.h:34
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::JetTagCollection
JetFloatAssociation::Container JetTagCollection
Definition: JetTag.h:17
edm::EventSetup
Definition: EventSetup.h:57
ScoutingParticle.h
ScoutingPFJetCollection
std::vector< ScoutingPFJet > ScoutingPFJetCollection
Definition: ScoutingPFJet.h:134
l1t::PFCandidateCollection
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:57
VertexFwd.h
HLTScoutingPFProducer::vertexCollection_
const edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
Definition: HLTScoutingPFProducer.cc:57
HLTScoutingPFProducer::doJetTags
const bool doJetTags
Definition: HLTScoutingPFProducer.cc:67
eostools.move
def move(src, dest)
Definition: eostools.py:511
reco::PFJetCollection
std::vector< PFJet > PFJetCollection
collection of PFJet objects
Definition: PFJetCollection.h:14
HLTScoutingPFProducer::rho_
const edm::EDGetTokenT< double > rho_
Definition: HLTScoutingPFProducer.cc:59
Frameworkfwd.h
spclusmultinvestigator_cfi.vertexCollection
vertexCollection
Definition: spclusmultinvestigator_cfi.py:4
metsig::jet
Definition: SignAlgoResolutions.h:47
HLTScoutingPFProducer::doMet
const bool doMet
Definition: HLTScoutingPFProducer.cc:69
HLT_2018_cff.pfJetTagCollection
pfJetTagCollection
Definition: HLT_2018_cff.py:88446
HLTScoutingPFProducer::pfCandidateCollection_
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollection_
Definition: HLTScoutingPFProducer.cc:56
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HLT_2018_cff.pfCandidateCollection
pfCandidateCollection
Definition: HLT_2018_cff.py:88455
ParameterSet.h
EDProducer.h
edm::Event
Definition: Event.h:73
HLT_2018_cff.doMet
doMet
Definition: HLT_2018_cff.py:88447
edm::InputTag
Definition: InputTag.h:15
PFCandidateFwd.h
ScoutingParticleCollection
std::vector< ScoutingParticle > ScoutingParticleCollection
Definition: ScoutingParticle.h:33
HLTScoutingPFProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTScoutingPFProducer.cc:254