CMS 3D CMS Logo

Typedefs | Functions

reco::JetTracksAssociation Namespace Reference

Typedefs

typedef edm::AssociationVector
< reco::JetRefBaseProd, Values
Container
typedef
math::PtEtaPhiELorentzVectorF 
LorentzVector
typedef edm::Ref< ContainerRef
typedef edm::RefProd< ContainerRefProd
typedef edm::RefVector< ContainerRefVector
typedef
Container::transient_vector_type 
transient_vector_type
typedef reco::TrackRefVector Value
typedef Container::value_type value_type
typedef std::vector< ValueValues

Functions

std::vector< reco::JetBaseRefallJets (const Container &)
 fill list of all jets associated with values. Return # of jets in the list
const reco::TrackRefVectorgetValue (const Container &, const reco::JetBaseRef &)
 get value for the association. Throw exception if no association found
const reco::TrackRefVectorgetValue (const Container &, const reco::Jet &)
 get value for the association. Throw exception if no association found
bool hasJet (const Container &, const reco::Jet &)
 check if jet is associated
bool hasJet (const Container &, const reco::JetBaseRef &)
 check if jet is associated
bool setValue (Container *, const reco::JetBaseRef &, reco::TrackRefVector)
 associate jet with value. Returns false and associate nothing if jet is already associated
bool setValue (Container &, const reco::JetBaseRef &, reco::TrackRefVector)
 associate jet with value. Returns false and associate nothing if jet is already associated
int tracksNumber (const Container &, const reco::JetBaseRef)
 Get number of tracks associated with jet.
int tracksNumber (const Container &, const reco::Jet &)
 Get number of tracks associated with jet.
LorentzVector tracksP4 (const Container &, const reco::JetBaseRef)
 Get LorentzVector as sum of all tracks associated with jet.
LorentzVector tracksP4 (const Container &, const reco::Jet &)
 Get LorentzVector as sum of all tracks associated with jet.

Typedef Documentation

Definition at line 29 of file JetTracksAssociation.h.

Definition at line 26 of file JetTracksAssociation.h.

Definition at line 32 of file JetTracksAssociation.h.

Definition at line 33 of file JetTracksAssociation.h.

Definition at line 34 of file JetTracksAssociation.h.

Definition at line 31 of file JetTracksAssociation.h.

Definition at line 27 of file JetTracksAssociation.h.

Definition at line 30 of file JetTracksAssociation.h.

Definition at line 28 of file JetTracksAssociation.h.


Function Documentation

std::vector< reco::JetBaseRef > reco::JetTracksAssociation::allJets ( const Container &  fContainer)

fill list of all jets associated with values. Return # of jets in the list

Definition at line 58 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::allJets().

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

get value for the association. Throw exception if no association found

Definition at line 48 of file JetTracksAssociation.cc.

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

get value for the association. Throw exception if no association found

Definition at line 53 of file JetTracksAssociation.cc.

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

check if jet is associated

Definition at line 67 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::hasJet().

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

check if jet is associated

Definition at line 62 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::hasJet().

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

associate jet with value. Returns false and associate nothing if jet is already associated

Definition at line 36 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::setValue().

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

associate jet with value. Returns false and associate nothing if jet is already associated

Definition at line 42 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::setValue().

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

Get number of tracks associated with jet.

Definition at line 7 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::getValue().

Referenced by JetValidation::analyze(), and JetExtender::produce().

                                                                                                  {
  return getValue (fContainer, fJet).size();
}
int reco::JetTracksAssociation::tracksNumber ( const Container &  fContainer,
const reco::Jet fJet 
)

Get number of tracks associated with jet.

Definition at line 10 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::getValue().

                                                                                            {
  return getValue (fContainer, fJet).size();
}
reco::JetTracksAssociation::LorentzVector reco::JetTracksAssociation::tracksP4 ( const Container &  fContainer,
const reco::JetBaseRef  fJet 
)

Get LorentzVector as sum of all tracks associated with jet.

Definition at line 15 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::getValue(), reco::TrackBase::p(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), query::result, edm::RefVector< C, T, F >::size(), lumiQTWidget::t, and testEve_cfg::tracks.

Referenced by JetValidation::analyze(), and JetExtender::produce().

                                                                                          {
  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 
)

Get LorentzVector as sum of all tracks associated with jet.

Definition at line 25 of file JetTracksAssociation.cc.

References reco::JetExtendedAssociation::getValue(), reco::TrackBase::p(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), query::result, edm::RefVector< C, T, F >::size(), lumiQTWidget::t, and testEve_cfg::tracks.

                                                                                    {
  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);
}