CMS 3D CMS Logo

GenJetFlavourInfoPreserver.cc
Go to the documentation of this file.
1 
14 
19 
22 
23 namespace pat {
24 
26  public:
27  explicit GenJetFlavourInfoPreserver(const edm::ParameterSet& iConfig);
29 
30  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
31 
32  private:
35 
38  };
39 
40 } // namespace pat
41 
43  : genJetsToken_(consumes<edm::View<reco::GenJet> >(iConfig.getParameter<edm::InputTag>("genJets"))),
44  slimmedGenJetsToken_(consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("slimmedGenJets"))),
45  genJetFlavourInfosToken_(
46  consumes<reco::JetFlavourInfoMatchingCollection>(iConfig.getParameter<edm::InputTag>("genJetFlavourInfos"))),
47  slimmedGenJetAssociationToken_(consumes<edm::Association<std::vector<reco::GenJet> > >(
48  iConfig.getParameter<edm::InputTag>("slimmedGenJetAssociation"))) {
49  produces<reco::JetFlavourInfoMatchingCollection>();
50 }
51 
53  using namespace edm;
54  using namespace std;
55 
57  iEvent.getByToken(genJetsToken_, genJets);
58 
60  iEvent.getByToken(slimmedGenJetsToken_, slimmedGenJets);
61 
63  iEvent.getByToken(genJetFlavourInfosToken_, genJetFlavourInfos);
64 
66  iEvent.getByToken(slimmedGenJetAssociationToken_, slimmedGenJetAssociation);
67 
68  auto jetFlavourInfos = std::make_unique<reco::JetFlavourInfoMatchingCollection>(reco::JetRefBaseProd(slimmedGenJets));
69  assert(genJets->size() == genJetFlavourInfos->size());
70 
71  edm::Ref<std::vector<reco::GenJet> > slimmedGenJetRef;
72 
73  for (unsigned int i = 0; i < genJets->size(); ++i) {
74  slimmedGenJetRef = (*slimmedGenJetAssociation)[genJets->refAt(i)];
75  if (!slimmedGenJetRef)
76  continue;
77  (*jetFlavourInfos)[reco::JetBaseRef(slimmedGenJetRef)] = (*genJetFlavourInfos)[i].second;
78  }
79 
81 }
82 
84 using namespace pat;
const edm::EDGetTokenT< edm::View< reco::Jet > > slimmedGenJetsToken_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Transfers the JetFlavourInfos from the original GenJets to the slimmedGenJets in MiniAOD.
assert(be >=bs)
Definition: HeavyIon.h:7
edm::RefToBaseProd< reco::Jet > JetRefBaseProd
Definition: JetCollection.h:13
int iEvent
Definition: GenABIO.cc:224
GenJetFlavourInfoPreserver(const edm::ParameterSet &iConfig)
Definition: Jet.py:1
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::RefToBase< Jet > JetBaseRef
Definition: JetCollection.h:12
const edm::EDGetTokenT< reco::JetFlavourInfoMatchingCollection > genJetFlavourInfosToken_
const edm::EDGetTokenT< edm::View< reco::GenJet > > genJetsToken_
fixed size matrix
HLT enums.
const edm::EDGetTokenT< edm::Association< std::vector< reco::GenJet > > > slimmedGenJetAssociationToken_
def move(src, dest)
Definition: eostools.py:511