Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66
#include "JetAssociationTemplate.icc"
#include "DataFormats/TrackReco/interface/Track.h"

#include "DataFormats/JetReco/interface/JetTracksAssociation.h"

/// Get number of tracks associated with jet
int reco::JetTracksAssociation::tracksNumber(const Container& fContainer, const reco::JetBaseRef fJet) {
  return getValue(fContainer, fJet).size();
}
int reco::JetTracksAssociation::tracksNumber(const Container& fContainer, const reco::Jet& fJet) {
  return getValue(fContainer, fJet).size();
}
/// Get LorentzVector as sum of all tracks associated with jet.
reco::JetTracksAssociation::LorentzVector reco::JetTracksAssociation::tracksP4(const Container& fContainer,
                                                                               const reco::JetBaseRef fJet) {
  const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
  math::XYZTLorentzVector result(0, 0, 0, 0);
  for (unsigned t = 0; t < tracks->size(); ++t) {
    const reco::Track& track = *((*tracks)[t]);
    result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p());  // massless hypothesis
  }
  return reco::JetTracksAssociation::LorentzVector(result);
}
reco::JetTracksAssociation::LorentzVector reco::JetTracksAssociation::tracksP4(const Container& fContainer,
                                                                               const reco::Jet& fJet) {
  const reco::TrackRefVector* tracks = &getValue(fContainer, fJet);
  math::XYZTLorentzVector result(0, 0, 0, 0);
  for (unsigned t = 0; t < tracks->size(); ++t) {
    const reco::Track& track = *((*tracks)[t]);
    result += math::XYZTLorentzVector(track.px(), track.py(), track.pz(), track.p());  // massless hypothesis
  }
  return reco::JetTracksAssociation::LorentzVector(result);
}

bool reco::JetTracksAssociation::setValue(Container* fContainer,
                                          const reco::JetBaseRef& fJet,
                                          reco::TrackRefVector fValue) {
  return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
}

bool reco::JetTracksAssociation::setValue(Container& fContainer,
                                          const reco::JetBaseRef& fJet,
                                          reco::TrackRefVector fValue) {
  return JetAssociationTemplate::setValue(fContainer, fJet, fValue);
}

const reco::TrackRefVector& reco::JetTracksAssociation::getValue(const Container& fContainer,
                                                                 const reco::JetBaseRef& fJet) {
  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
}

const reco::TrackRefVector& reco::JetTracksAssociation::getValue(const Container& fContainer, const reco::Jet& fJet) {
  return JetAssociationTemplate::getValue<Container, reco::TrackRefVector>(fContainer, fJet);
}

std::vector<reco::JetBaseRef> reco::JetTracksAssociation::allJets(const Container& fContainer) {
  return JetAssociationTemplate::allJets(fContainer);
}

bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::JetBaseRef& fJet) {
  return JetAssociationTemplate::hasJet(fContainer, fJet);
}

bool reco::JetTracksAssociation::hasJet(const Container& fContainer, const reco::Jet& fJet) {
  return JetAssociationTemplate::hasJet(fContainer, fJet);
}