Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:25

0001 #include "JetAssociationTemplate.icc"
0002 #include "DataFormats/TrackReco/interface/Track.h"
0003 
0004 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0005 
0006 /// Get number of tracks associated with jet
0007 int reco::JetTracksAssociation::tracksNumber(const Container& fContainer, const reco::JetBaseRef fJet) {
0008   return getValue(fContainer, fJet).size();
0009 }
0010 int reco::JetTracksAssociation::tracksNumber(const Container& fContainer, const reco::Jet& fJet) {
0011   return getValue(fContainer, fJet).size();
0012 }
0013 /// Get LorentzVector as sum of all tracks associated with jet.
0014 reco::JetTracksAssociation::LorentzVector reco::JetTracksAssociation::tracksP4(const Container& fContainer,
0015                                                                                const reco::JetBaseRef fJet) {
0016   const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
0017   math::XYZTLorentzVector result(0, 0, 0, 0);
0018   for (unsigned t = 0; t < tracks->size(); ++t) {
0019     const reco::Track& track = *((*tracks)[t]);
0020     result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p());  // massless hypothesis
0021   }
0022   return reco::JetTracksAssociation::LorentzVector(result);
0023 }
0024 reco::JetTracksAssociation::LorentzVector reco::JetTracksAssociation::tracksP4(const Container& fContainer,
0025                                                                                const reco::Jet& fJet) {
0026   const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
0027   math::XYZTLorentzVector result(0, 0, 0, 0);
0028   for (unsigned t = 0; t < tracks->size(); ++t) {
0029     const reco::Track& track = *((*tracks)[t]);
0030     result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p());  // massless hypothesis
0031   }
0032   return reco::JetTracksAssociation::LorentzVector(result);
0033 }
0034 
0035 bool reco::JetTracksAssociation::setValue(Container* fContainer,
0036                                           const reco::JetBaseRef& fJet,
0037                                           reco::TrackRefVector fValue) {
0038   return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
0039 }
0040 
0041 bool reco::JetTracksAssociation::setValue(Container& fContainer,
0042                                           const reco::JetBaseRef& fJet,
0043                                           reco::TrackRefVector fValue) {
0044   return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
0045 }
0046 
0047 const reco::TrackRefVector& reco::JetTracksAssociation::getValue(const Container& fContainer,
0048                                                                  const reco::JetBaseRef& fJet) {
0049   return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
0050 }
0051 
0052 const reco::TrackRefVector& reco::JetTracksAssociation::getValue(const Container& fContainer, const reco::Jet& fJet) {
0053   return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
0054 }
0055 
0056 std::vector<reco::JetBaseRef> reco::JetTracksAssociation::allJets(const Container& fContainer) {
0057   return JetAssociationTemplate::allJets(fContainer);
0058 }
0059 
0060 bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::JetBaseRef& fJet) {
0061   return JetAssociationTemplate::hasJet(fContainer, fJet);
0062 }
0063 
0064 bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::Jet& fJet) {
0065   return JetAssociationTemplate::hasJet(fContainer, fJet);
0066 }