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 //
8 
9 // user include files
15 
17 
19  : mAssociator (fConfig.getParameter<double> ("coneSize")),
20  mAssociatorAssigned (fConfig.getParameter<double> ("coneSize")),
22 {
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  std::auto_ptr<reco::JetTracksAssociation::Container> jetTracks (new 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) allJets.push_back (jets_h->refAt(i));
47  std::vector <reco::TrackRef> allTracks;
48  allTracks.reserve (tracks_h->size());
49  // run algo
50  for (unsigned i = 0; i < tracks_h->size(); ++i) {
51  allTracks.push_back (reco::TrackRef (tracks_h, i));
52  }
53  if ( !useAssigned ) {
54  mAssociator.produce (&*jetTracks, allJets, allTracks);
55  } else {
57  fEvent.getByToken(pvSrc,pvHandle);
58  const reco::VertexCollection & vertices = *pvHandle.product();
59  mAssociatorAssigned.produce (&*jetTracks, allJets, allTracks,vertices);
60  }
61 
62 
63  // store output
64  fEvent.put (jetTracks);
65 }
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
edm::EDGetTokenT< reco::TrackCollection > mTracks
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
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 &)
edm::EDGetTokenT< edm::View< reco::Jet > > mJets
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:113
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:81
JetTracksAssociationDRVertex mAssociator
volatile std::atomic< bool > shutdown_flag false