CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 // $Id: JetTracksAssociatorAtVertex.cc,v 1.8 2011/11/11 19:44:39 srappocc Exp $
7 //
9 
10 // user include files
16 
18 
20  : mJets (fConfig.getParameter<edm::InputTag> ("jets")),
21  mTracks (fConfig.getParameter<edm::InputTag> ("tracks")),
22  mAssociator (fConfig.getParameter<double> ("coneSize")),
23  mAssociatorAssigned (fConfig.getParameter<double> ("coneSize")),
24  useAssigned(false), pvSrc()
25 {
26 
27  if ( fConfig.exists("useAssigned") ) {
28  useAssigned = fConfig.getParameter<bool> ("useAssigned");
29  pvSrc = fConfig.getParameter<edm::InputTag> ("pvSrc");
30  }
31 
32  produces<reco::JetTracksAssociation::Container> ();
33 }
34 
36 
39  fEvent.getByLabel (mJets, jets_h);
41  fEvent.getByLabel (mTracks, tracks_h);
42 
43  std::auto_ptr<reco::JetTracksAssociation::Container> jetTracks (new reco::JetTracksAssociation::Container (reco::JetRefBaseProd(jets_h)));
44 
45  // format inputs
46  std::vector <edm::RefToBase<reco::Jet> > allJets;
47  allJets.reserve (jets_h->size());
48  for (unsigned i = 0; i < jets_h->size(); ++i) allJets.push_back (jets_h->refAt(i));
49  std::vector <reco::TrackRef> allTracks;
50  allTracks.reserve (tracks_h->size());
51  // run algo
52  for (unsigned i = 0; i < tracks_h->size(); ++i) {
53  allTracks.push_back (reco::TrackRef (tracks_h, i));
54  }
55  if ( !useAssigned ) {
56  mAssociator.produce (&*jetTracks, allJets, allTracks);
57  } else {
59  fEvent.getByLabel(pvSrc,pvHandle);
60  const reco::VertexCollection & vertices = *pvHandle.product();
61  mAssociatorAssigned.produce (&*jetTracks, allJets, allTracks,vertices);
62  }
63 
64 
65  // store output
66  fEvent.put (jetTracks);
67 }
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
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 &)
virtual void produce(edm::Event &, const edm::EventSetup &)
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
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks) const
edm::InputTag pvSrc
if true, use the track/jet association with vertex assignment to tracks
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
JetTracksAssociationDRVertexAssigned mAssociatorAssigned
JetTracksAssociationDRVertex mAssociator
T const * product() const
Definition: Handle.h:74