CMS 3D CMS Logo

Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Attributes

reco::Jet Class Reference

Base class for all types of Jets. More...

#include <Jet.h>

Inheritance diagram for reco::Jet:
reco::CompositePtrCandidate reco::LeafCandidate reco::Candidate pat::PATObject< reco::Jet > reco::BasicJet reco::CaloJet reco::GenJet reco::JPTJet reco::PFClusterJet reco::PFJet reco::TrackJet

List of all members.

Classes

class  EtaPhiMoments
 record to store eta-phi first and second moments More...

Public Types

typedef edm::Ptr< CandidateConstituent
typedef std::vector< ConstituentConstituents

Public Member Functions

float constituentEtaPhiSpread () const
float constituentPtDistribution () const
float etaetaMoment () const
 eta-eta second moment, ET weighted
float etaphiMoment () const
 eta-phi second moment, ET weighted
EtaPhiMoments etaPhiStatistics () const
 eta-phi statistics, ET weighted
float etInAnnulus (float fRmin, float fRmax) const
 ET in annulus between rmin and rmax around jet direction.
virtual Constituents getJetConstituents () const
 list of constituents
virtual std::vector< const
reco::Candidate * > 
getJetConstituentsQuick () const
 quick list of constituents
bool isJet () const
 Jet ()
 Default constructor.
 Jet (const LorentzVector &fP4, const Point &fVertex)
 Initiator.
 Jet (const LorentzVector &fP4, const Point &fVertex, const Constituents &fConstituents)
virtual float jetArea () const
 get jet area
float maxDistance () const
 maximum distance from jet to constituent
int nCarrying (float fFraction) const
 return # of constituent carrying fraction of energy
virtual int nConstituents () const
 # of constituents
virtual int nPasses () const
 number of passes taken by algorithm
float phiphiMoment () const
 phi-phi second moment, ET weighted
virtual float pileup () const
 pileup energy contribution as calculated by algorithm
virtual std::string print () const
 Print object.
virtual void scaleEnergy (double fScale)
 scale energy of the jet
virtual void setJetArea (float fArea)
 set jet area
virtual void setNPasses (int fPasses)
 Set number of passes taken by algorithm.
virtual void setPileup (float fEnergy)
 Set pileup energy contribution as calculated by algorithm.
virtual ~Jet ()
 Destructor.

Static Public Member Functions

static float detectorEta (float fZVertex, float fPhysicsEta)
 static function to convert physics eta to detector eta
static Candidate::LorentzVector detectorP4 (const Candidate::Point &vertex, const Candidate &inParticle)
static float physicsEta (float fZVertex, float fDetectorEta)
 static function to convert detector eta to physics eta
static Candidate::LorentzVector physicsP4 (const Candidate::Point &newVertex, const Candidate &inParticle, const Candidate::Point &oldVertex=Candidate::Point(0, 0, 0))

Private Attributes

float mJetArea
int mPassNumber
float mPileupEnergy

Detailed Description

Base class for all types of Jets.

GenericJet describes jets made from arbitrary constituents, No direct constituents reference is stored for now

Author:
Fedor Ratnikov, UMd
Version:
Mar 23, 2007 by F.R.
Id:
GenericJet.h,v 1.8 2007/08/24 17:35:23 fedor Exp

Jet describes properties common for all kinds of jets, essentially kinematics. Base class for all types of Jets.

Author:
Fedor Ratnikov, UMd
Version:
Original: April 22, 2005 by Fernando Varela Rodriguez.
May 23, 2006 by F.R.
Id:
Jet.h,v 1.32 2012/02/01 14:51:10 pandolf Exp

Definition at line 21 of file Jet.h.


Member Typedef Documentation

Definition at line 23 of file Jet.h.

typedef std::vector<Constituent> reco::Jet::Constituents

Definition at line 24 of file Jet.h.


Constructor & Destructor Documentation

reco::Jet::Jet ( ) [inline]

Default constructor.

Reimplemented in pat::Jet.

Definition at line 37 of file Jet.h.

: mJetArea (0), mPileupEnergy (0), mPassNumber (0) {}
Jet::Jet ( const LorentzVector fP4,
const Point fVertex 
)

Initiator.

Definition at line 175 of file Jet.cc.

  :  CompositePtrCandidate (0, fP4, fVertex),
     mJetArea (0),
     mPileupEnergy (0),
     mPassNumber (0)
{}
Jet::Jet ( const LorentzVector fP4,
const Point fVertex,
const Constituents fConstituents 
)

Definition at line 162 of file Jet.cc.

References reco::CompositePtrCandidate::addDaughter(), and i.

  :  CompositePtrCandidate (0, fP4, fVertex),
     mJetArea (0),
     mPileupEnergy (0),
     mPassNumber (0)
{
  for (unsigned i = 0; i < fConstituents.size (); i++) {
    addDaughter (fConstituents [i]);
  }
}  
virtual reco::Jet::~Jet ( ) [inline, virtual]

Destructor.

Reimplemented in pat::Jet.

Definition at line 42 of file Jet.h.

{}

Member Function Documentation

float Jet::constituentEtaPhiSpread ( ) const

Definition at line 395 of file Jet.cc.

References reco::deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, getJetConstituents(), and reco::LeafCandidate::pt().

Referenced by FFTGenericScaleCalculator::mapFFTJet().

                                         {

  Jet::Constituents constituents = this->getJetConstituents();


  float sum_pt2 = 0.;
  float sum_pt2deltaR2 = 0.;

  for( unsigned iConst=0; iConst<constituents.size(); ++iConst ) {

    LorentzVector thisConstituent = constituents[iConst]->p4();

    float pt = thisConstituent.Pt();
    float pt2 = pt*pt;
    double dR = deltaR (*this, *(constituents[iConst]));
    float pt2deltaR2 = pt*pt*dR*dR;

    sum_pt2 += pt2;
    sum_pt2deltaR2 += pt2deltaR2;

  } //for constituents

  float rmsCand_value = (sum_pt2>0.) ? sum_pt2deltaR2/sum_pt2 : 0.;

  return rmsCand_value;

} //constituentEtaPhiSpread
float Jet::constituentPtDistribution ( ) const

Definition at line 370 of file Jet.cc.

References getJetConstituents(), reco::LeafCandidate::pt(), and mathSSE::sqrt().

Referenced by FFTGenericScaleCalculator::mapFFTJet().

                                           {

  Jet::Constituents constituents = this->getJetConstituents();

  float sum_pt2 = 0.;
  float sum_pt  = 0.;

  for( unsigned iConst=0; iConst<constituents.size(); ++iConst ) {

    float pt = constituents[iConst]->p4().Pt();
    float pt2 = pt*pt;

    sum_pt += pt;
    sum_pt2 += pt2;

  } //for constituents

  float ptD_value = (sum_pt>0.) ? sqrt( sum_pt2 / (sum_pt*sum_pt) ) : 0.;

  return ptD_value;

} //constituentPtDistribution
float Jet::detectorEta ( float  fZVertex,
float  fPhysicsEta 
) [static]

static function to convert physics eta to detector eta

Definition at line 326 of file Jet.cc.

                                                         {
  CaloPointZ refPoint (fZVertex, fPhysicsEta);
  return refPoint.etaReference (0.);
}
Candidate::LorentzVector Jet::detectorP4 ( const Candidate::Point vertex,
const Candidate inParticle 
) [static]

Definition at line 340 of file Jet.cc.

References reco::Candidate::energy(), reco::Candidate::momentum(), np, reco::LeafCandidate::p(), and p3.

                                                                                                 {
  CaloPoint3D<Point> caloPoint(vertex,inParticle.momentum()); // Jet position in Calo.
  static const Point np(0,0,0);
  Vector detectorDir = caloPoint.caloPoint() - np;
  double p = inParticle.momentum().r();
  Vector p3 = p * detectorDir.unit();
  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
  return returnVector;
}
float Jet::etaetaMoment ( ) const

eta-eta second moment, ET weighted

eta-eta second moment

Definition at line 224 of file Jet.cc.

References getJetConstituentsQuick(), i, relativeConstraints::value, and CommonMethods::weight().

Referenced by fireworks::makeEveJetCone(), and JetIDSelectionFunctor::operator()().

                               {
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  double sumw = 0;
  double sum = 0;
  double sum2 = 0;
  int i = towers.size();
  while (--i >= 0) {
    double value = towers[i]->eta();
    double weight = towers[i]->et();
    sumw += weight;
    sum += value * weight;
    sum2 += value * value * weight;
  }
  return sumw > 0 ? (sum2 - sum*sum/sumw ) / sumw : 0;
}
float Jet::etaphiMoment ( ) const

eta-phi second moment, ET weighted

eta-phi second moment

Definition at line 259 of file Jet.cc.

References reco::deltaPhi(), getJetConstituentsQuick(), i, reco::LeafCandidate::phi(), and CommonMethods::weight().

                               {
  double phiRef = phi();
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  double sumw = 0;
  double sumA = 0;
  double sumB = 0;
  double sumAB = 0;
  int i = towers.size();
  while (--i >= 0) {
    double valueA = towers[i]->eta();
    double valueB = deltaPhi (towers[i]->phi(), phiRef);
    double weight = towers[i]->et();
    sumw += weight;
    sumA += valueA * weight;
    sumB += valueB * weight;
    sumAB += valueA * valueB * weight;
  }
  return sumw > 0 ? (sumAB - sumA*sumB/sumw ) / sumw : 0;
}
Jet::EtaPhiMoments Jet::etaPhiStatistics ( ) const

eta-phi statistics, ET weighted

eta-phi statistics

Definition at line 184 of file Jet.cc.

References reco::deltaPhi(), reco::LeafCandidate::eta(), reco::Jet::EtaPhiMoments::etaEtaMoment, reco::Jet::EtaPhiMoments::etaMean, reco::Jet::EtaPhiMoments::etaPhiMoment, getJetConstituentsQuick(), i, reco::LeafCandidate::phi(), reco::Jet::EtaPhiMoments::phiMean, reco::Jet::EtaPhiMoments::phiPhiMoment, query::result, and CommonMethods::weight().

                                              {
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  double phiRef = phi();
  double sumw = 0;
  double sumEta = 0;
  double sumEta2 = 0;
  double sumPhi = 0;
  double sumPhi2 = 0;
  double sumEtaPhi = 0;
  int i = towers.size();
  while (--i >= 0) {
    double eta = towers[i]->eta();
    double phi = deltaPhi (towers[i]->phi(), phiRef);
    double weight = towers[i]->et();
    sumw += weight;
    sumEta += eta * weight;
    sumEta2 += eta * eta * weight;
    sumPhi += phi * weight;
    sumPhi2 += phi * phi * weight;
    sumEtaPhi += eta * phi * weight;
  }
  Jet::EtaPhiMoments result;
  if (sumw > 0) {
    result.etaMean = sumEta / sumw;
    result.phiMean = deltaPhi (phiRef + sumPhi, 0.);
    result.etaEtaMoment = (sumEta2 - sumEta * sumEta / sumw) / sumw;
    result.phiPhiMoment = (sumPhi2 - sumPhi * sumPhi / sumw) / sumw;
    result.etaPhiMoment = (sumEtaPhi - sumEta * sumPhi / sumw) / sumw;
  }
  else {
    result.etaMean = 0;
    result.phiMean = 0;
    result.etaEtaMoment = 0;
    result.phiPhiMoment = 0;
    result.etaPhiMoment = 0;
  }
  return result;
}
float Jet::etInAnnulus ( float  fRmin,
float  fRmax 
) const

ET in annulus between rmin and rmax around jet direction.

et in annulus between rmin and rmax around jet direction

Definition at line 280 of file Jet.cc.

References reco::deltaR(), getJetConstituentsQuick(), i, and query::result.

                                                      {
  float result = 0;
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  int i = towers.size ();
  while (--i >= 0) {
    double r = deltaR (*this, *(towers[i]));
    if (r >= fRmin && r < fRmax) result += towers[i]->et ();
  }
  return result;
}
Jet::Constituents Jet::getJetConstituents ( ) const [virtual]
std::vector< const Candidate * > Jet::getJetConstituentsQuick ( ) const [virtual]

quick list of constituents

Definition at line 359 of file Jet.cc.

References reco::CompositePtrCandidate::daughter(), i, reco::CompositePtrCandidate::numberOfDaughters(), and query::result.

Referenced by FWPFTauProxyBuilder::buildViewType(), etaetaMoment(), etaphiMoment(), etaPhiStatistics(), etInAnnulus(), maxDistance(), nCarrying(), and phiphiMoment().

                                                                {
  std::vector<const Candidate*> result;
  int nDaughters = numberOfDaughters();
  for (int i = 0; i < nDaughters; ++i) { 
    result.push_back (daughter (i));
  }
  return result;
}
bool Jet::isJet ( ) const [virtual]

Reimplemented from reco::LeafCandidate.

Definition at line 449 of file Jet.cc.

                      {
  return true;
}
virtual float reco::Jet::jetArea ( ) const [inline, virtual]

get jet area

Definition at line 106 of file Jet.h.

References mJetArea.

Referenced by L1FastjetCorrector::correction(), PileupJetIdProducer::produce(), and HiL1Subtractor::produce().

{return mJetArea;}
float Jet::maxDistance ( ) const

maximum distance from jet to constituent

Definition at line 307 of file Jet.cc.

References reco::deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, getJetConstituentsQuick(), i, and query::result.

                              {
  float result = 0;
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  for (unsigned i = 0; i < towers.size(); ++i) {
    float dR = deltaR (*(towers[i]), *this);
    if (dR > result)  result = dR;
  }
  return result;
}
int Jet::nCarrying ( float  fFraction) const

return # of constituent carrying fraction of energy

return # of constituent carring fraction of energy. Assume ordered towers

Definition at line 292 of file Jet.cc.

References reco::LeafCandidate::et(), getJetConstituentsQuick(), i, and query::result.

Referenced by pat::Jet::n60(), reco::CaloJet::n60(), pat::Jet::n90(), and reco::CaloJet::n90().

                                         {
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  if (fFraction >= 1) return towers.size();
  double totalEt = 0;
  for (unsigned i = 0; i < towers.size(); ++i) totalEt += towers[i]->et();
  double fractionEnergy = totalEt * fFraction;
  unsigned result = 0;
  for (; result < towers.size(); ++result) {
    fractionEnergy -= towers[result]->et();
    if (fractionEnergy <= 0) return result+1;
  }
  return 0;
}
virtual int reco::Jet::nConstituents ( ) const [inline, virtual]
virtual int reco::Jet::nPasses ( ) const [inline, virtual]

number of passes taken by algorithm

Definition at line 116 of file Jet.h.

References mPassNumber.

{return mPassNumber;}
float Jet::phiphiMoment ( ) const

phi-phi second moment, ET weighted

phi-phi second moment

Definition at line 241 of file Jet.cc.

References reco::deltaPhi(), getJetConstituentsQuick(), i, reco::LeafCandidate::phi(), relativeConstraints::value, and CommonMethods::weight().

Referenced by fireworks::makeEveJetCone(), and JetIDSelectionFunctor::operator()().

                               {
  double phiRef = phi();
  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
  double sumw = 0;
  double sum = 0;
  double sum2 = 0;
  int i = towers.size();
  while (--i >= 0) {
    double value = deltaPhi (towers[i]->phi(), phiRef);
    double weight = towers[i]->et();
    sumw += weight;
    sum += value * weight;
    sum2 += value * value * weight;
  }
  return sumw > 0 ? (sum2 - sum*sum/sumw ) / sumw : 0;
}
float Jet::physicsEta ( float  fZVertex,
float  fDetectorEta 
) [static]

static function to convert detector eta to physics eta

Definition at line 319 of file Jet.cc.

Referenced by JetTracksAssociationDRVertexAssigned::produce().

                                                         {
  CaloPointZ refPoint (0., fDetectorEta);
  return refPoint.etaReference (fZVertex);
}
Candidate::LorentzVector Jet::physicsP4 ( const Candidate::Point newVertex,
const Candidate inParticle,
const Candidate::Point oldVertex = Candidate::Point(0,0,0) 
) [static]

Definition at line 331 of file Jet.cc.

References reco::Candidate::energy(), reco::Candidate::momentum(), reco::LeafCandidate::p(), and p3.

                                                                                                                                   {
  CaloPoint3D<Point> caloPoint(oldVertex,inParticle.momentum()); // Jet position in Calo.
  Vector physicsDir = caloPoint.caloPoint() - newVertex;
  double p = inParticle.momentum().r();
  Vector p3 = p * physicsDir.unit();
  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
  return returnVector;
}
virtual float reco::Jet::pileup ( ) const [inline, virtual]

pileup energy contribution as calculated by algorithm

Definition at line 111 of file Jet.h.

References mPileupEnergy.

{return mPileupEnergy;}
std::string Jet::print ( void  ) const [virtual]

Print object.

Reimplemented in reco::BasicJet, reco::CaloJet, reco::GenJet, reco::JPTJet, reco::PFClusterJet, reco::PFJet, and reco::TrackJet.

Definition at line 425 of file Jet.cc.

References reco::CompositePtrCandidate::daughterPtr(), reco::LeafCandidate::eta(), edm::Ptr< T >::id(), edm::Ptr< T >::isNonnull(), edm::Ptr< T >::key(), nConstituents(), reco::CompositePtrCandidate::numberOfDaughters(), dbtoconf::out, reco::LeafCandidate::p(), reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::LeafCandidate::pz().

Referenced by reco::TrackJet::print(), and reco::PFClusterJet::print().

                            {
  std::ostringstream out;
  out << "Jet p/px/py/pz/pt: " << p() << '/' << px () << '/' << py() << '/' << pz() << '/' << pt() << std::endl
      << "    eta/phi: " << eta () << '/' << phi () << std::endl
       << "    # of constituents: " << nConstituents () << std::endl;
  out << "    Constituents:" << std::endl;
  for (unsigned index = 0; index < numberOfDaughters(); index++) {
    Constituent constituent = daughterPtr (index); // deref
    if (constituent.isNonnull()) {
      out << "      #" << index << " p/pt/eta/phi: " 
          << constituent->p() << '/' << constituent->pt() << '/' << constituent->eta() << '/' << constituent->phi() 
          << "    productId/index: " << constituent.id() << '/' << constituent.key() << std::endl; 
    }
    else {
      out << "      #" << index << " constituent is not available in the event"  << std::endl;
    }
  }
  return out.str ();
}
void Jet::scaleEnergy ( double  fScale) [virtual]
virtual void reco::Jet::setJetArea ( float  fArea) [inline, virtual]
virtual void reco::Jet::setNPasses ( int  fPasses) [inline, virtual]

Set number of passes taken by algorithm.

Definition at line 114 of file Jet.h.

References mPassNumber.

{mPassNumber = fPasses;}
virtual void reco::Jet::setPileup ( float  fEnergy) [inline, virtual]

Set pileup energy contribution as calculated by algorithm.

Definition at line 109 of file Jet.h.

References mPileupEnergy.

Referenced by HiL1Subtractor::produce(), FastjetJetProducer::produceTrackJets(), and cms::SubEventGenJetProducer::runAlgorithm().

{mPileupEnergy = fEnergy;}

Member Data Documentation

float reco::Jet::mJetArea [private]

Definition at line 121 of file Jet.h.

Referenced by jetArea(), and setJetArea().

int reco::Jet::mPassNumber [private]

Definition at line 123 of file Jet.h.

Referenced by nPasses(), and setNPasses().

float reco::Jet::mPileupEnergy [private]

Definition at line 122 of file Jet.h.

Referenced by pileup(), and setPileup().