CMS 3D CMS Logo

JetTracksAssociation.cc
Go to the documentation of this file.
1 #include "JetAssociationTemplate.icc"
3 
5 
8  return getValue(fContainer, fJet).size();
9 }
10 int reco::JetTracksAssociation::tracksNumber(const Container& fContainer, const reco::Jet& fJet) {
11  return getValue(fContainer, fJet).size();
12 }
15  const reco::JetBaseRef fJet) {
16  const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
17  math::XYZTLorentzVector result(0, 0, 0, 0);
18  for (unsigned t = 0; t < tracks->size(); ++t) {
19  const reco::Track& track = *((*tracks)[t]);
20  result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p()); // massless hypothesis
21  }
23 }
25  const reco::Jet& fJet) {
26  const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
27  math::XYZTLorentzVector result(0, 0, 0, 0);
28  for (unsigned t = 0; t < tracks->size(); ++t) {
29  const reco::Track& track = *((*tracks)[t]);
30  result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p()); // massless hypothesis
31  }
33 }
34 
36  const reco::JetBaseRef& fJet,
37  reco::TrackRefVector fValue) {
38  return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
39 }
40 
42  const reco::JetBaseRef& fJet,
43  reco::TrackRefVector fValue) {
44  return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
45 }
46 
48  const reco::JetBaseRef& fJet) {
49  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
50 }
51 
53  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
54 }
55 
56 std::vector<reco::JetBaseRef> reco::JetTracksAssociation::allJets(const Container& fContainer) {
57  return JetAssociationTemplate::allJets(fContainer);
58 }
59 
60 bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::JetBaseRef& fJet) {
61  return JetAssociationTemplate::hasJet(fContainer, fJet);
62 }
63 
64 bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::Jet& fJet) {
65  return JetAssociationTemplate::hasJet(fContainer, fJet);
66 }
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
Base class for all types of Jets.
Definition: Jet.h:20
math::PtEtaPhiELorentzVectorF LorentzVector
bool setValue(Container &, const reco::JetBaseRef &, const JetExtendedData &)
associate jet with value. Returns false and associate nothing if jet is already associated ...
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool setValue(Container &, const reco::JetBaseRef &, reco::TrackRefVector)
associate jet with value. Returns false and associate nothing if jet is already associated ...
bool hasJet(const Container &, const reco::JetBaseRef &)
check if jet is associated
const reco::TrackRefVector & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
auto const & tracks
cannot be loose
LorentzVector tracksP4(const Container &, const reco::JetBaseRef)
Get LorentzVector as sum of all tracks associated with jet.
int tracksNumber(const Container &, const reco::JetBaseRef)
Get number of tracks associated with jet.
const JetExtendedData & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
bool hasJet(const Container &, const reco::JetBaseRef &)
check if jet is associated