File indexing completed on 2024-04-06 12:25:29
0001
0002
0003
0004
0005
0006
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "DataFormats/JetReco/interface/JetExtendedAssociation.h"
0012
0013 #include "JetExtender.h"
0014
0015 JetExtender::JetExtender(const edm::ParameterSet& fConfig)
0016 : mJets(fConfig.getParameter<edm::InputTag>("jets")),
0017 mJet2TracksAtVX(fConfig.getParameter<edm::InputTag>("jet2TracksAtVX")),
0018 mJet2TracksAtCALO(fConfig.getParameter<edm::InputTag>("jet2TracksAtCALO")) {
0019 token_mJets = consumes<edm::View<reco::Jet> >(mJets);
0020 if (!(mJet2TracksAtVX.label().empty()))
0021 token_mJet2TracksAtVX = consumes<reco::JetTracksAssociation::Container>(mJet2TracksAtVX);
0022 if (!(mJet2TracksAtCALO.label().empty()))
0023 token_mJet2TracksAtCALO = consumes<reco::JetTracksAssociation::Container>(mJet2TracksAtCALO);
0024
0025 produces<reco::JetExtendedAssociation::Container>();
0026 }
0027
0028 JetExtender::~JetExtender() {}
0029
0030 void JetExtender::produce(edm::Event& fEvent, const edm::EventSetup& fSetup) {
0031
0032 edm::Handle<edm::View<reco::Jet> > jets_h;
0033 fEvent.getByToken(token_mJets, jets_h);
0034 edm::Handle<reco::JetTracksAssociation::Container> j2tVX_h;
0035 if (!(mJet2TracksAtVX.label().empty()))
0036 fEvent.getByToken(token_mJet2TracksAtVX, j2tVX_h);
0037 edm::Handle<reco::JetTracksAssociation::Container> j2tCALO_h;
0038 if (!(mJet2TracksAtCALO.label().empty()))
0039 fEvent.getByToken(token_mJet2TracksAtCALO, j2tCALO_h);
0040
0041 auto jetExtender = std::make_unique<reco::JetExtendedAssociation::Container>(reco::JetRefBaseProd(jets_h));
0042
0043
0044
0045 for (unsigned j = 0; j < jets_h->size(); ++j) {
0046 edm::RefToBase<reco::Jet> jet = jets_h->refAt(j);
0047 reco::JetExtendedAssociation::JetExtendedData extendedData;
0048 if (j2tVX_h.isValid()) {
0049 try {
0050 extendedData.mTracksAtVertexNumber = reco::JetTracksAssociation::tracksNumber(*j2tVX_h, jet);
0051 extendedData.mTracksAtVertexP4 = reco::JetTracksAssociation::tracksP4(*j2tVX_h, jet);
0052 } catch (cms::Exception const&) {
0053 edm::LogError("MismatchedJets") << "Jets in original collection " << mJets
0054 << " mismatch jets in j2t VX association " << mJet2TracksAtVX
0055 << ". Wrong collections?";
0056 throw;
0057 }
0058 }
0059 if (j2tCALO_h.isValid()) {
0060 try {
0061 extendedData.mTracksAtCaloNumber = reco::JetTracksAssociation::tracksNumber(*j2tCALO_h, jet);
0062 extendedData.mTracksAtCaloP4 = reco::JetTracksAssociation::tracksP4(*j2tCALO_h, jet);
0063 } catch (cms::Exception const&) {
0064 edm::LogError("MismatchedJets") << "Jets in original collection " << mJets
0065 << " mismatch jets in j2t CALO association " << mJet2TracksAtCALO
0066 << ". Wrong collections?";
0067 throw;
0068 }
0069 }
0070 reco::JetExtendedAssociation::setValue(&*jetExtender, jet, extendedData);
0071 }
0072 fEvent.put(std::move(jetExtender));
0073 }