CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 }
16  const reco::TrackRefVector* tracks = &getValue (fContainer, fJet);
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 reco::JetTracksAssociation::tracksP4 (const Container& fContainer, const reco::Jet& fJet) {
26  const reco::TrackRefVector* tracks = &getValue (fContainer, fJet);
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 
35 
37  const reco::JetBaseRef& fJet,
38  reco::TrackRefVector fValue) {
39  return JetAssociationTemplate::setValue (fContainer, fJet,fValue);
40 }
41 
43  const reco::JetBaseRef& fJet,
44  reco::TrackRefVector fValue) {
45  return JetAssociationTemplate::setValue (fContainer, fJet,fValue);
46 }
47 
49  const reco::JetBaseRef& fJet) {
50  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector> (fContainer, fJet);
51 }
52 
54  const reco::Jet& fJet) {
55  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector> (fContainer, fJet);
56 }
57 
58 std::vector<reco::JetBaseRef > reco::JetTracksAssociation::allJets (const Container& fContainer) {
59  return JetAssociationTemplate::allJets (fContainer);
60 }
61 
63  const reco::JetBaseRef& fJet) {
64  return JetAssociationTemplate::hasJet (fContainer, fJet);
65 }
66 
68  const reco::Jet& fJet) {
69  return JetAssociationTemplate::hasJet (fContainer, fJet);
70 }
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
double p() const
momentum vector magnitude
Definition: TrackBase.h:602
Base class for all types of Jets.
Definition: Jet.h:20
math::PtEtaPhiELorentzVectorF LorentzVector
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:614
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
tuple result
Definition: query.py:137
const reco::TrackRefVector & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:626
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
tuple tracks
Definition: testEve_cfg.py:39
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.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
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
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:620