CMS 3D CMS Logo

GenJetFlavourTableProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
7 
9 
13 
18 
20 
23 
25 public:
27  : name_(iConfig.getParameter<std::string>("name")),
28  src_(consumes<std::vector<reco::GenJet> >(iConfig.getParameter<edm::InputTag>("src"))),
29  cut_(iConfig.getParameter<std::string>("cut"), true),
30  deltaR_(iConfig.getParameter<double>("deltaR")),
32  consumes<reco::JetFlavourInfoMatchingCollection>(iConfig.getParameter<edm::InputTag>("jetFlavourInfos"))) {
33  produces<nanoaod::FlatTable>();
34  }
35 
37 
38  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
40  desc.add<edm::InputTag>("src")->setComment("input genJet collection");
41  desc.add<edm::InputTag>("jetFlavourInfos")->setComment("input flavour info collection");
42  desc.add<std::string>("name")->setComment("name of the genJet FlatTable we are extending with flavour information");
43  desc.add<std::string>("cut")->setComment("cut on input genJet collection");
44  desc.add<double>("deltaR")->setComment("deltaR to match genjets");
45  descriptions.add("genJetFlavourTable", desc);
46  }
47 
48 private:
49  void produce(edm::Event&, edm::EventSetup const&) override;
50 
54  const double deltaR_;
56 };
57 
58 // ------------ method called to produce the data ------------
60  const auto& jetsProd = iEvent.get(src_);
61  const auto& jetFlavourInfosProd = iEvent.get(jetFlavourInfosToken_);
62 
63  unsigned int ncand = 0;
64  std::vector<int> partonFlavour;
65  std::vector<uint8_t> hadronFlavour;
66 
67  for (const reco::GenJet& jet : jetsProd) {
68  if (!cut_(jet))
69  continue;
70  ++ncand;
71  bool matched = false;
72  for (const reco::JetFlavourInfoMatching& jetFlavourInfoMatching : jetFlavourInfosProd) {
73  if (deltaR(jet.p4(), jetFlavourInfoMatching.first->p4()) < deltaR_) {
74  partonFlavour.push_back(jetFlavourInfoMatching.second.getPartonFlavour());
75  hadronFlavour.push_back(jetFlavourInfoMatching.second.getHadronFlavour());
76  matched = true;
77  break;
78  }
79  }
80  if (!matched) {
81  partonFlavour.push_back(0);
82  hadronFlavour.push_back(0);
83  }
84  }
85 
86  auto tab = std::make_unique<nanoaod::FlatTable>(ncand, name_, false, true);
87  tab->addColumn<int>("partonFlavour", partonFlavour, "flavour from parton matching");
88  tab->addColumn<uint8_t>("hadronFlavour", hadronFlavour, "flavour from hadron ghost clustering");
89 
90  iEvent.put(std::move(tab));
91 }
92 
94 //define this as a plug-in
GenJetFlavourTableProducer(const edm::ParameterSet &iConfig)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::JetFlavourInfoMatchingCollection > jetFlavourInfosToken_
edm::EDGetTokenT< std::vector< reco::GenJet > > src_
Jets made from MC generator particles.
Definition: GenJet.h:23
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, edm::EventSetup const &) override
JetFlavourInfoMatchingCollection::value_type JetFlavourInfoMatching
const StringCutObjectSelector< reco::GenJet > cut_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
def move(src, dest)
Definition: eostools.py:511