CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
JetTracksAssociatorAtVertex.cc
Go to the documentation of this file.
1 // \class JetTracksAssociatorAtVertex JetTracksAssociatorAtVertex.cc
2 //
3 // Original Author: Andrea Rizzi
4 // Created: Wed Apr 12 11:12:49 CEST 2006
5 // Accommodated for Jet Package by: Fedor Ratnikov Jul. 30, 2007
6 //
8 
9 // user include files
15 
17 
19  : mAssociator(fConfig.getParameter<double>("coneSize")),
20  mAssociatorAssigned(fConfig.getParameter<double>("coneSize")),
22  pvSrc() {
23  mJets = consumes<edm::View<reco::Jet> >(fConfig.getParameter<edm::InputTag>("jets"));
24  mTracks = consumes<reco::TrackCollection>(fConfig.getParameter<edm::InputTag>("tracks"));
25  if (fConfig.exists("useAssigned")) {
26  useAssigned = fConfig.getParameter<bool>("useAssigned");
27  pvSrc = consumes<reco::VertexCollection>(fConfig.getParameter<edm::InputTag>("pvSrc"));
28  }
29 
30  produces<reco::JetTracksAssociation::Container>();
31 }
32 
34 
37  fEvent.getByToken(mJets, jets_h);
39  fEvent.getByToken(mTracks, tracks_h);
40 
41  auto jetTracks = std::make_unique<reco::JetTracksAssociation::Container>(reco::JetRefBaseProd(jets_h));
42 
43  // format inputs
44  std::vector<edm::RefToBase<reco::Jet> > allJets;
45  allJets.reserve(jets_h->size());
46  for (unsigned i = 0; i < jets_h->size(); ++i)
47  allJets.push_back(jets_h->refAt(i));
48  std::vector<reco::TrackRef> allTracks;
49  allTracks.reserve(tracks_h->size());
50  // run algo
51  for (unsigned i = 0; i < tracks_h->size(); ++i) {
52  allTracks.push_back(reco::TrackRef(tracks_h, i));
53  }
54  if (!useAssigned) {
55  mAssociator.produce(&*jetTracks, allJets, allTracks);
56  } else {
58  fEvent.getByToken(pvSrc, pvHandle);
59  const reco::VertexCollection& vertices = *pvHandle.product();
60  mAssociatorAssigned.produce(&*jetTracks, allJets, allTracks, vertices);
61  }
62 
63  // store output
64  fEvent.put(std::move(jetTracks));
65 }
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
edm::EDGetTokenT< reco::TrackCollection > mTracks
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
JetTracksAssociatorAtVertex(const edm::ParameterSet &)
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks, const reco::VertexCollection &vertices) const
edm::RefToBaseProd< reco::Jet > JetRefBaseProd
Definition: JetCollection.h:13
def move
Definition: eostools.py:511
edm::EDGetTokenT< reco::VertexCollection > pvSrc
if true, use the track/jet association with vertex assignment to tracks
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks) const
JetTracksAssociationDRVertexAssigned mAssociatorAssigned
T const * product() const
Definition: Handle.h:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
JetTracksAssociationDRVertex mAssociator
edm::EDGetTokenT< edm::View< reco::Jet > > mJets
void produce(edm::Event &, const edm::EventSetup &) override