Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:29

0001 // \class JetTracksAssociatorExplicit JetTracksAssociatorExplicit.cc
0002 //
0003 // Original Author:  Andrea Rizzi
0004 //         Created:  Wed Apr 12 11:12:49 CEST 2006
0005 // Accommodated for Jet Package by: Fedor Ratnikov Jul. 30, 2007
0006 //
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 
0009 // user include files
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "DataFormats/Common/interface/View.h"
0012 #include "DataFormats/JetReco/interface/Jet.h"
0013 #include "DataFormats/TrackReco/interface/Track.h"
0014 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0015 
0016 #include "JetTracksAssociatorExplicit.h"
0017 
0018 JetTracksAssociatorExplicit::JetTracksAssociatorExplicit(const edm::ParameterSet& fConfig) : mAssociatorExplicit() {
0019   mJets = consumes<edm::View<reco::Jet> >(fConfig.getParameter<edm::InputTag>("jets"));
0020   mTracks = consumes<reco::TrackCollection>(fConfig.getParameter<edm::InputTag>("tracks"));
0021 
0022   produces<reco::JetTracksAssociation::Container>();
0023 }
0024 
0025 JetTracksAssociatorExplicit::~JetTracksAssociatorExplicit() {}
0026 
0027 void JetTracksAssociatorExplicit::produce(edm::Event& fEvent, const edm::EventSetup& fSetup) {
0028   edm::Handle<edm::View<reco::Jet> > jets_h;
0029   fEvent.getByToken(mJets, jets_h);
0030   edm::Handle<reco::TrackCollection> tracks_h;
0031   fEvent.getByToken(mTracks, tracks_h);
0032 
0033   auto jetTracks = std::make_unique<reco::JetTracksAssociation::Container>(reco::JetRefBaseProd(jets_h));
0034 
0035   // format inputs
0036   std::vector<edm::RefToBase<reco::Jet> > allJets;
0037   allJets.reserve(jets_h->size());
0038   for (unsigned i = 0; i < jets_h->size(); ++i)
0039     allJets.push_back(jets_h->refAt(i));
0040   std::vector<reco::TrackRef> allTracks;
0041   allTracks.reserve(tracks_h->size());
0042   // run algo
0043   for (unsigned i = 0; i < tracks_h->size(); ++i) {
0044     allTracks.push_back(reco::TrackRef(tracks_h, i));
0045   }
0046 
0047   mAssociatorExplicit.produce(&*jetTracks, allJets, allTracks);
0048 
0049   // store output
0050   fEvent.put(std::move(jetTracks));
0051 }