File indexing completed on 2024-04-06 12:24:04
0001 #include "PhysicsTools/PatUtils/interface/BoostedJetMerger.h"
0002
0003 BoostedJetMerger::BoostedJetMerger(const edm::ParameterSet& iConfig)
0004 : jetToken_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("jetSrc"))),
0005 subjetToken_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("subjetSrc"))) {
0006
0007 produces<std::vector<pat::Jet>>();
0008 produces<std::vector<pat::Jet>>("SubJets");
0009 }
0010
0011 BoostedJetMerger::~BoostedJetMerger() {}
0012
0013
0014 void BoostedJetMerger::produce(edm::Event& iEvent, const edm::EventSetup&) {
0015 auto outputs = std::make_unique<std::vector<pat::Jet>>();
0016 auto outputSubjets = std::make_unique<std::vector<pat::Jet>>();
0017
0018 edm::RefProd<std::vector<pat::Jet>> h_subJetsOut = iEvent.getRefBeforePut<std::vector<pat::Jet>>("SubJets");
0019
0020 edm::Handle<edm::View<pat::Jet>> jetHandle;
0021 edm::Handle<edm::View<pat::Jet>> subjetHandle;
0022
0023 iEvent.getByToken(jetToken_, jetHandle);
0024 iEvent.getByToken(subjetToken_, subjetHandle);
0025
0026 for (edm::View<pat::Jet>::const_iterator ijetBegin = jetHandle->begin(), ijetEnd = jetHandle->end(), ijet = ijetBegin;
0027 ijet != ijetEnd;
0028 ++ijet) {
0029 outputs->push_back(*ijet);
0030 std::vector<edm::Ptr<reco::Candidate>> nextSubjets;
0031
0032 for (unsigned int isubjet = 0; isubjet < ijet->numberOfDaughters(); ++isubjet) {
0033 edm::Ptr<reco::Candidate> const& subjet = ijet->daughterPtr(isubjet);
0034 edm::View<pat::Jet>::const_iterator ifound =
0035 find_if(subjetHandle->begin(), subjetHandle->end(), FindCorrectedSubjet(subjet));
0036 if (ifound != subjetHandle->end()) {
0037 outputSubjets->push_back(*ifound);
0038
0039 edm::Ref<std::vector<pat::Jet>> subjetRef(h_subJetsOut, outputSubjets->size() - 1);
0040 edm::Ptr<pat::Jet> subjetPtr(h_subJetsOut.id(), subjetRef.key(), h_subJetsOut.productGetter());
0041 nextSubjets.push_back(subjetPtr);
0042 }
0043 }
0044 outputs->back().clearDaughters();
0045 for (std::vector<edm::Ptr<reco::Candidate>>::const_iterator nextSubjet = nextSubjets.begin(),
0046 nextSubjetEnd = nextSubjets.end();
0047 nextSubjet != nextSubjetEnd;
0048 ++nextSubjet) {
0049 outputs->back().addDaughter(*nextSubjet);
0050 }
0051 }
0052
0053 iEvent.put(std::move(outputs));
0054 iEvent.put(std::move(outputSubjets), "SubJets");
0055 }
0056
0057
0058 DEFINE_FWK_MODULE(BoostedJetMerger);