CMS 3D CMS Logo

BoostedJetMerger.cc
Go to the documentation of this file.
2 
4  : jetToken_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("jetSrc"))),
5  subjetToken_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("subjetSrc"))) {
6  //register products
7  produces<std::vector<pat::Jet>>();
8  produces<std::vector<pat::Jet>>("SubJets");
9 }
10 
12 
13 // ------------ method called to produce the data ------------
15  auto outputs = std::make_unique<std::vector<pat::Jet>>();
16  auto outputSubjets = std::make_unique<std::vector<pat::Jet>>();
17 
18  edm::RefProd<std::vector<pat::Jet>> h_subJetsOut = iEvent.getRefBeforePut<std::vector<pat::Jet>>("SubJets");
19 
21  edm::Handle<edm::View<pat::Jet>> subjetHandle;
22 
23  iEvent.getByToken(jetToken_, jetHandle);
24  iEvent.getByToken(subjetToken_, subjetHandle);
25 
26  for (edm::View<pat::Jet>::const_iterator ijetBegin = jetHandle->begin(), ijetEnd = jetHandle->end(), ijet = ijetBegin;
27  ijet != ijetEnd;
28  ++ijet) {
29  outputs->push_back(*ijet);
30  std::vector<edm::Ptr<reco::Candidate>> nextSubjets;
31 
32  for (unsigned int isubjet = 0; isubjet < ijet->numberOfDaughters(); ++isubjet) {
33  edm::Ptr<reco::Candidate> const& subjet = ijet->daughterPtr(isubjet);
35  find_if(subjetHandle->begin(), subjetHandle->end(), FindCorrectedSubjet(subjet));
36  if (ifound != subjetHandle->end()) {
37  outputSubjets->push_back(*ifound);
38 
39  edm::Ref<std::vector<pat::Jet>> subjetRef(h_subJetsOut, outputSubjets->size() - 1);
40  edm::Ptr<pat::Jet> subjetPtr(h_subJetsOut.id(), subjetRef.key(), h_subJetsOut.productGetter());
41  nextSubjets.push_back(subjetPtr);
42  }
43  }
44  outputs->back().clearDaughters();
45  for (std::vector<edm::Ptr<reco::Candidate>>::const_iterator nextSubjet = nextSubjets.begin(),
46  nextSubjetEnd = nextSubjets.end();
47  nextSubjet != nextSubjetEnd;
48  ++nextSubjet) {
49  outputs->back().addDaughter(*nextSubjet);
50  }
51  }
52 
53  iEvent.put(std::move(outputs));
54  iEvent.put(std::move(outputSubjets), "SubJets");
55 }
56 
57 //define this as a plug-in
~BoostedJetMerger() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
BoostedJetMerger(const edm::ParameterSet &)
Definition: Jet.py:1
edm::EDGetTokenT< edm::View< pat::Jet > > jetToken_
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< edm::View< pat::Jet > > subjetToken_
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511