File indexing completed on 2024-04-06 12:31:05
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "DataFormats/Common/interface/EDProductfwd.h"
0010 #include "DataFormats/Common/interface/View.h"
0011 #include "DataFormats/JetReco/interface/Jet.h"
0012 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0013 #include "DataFormats/TrackReco/interface/Track.h"
0014
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 #include "FWCore/Framework/interface/stream/EDProducer.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "FWCore/Utilities/interface/InputTag.h"
0021
0022 #include "SimTracker/TrackHistory/interface/JetVetoedTracksAssociatorDRVertex.h"
0023
0024 class JetVetoedTracksAssociatorAtVertex : public edm::stream::EDProducer<> {
0025 public:
0026 JetVetoedTracksAssociatorAtVertex(const edm::ParameterSet &);
0027 ~JetVetoedTracksAssociatorAtVertex() override;
0028 void produce(edm::Event &, const edm::EventSetup &) override;
0029
0030 private:
0031 edm::EDGetTokenT<edm::View<reco::Jet>> mJets;
0032 edm::EDGetTokenT<reco::TrackCollection> mTracks;
0033 JetVetoedTracksAssociationDRVertex mAssociator;
0034 TrackClassifier classifier;
0035 };
0036
0037 JetVetoedTracksAssociatorAtVertex::JetVetoedTracksAssociatorAtVertex(const edm::ParameterSet &fConfig)
0038 : mJets(consumes<edm::View<reco::Jet>>(fConfig.getParameter<edm::InputTag>("jets"))),
0039 mTracks(consumes<reco::TrackCollection>(fConfig.getParameter<edm::InputTag>("tracks"))),
0040 mAssociator(fConfig.getParameter<double>("coneSize")),
0041 classifier(fConfig, consumesCollector()) {
0042 produces<reco::JetTracksAssociation::Container>();
0043 }
0044
0045 JetVetoedTracksAssociatorAtVertex::~JetVetoedTracksAssociatorAtVertex() {}
0046
0047 void JetVetoedTracksAssociatorAtVertex::produce(edm::Event &fEvent, const edm::EventSetup &fSetup) {
0048
0049 classifier.newEvent(fEvent, fSetup);
0050
0051 edm::Handle<edm::View<reco::Jet>> jets_h;
0052 fEvent.getByToken(mJets, jets_h);
0053 edm::Handle<reco::TrackCollection> tracks_h;
0054 fEvent.getByToken(mTracks, tracks_h);
0055
0056 std::unique_ptr<reco::JetTracksAssociation::Container> jetTracks(
0057 new reco::JetTracksAssociation::Container(reco::JetRefBaseProd(jets_h)));
0058
0059
0060 std::vector<edm::RefToBase<reco::Jet>> allJets;
0061 allJets.reserve(jets_h->size());
0062 for (unsigned i = 0; i < jets_h->size(); ++i)
0063 allJets.push_back(jets_h->refAt(i));
0064 std::vector<reco::TrackRef> allTracks;
0065 allTracks.reserve(tracks_h->size());
0066 for (unsigned i = 0; i < tracks_h->size(); ++i)
0067 allTracks.push_back(reco::TrackRef(tracks_h, i));
0068
0069 mAssociator.produce(&*jetTracks, allJets, allTracks, classifier);
0070
0071 fEvent.put(std::move(jetTracks));
0072 }
0073
0074 DEFINE_FWK_MODULE(JetVetoedTracksAssociatorAtVertex);