Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // \class JetTracksAssociatorAtVertex JetTracksAssociatorAtVertex.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 "JetTracksAssociatorAtVertex.h"
0017 
0018 JetTracksAssociatorAtVertex::JetTracksAssociatorAtVertex(const edm::ParameterSet& fConfig)
0019     : mAssociator(fConfig.getParameter<double>("coneSize")),
0020       mAssociatorAssigned(fConfig.getParameter<double>("coneSize")),
0021       useAssigned(false),
0022       pvSrc() {
0023   mJets = consumes<edm::View<reco::Jet> >(fConfig.getParameter<edm::InputTag>("jets"));
0024   mTracks = consumes<reco::TrackCollection>(fConfig.getParameter<edm::InputTag>("tracks"));
0025   if (fConfig.exists("useAssigned")) {
0026     useAssigned = fConfig.getParameter<bool>("useAssigned");
0027     pvSrc = consumes<reco::VertexCollection>(fConfig.getParameter<edm::InputTag>("pvSrc"));
0028   }
0029 
0030   produces<reco::JetTracksAssociation::Container>();
0031 }
0032 
0033 JetTracksAssociatorAtVertex::~JetTracksAssociatorAtVertex() {}
0034 
0035 void JetTracksAssociatorAtVertex::produce(edm::Event& fEvent, const edm::EventSetup& fSetup) {
0036   edm::Handle<edm::View<reco::Jet> > jets_h;
0037   fEvent.getByToken(mJets, jets_h);
0038   edm::Handle<reco::TrackCollection> tracks_h;
0039   fEvent.getByToken(mTracks, tracks_h);
0040 
0041   auto jetTracks = std::make_unique<reco::JetTracksAssociation::Container>(reco::JetRefBaseProd(jets_h));
0042 
0043   // format inputs
0044   std::vector<edm::RefToBase<reco::Jet> > allJets;
0045   allJets.reserve(jets_h->size());
0046   for (unsigned i = 0; i < jets_h->size(); ++i)
0047     allJets.push_back(jets_h->refAt(i));
0048   std::vector<reco::TrackRef> allTracks;
0049   allTracks.reserve(tracks_h->size());
0050   // run algo
0051   for (unsigned i = 0; i < tracks_h->size(); ++i) {
0052     allTracks.push_back(reco::TrackRef(tracks_h, i));
0053   }
0054   if (!useAssigned) {
0055     mAssociator.produce(&*jetTracks, allJets, allTracks);
0056   } else {
0057     edm::Handle<reco::VertexCollection> pvHandle;
0058     fEvent.getByToken(pvSrc, pvHandle);
0059     const reco::VertexCollection& vertices = *pvHandle.product();
0060     mAssociatorAssigned.produce(&*jetTracks, allJets, allTracks, vertices);
0061   }
0062 
0063   // store output
0064   fEvent.put(std::move(jetTracks));
0065 }