GeneratorInterface
GenFilters
plugins
NJetsMC.cc
Go to the documentation of this file.
1
// -*- C++ -*-
2
//
3
// Package: NJetsMC
4
// Class: NJetsMC
5
//
13
//
14
// Original Author: "Nathaniel Odell"
15
// Created: Thu Aug 12 09:24:46 CDT 2010
16
// then moved to more general N-jets purpose in GeneratorInterface/GenFilters
17
//
18
19
#include "
DataFormats/Common/interface/Handle.h
"
20
#include "
FWCore/Framework/interface/Event.h
"
21
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
22
#include "
FWCore/Framework/interface/global/EDFilter.h
"
23
#include "
FWCore/Framework/interface/MakerMacros.h
"
24
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
25
#include "
FWCore/Utilities/interface/EDGetToken.h
"
26
#include "
FWCore/Utilities/interface/InputTag.h
"
27
28
#include "
DataFormats/JetReco/interface/GenJetCollection.h
"
29
#include "
DataFormats/JetReco/interface/GenJet.h
"
30
31
#include <cstdint>
32
33
class
NJetsMC
:
public
edm::global::EDFilter
<> {
34
public
:
35
explicit
NJetsMC
(
const
edm::ParameterSet
&);
36
37
private
:
38
bool
filter
(
edm::StreamID
,
edm::Event
&,
const
edm::EventSetup
&)
const override
;
39
40
// ----------member data ---------------------------
41
42
const
edm::EDGetTokenT<reco::GenJetCollection>
genToken_
;
43
const
int
njets_
;
44
const
double
minpt_
;
45
};
46
47
NJetsMC::NJetsMC
(
const
edm::ParameterSet
& iConfig)
48
: genToken_(consumes<
reco
::
GenJetCollection
>(iConfig.getUntrackedParameter<
edm
::
InputTag
>(
"GenTag"
))),
49
njets_(iConfig.getParameter<int32_t>(
"Njets"
)),
50
minpt_(iConfig.getParameter<double>(
"MinPt"
)) {}
51
52
bool
NJetsMC::filter
(
edm::StreamID
,
edm::Event
&
iEvent
,
const
edm::EventSetup
&)
const
{
53
edm::Handle<reco::GenJetCollection>
genJets
;
54
iEvent
.getByToken(
genToken_
,
genJets
);
55
56
int
count
= 0;
57
bool
result
=
false
;
58
59
for
(reco::GenJetCollection::const_iterator iJet =
genJets
->begin(); iJet !=
genJets
->end(); ++iJet) {
60
reco::GenJet
myJet =
reco::GenJet
(*iJet);
61
62
if
(myJet.
pt
() >
minpt_
)
63
++
count
;
64
}
65
66
if
(
count
>=
njets_
)
67
result
=
true
;
68
69
return
result
;
70
}
71
//define this as a plug-in
72
DEFINE_FWK_MODULE
(
NJetsMC
);
edm::StreamID
Definition:
StreamID.h:30
GenJetCollection.h
NJetsMC::njets_
const int njets_
Definition:
NJetsMC.cc:43
Handle.h
reco::GenJet
Jets made from MC generator particles.
Definition:
GenJet.h:23
reco::GenJetCollection
std::vector< GenJet > GenJetCollection
collection of GenJet objects
Definition:
GenJetCollection.h:14
edm::EDGetTokenT< reco::GenJetCollection >
edm
HLT enums.
Definition:
AlignableModifier.h:19
NJetsMC::minpt_
const double minpt_
Definition:
NJetsMC.cc:44
HLT_FULL_cff.InputTag
InputTag
Definition:
HLT_FULL_cff.py:89285
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition:
LeafCandidate.h:146
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
edm::Handle< reco::GenJetCollection >
NJetsMC
Definition:
NJetsMC.cc:33
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
NJetsMC::genToken_
const edm::EDGetTokenT< reco::GenJetCollection > genToken_
Definition:
NJetsMC.cc:42
submitPVResolutionJobs.count
count
Definition:
submitPVResolutionJobs.py:352
EDGetToken.h
NJetsMC::filter
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition:
NJetsMC.cc:52
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
iEvent
int iEvent
Definition:
GenABIO.cc:224
ttbarCategorization_cff.genJets
genJets
Definition:
ttbarCategorization_cff.py:29
edm::EventSetup
Definition:
EventSetup.h:58
InputTag.h
Frameworkfwd.h
NJetsMC::NJetsMC
NJetsMC(const edm::ParameterSet &)
Definition:
NJetsMC.cc:47
edm::global::EDFilter
Definition:
EDFilter.h:32
mps_fire.result
result
Definition:
mps_fire.py:311
GenJet.h
nanoDQM_cfi.GenJet
GenJet
Definition:
nanoDQM_cfi.py:262
ParameterSet.h
edm::Event
Definition:
Event.h:73
EDFilter.h
Generated for CMSSW Reference Manual by
1.8.16