A minimal class representing a jet of particles. More...
#include <Jet.hh>
Public Types | |
typedef vector< FourMomentum > ::iterator | iterator |
Define a Jet::iterator via a typedef. | |
typedef vector< FourMomentum > ::const_iterator | const_iterator |
Define a Jet::const_iterator via a typedef. | |
Public Member Functions | |
Jet () | |
Constructor. | |
iterator | begin () |
Get a begin iterator over the particle/track four-momenta in this jet. | |
iterator | end () |
Get an end iterator over the particle/track four-momenta in this jet. | |
const_iterator | begin () const |
Get a const begin iterator over the particle/track four-momenta in this jet. | |
const_iterator | end () const |
Get a const end iterator over the particle/track four-momenta in this jet. | |
vector< FourMomentum > & | momenta () |
Get the track momenta in this jet. | |
const vector< FourMomentum > & | momenta () const |
Get the track momenta in this jet (const version). | |
vector< Particle > & | particles () |
Get the Rivet::Particles (full information) in this jet. | |
const vector< Particle > & | particles () const |
Get the Rivet::Particles (full information) in this jet (const version). | |
size_t | size () const |
Number of particles (tracks) in this jet. | |
Jet & | setParticles (const vector< FourMomentum > &particles) |
Set the particles/tracks collection. | |
Jet & | addParticle (const FourMomentum &particle) |
Add a particle/track to this jet. | |
Jet & | addParticle (const Particle &particle) |
Add a particle/track to this jet. | |
bool | containsParticle (const Particle &particle) const |
Check whether this jet contains a particular particle. | |
bool | containsParticleId (PdgId pid) const |
Check whether this jet contains a certain particle type. | |
bool | containsParticleId (const vector< PdgId > &pids) const |
Check whether this jet contains at least one of certain particle types. | |
bool | containsCharm () const |
Check whether this jet contains a charm-flavoured hadron (or decay products from one). | |
bool | containsBottom () const |
Check whether this jet contains a bottom-flavoured hadron (or decay products from one). | |
Jet & | clear () |
Reset this jet as empty. | |
double | ptWeightedEta () const |
double | ptWeightedPhi () const |
double | eta () const |
Get the unweighted average ![]() | |
double | phi () const |
Get the unweighted average ![]() | |
const FourMomentum & | momentum () const |
Get equivalent single momentum four-vector. (caches). | |
double | totalEnergy () const |
Get the total energy of this jet. | |
double | neutralEnergy () const |
Get the energy carried in this jet by neutral particles. | |
double | hadronicEnergy () const |
Get the energy carried in this jet by hadrons. | |
double | ptSum () const |
Get the sum of the ![]() | |
double | EtSum () const |
Get the sum of the ![]() | |
Private Member Functions | |
void | _resetCaches () const |
Clear the internal cached values. Const because cache variables are mutable. | |
void | _calcMomVector () const |
Calculate cached equivalent momentum vector. Const because cache variables are mutable. | |
void | _calcPtAvgs () const |
Private Attributes | |
std::vector< FourMomentum > | _particles |
The particle tracks. | |
ParticleVector | _fullParticles |
Full particle information including tracks, ID etc. | |
double | _ptWeightedPhi |
Cached values of the ![]() ![]() ![]() | |
double | _ptWeightedEta |
bool | _okPtWeightedPhi |
bool | _okPtWeightedEta |
FourMomentum | _momentum |
Cached effective jet 4-vector. | |
bool | _okMomentum |
A minimal class representing a jet of particles.
Definition at line 12 of file Jet.hh.
typedef vector<FourMomentum>::const_iterator const_iterator |
Define a Jet::const_iterator via a typedef.
typedef vector<FourMomentum>::iterator iterator |
Define a Jet::iterator via a typedef.
Jet | ( | ) |
Constructor.
Definition at line 10 of file Jet.cc.
References Jet::clear().
00011 : ParticleBase() 00012 { 00013 clear(); 00014 }
void _calcMomVector | ( | ) | const [private] |
Calculate cached equivalent momentum vector. Const because cache variables are mutable.
Definition at line 194 of file Jet.cc.
References Jet::_momentum, Jet::_okMomentum, Jet::begin(), and Jet::end().
Referenced by Jet::momentum().
00194 { 00195 if (!_okMomentum) { 00196 _momentum = accumulate(begin(), end(), FourMomentum()); 00197 _okMomentum = true; 00198 } 00199 }
void _calcPtAvgs | ( | ) | const [private] |
Internal caching method to calculate the average and
for this jet, weighted by the
values of the constituent tracks. Const because cache variables are mutable.
Definition at line 202 of file Jet.cc.
References Jet::_okPtWeightedEta, Jet::_okPtWeightedPhi, Jet::_ptWeightedEta, Jet::_ptWeightedPhi, FourVector::azimuthalAngle(), Rivet::mapAngleMPiToPi(), Jet::momenta(), Jet::phi(), FourVector::pseudorapidity(), and FourMomentum::pT().
Referenced by Jet::ptWeightedEta(), and Jet::ptWeightedPhi().
00202 { 00203 if (!_okPtWeightedEta || !_okPtWeightedPhi) { 00204 double ptwetasum(0.0), ptwdphisum(0.0), ptsum(0.0); 00205 double phi0 = phi(); 00206 foreach (const FourMomentum& p, momenta()) { 00207 double pt = p.pT(); 00208 ptsum += pt; 00209 ptwetasum += pt * p.pseudorapidity(); 00210 ptwdphisum += pt * mapAngleMPiToPi(phi0 - p.azimuthalAngle()); 00211 } 00212 _ptWeightedEta = ptwetasum/ptsum; 00213 _okPtWeightedEta = true; 00214 _ptWeightedPhi = mapAngleMPiToPi(phi0 + ptwdphisum/ptsum); 00215 _okPtWeightedPhi = true; 00216 } 00217 }
void _resetCaches | ( | ) | const [private] |
Clear the internal cached values. Const because cache variables are mutable.
Definition at line 187 of file Jet.cc.
References Jet::_okMomentum, Jet::_okPtWeightedEta, and Jet::_okPtWeightedPhi.
Referenced by Jet::addParticle(), Jet::clear(), and Jet::setParticles().
00187 { 00188 _okPtWeightedPhi = false; 00189 _okPtWeightedEta = false; 00190 _okMomentum = false; 00191 }
Add a particle/track to this jet.
Definition at line 31 of file Jet.cc.
References Jet::_fullParticles, Jet::_particles, Jet::_resetCaches(), and Particle::momentum().
00031 { 00032 _fullParticles.push_back(particle); 00033 _particles.push_back(particle.momentum()); 00034 _resetCaches(); 00035 return *this; 00036 }
Jet & addParticle | ( | const FourMomentum & | particle | ) |
Add a particle/track to this jet.
Definition at line 24 of file Jet.cc.
References Jet::_particles, and Jet::_resetCaches().
Referenced by FastJets::_pseudojetsToJets().
00024 { 00025 _particles.push_back(particle); 00026 _resetCaches(); 00027 return *this; 00028 }
const_iterator begin | ( | ) | const [inline] |
Get a const begin iterator over the particle/track four-momenta in this jet.
Definition at line 35 of file Jet.hh.
References Jet::_particles.
00035 { 00036 return _particles.begin(); 00037 }
iterator begin | ( | ) | [inline] |
Get a begin iterator over the particle/track four-momenta in this jet.
Definition at line 25 of file Jet.hh.
References Jet::_particles.
Referenced by Jet::_calcMomVector().
00025 { 00026 return _particles.begin(); 00027 }
Jet & clear | ( | ) |
Reset this jet as empty.
Definition at line 132 of file Jet.cc.
References Jet::_fullParticles, Jet::_particles, and Jet::_resetCaches().
Referenced by Jet::Jet().
00132 { 00133 _particles.clear(); 00134 _fullParticles.clear(); 00135 _resetCaches(); 00136 return *this; 00137 }
bool containsBottom | ( | ) | const |
Check whether this jet contains a bottom-flavoured hadron (or decay products from one).
Definition at line 115 of file Jet.cc.
References Rivet::BQUARK, Particle::genParticle(), Rivet::PID::hasBottom(), Rivet::PID::isHadron(), Rivet::particles(), Jet::particles(), Particle::pdgId(), and Rivet::pi.
Referenced by MC_TTBAR::analyze(), ExampleAnalysis::analyze(), and CDF_2008_S7782535::analyze().
00115 { 00116 foreach (const Particle& p, particles()) { 00117 const PdgId pid = p.pdgId(); 00118 if (abs(pid) == BQUARK) return true; 00119 if (PID::isHadron(pid) && PID::hasBottom(pid)) return true; 00120 HepMC::GenVertex* gv = p.genParticle().production_vertex(); 00121 if (gv) { 00122 foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { 00123 const PdgId pid2 = pi->pdg_id(); 00124 if (PID::isHadron(pid2) && PID::hasBottom(pid2)) return true; 00125 } 00126 } 00127 } 00128 return false; 00129 }
bool containsCharm | ( | ) | const |
Check whether this jet contains a charm-flavoured hadron (or decay products from one).
Definition at line 98 of file Jet.cc.
References Rivet::CQUARK, Particle::genParticle(), Rivet::PID::hasCharm(), Rivet::PID::isHadron(), Rivet::particles(), Jet::particles(), Particle::pdgId(), and Rivet::pi.
00098 { 00099 foreach (const Particle& p, particles()) { 00100 const PdgId pid = p.pdgId(); 00101 if (abs(pid) == CQUARK) return true; 00102 if (PID::isHadron(pid) && PID::hasCharm(pid)) return true; 00103 HepMC::GenVertex* gv = p.genParticle().production_vertex(); 00104 if (gv) { 00105 foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { 00106 const PdgId pid2 = pi->pdg_id(); 00107 if (PID::isHadron(pid2) && PID::hasCharm(pid2)) return true; 00108 } 00109 } 00110 } 00111 return false; 00112 }
bool containsParticle | ( | const Particle & | particle | ) | const |
Check whether this jet contains a particular particle.
Definition at line 39 of file Jet.cc.
References Particle::genParticle(), and Jet::particles().
00039 { 00040 const int barcode = particle.genParticle().barcode(); 00041 foreach (const Particle& p, particles()) { 00042 if (p.genParticle().barcode() == barcode) return true; 00043 } 00044 return false; 00045 }
bool containsParticleId | ( | const vector< PdgId > & | pids | ) | const |
Check whether this jet contains at least one of certain particle types.
Definition at line 56 of file Jet.cc.
References Jet::particles(), and Particle::pdgId().
bool containsParticleId | ( | PdgId | pid | ) | const |
Check whether this jet contains a certain particle type.
Definition at line 48 of file Jet.cc.
References Jet::particles(), and Particle::pdgId().
00048 { 00049 foreach (const Particle& p, particles()) { 00050 if (p.pdgId() == pid) return true; 00051 } 00052 return false; 00053 }
const_iterator end | ( | ) | const [inline] |
Get a const end iterator over the particle/track four-momenta in this jet.
Definition at line 40 of file Jet.hh.
References Jet::_particles.
00040 { 00041 return _particles.end(); 00042 }
iterator end | ( | ) | [inline] |
Get an end iterator over the particle/track four-momenta in this jet.
Definition at line 30 of file Jet.hh.
References Jet::_particles.
Referenced by Jet::_calcMomVector().
00030 { 00031 return _particles.end(); 00032 }
double eta | ( | ) | const |
Get the unweighted average for this jet. (caches).
Definition at line 154 of file Jet.cc.
References FourVector::eta(), and Jet::momentum().
Referenced by STAR_2006_S6870392::analyze(), and ATLAS_2010_S8919674::analyze().
00154 { 00155 return momentum().eta(); 00156 00157 }
double EtSum | ( | ) | const |
Get the sum of the values of the constituent tracks. (caches).
Definition at line 182 of file Jet.cc.
References FourMomentum::Et(), and Jet::momentum().
Referenced by CDF_2004_S5839831::analyze(), Rivet::cmpJetsByEt(), and Rivet::cmpJetsByEtDesc().
00182 { 00183 return momentum().Et(); 00184 }
double hadronicEnergy | ( | ) | const |
Get the energy carried in this jet by hadrons.
Definition at line 86 of file Jet.cc.
References FourMomentum::E(), Rivet::PID::isHadron(), Particle::momentum(), Jet::particles(), and Particle::pdgId().
00086 { 00087 double e_hadr = 0.0; 00088 foreach (const Particle& p, particles()) { 00089 const PdgId pid = p.pdgId(); 00090 if (PID::isHadron(pid)) { 00091 e_hadr += p.momentum().E(); 00092 } 00093 } 00094 return e_hadr; 00095 }
const vector<FourMomentum>& momenta | ( | ) | const [inline] |
Get the track momenta in this jet (const version).
Definition at line 50 of file Jet.hh.
References Jet::_particles.
00050 { 00051 return _particles; 00052 }
vector<FourMomentum>& momenta | ( | ) | [inline] |
Get the track momenta in this jet.
Definition at line 45 of file Jet.hh.
References Jet::_particles.
Referenced by Jet::_calcPtAvgs().
00045 { 00046 return _particles; 00047 }
const FourMomentum & momentum | ( | ) | const [virtual] |
Get equivalent single momentum four-vector. (caches).
Implements ParticleBase.
Definition at line 165 of file Jet.cc.
References Jet::_calcMomVector(), and Jet::_momentum.
Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), ZEUS_2001_S4815815::analyze(), STAR_2009_UE_HELEN::analyze(), STAR_2006_S6870392::analyze(), MC_ZZJETS::analyze(), MC_WWJETS::analyze(), MC_TTBAR::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_JetAnalysis::analyze(), MC_DIJET::analyze(), D0_2009_S8349509::analyze(), D0_2009_S8202443::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7662670::analyze(), D0_2008_S6879055::analyze(), D0_1996_S3324664::analyze(), D0_1996_S3214044::analyze(), CMS_2011_S8957746::analyze(), CDF_2008_S7828950::analyze(), CDF_2008_S7782535::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_S7540469::analyze(), CDF_2007_S7057202::analyze(), CDF_2006_S6450792::analyze(), CDF_2004_S5839831::analyze(), CDF_2001_S4563131::analyze(), CDF_1998_S3618439::analyze(), CDF_1997_S3541940::analyze(), CDF_1996_S3349578::analyze(), CDF_1996_S3108457::analyze(), CDF_1993_S2742446::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2011_S8971293::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2010_S8817804::analyze(), ATLAS_2010_CONF_2010_049::analyze(), JetShape::calc(), Rivet::cmpJetsByAscAbsPseudorapidity(), Rivet::cmpJetsByAscAbsRapidity(), Rivet::cmpJetsByAscE(), Rivet::cmpJetsByAscP(), Rivet::cmpJetsByAscPseudorapidity(), Rivet::cmpJetsByAscRapidity(), Rivet::cmpJetsByDescAbsPseudorapidity(), Rivet::cmpJetsByDescAbsRapidity(), Rivet::cmpJetsByDescPseudorapidity(), Rivet::cmpJetsByDescRapidity(), Rivet::cmpJetsByE(), Rivet::cmpJetsByP(), Rivet::deltaEta(), Rivet::deltaPhi(), Rivet::deltaR(), Jet::eta(), Jet::EtSum(), JetAlg::jets(), Jet::phi(), Jet::ptSum(), and Jet::totalEnergy().
00165 { 00166 _calcMomVector(); 00167 return _momentum; 00168 }
double neutralEnergy | ( | ) | const |
Get the energy carried in this jet by neutral particles.
Definition at line 74 of file Jet.cc.
References FourMomentum::E(), Particle::momentum(), Jet::particles(), Particle::pdgId(), and Rivet::PID::threeCharge().
Referenced by STAR_2009_UE_HELEN::analyze().
00074 { 00075 double e_neutral = 0.0; 00076 foreach (const Particle& p, particles()) { 00077 const PdgId pid = p.pdgId(); 00078 if (PID::threeCharge(pid) == 0) { 00079 e_neutral += p.momentum().E(); 00080 } 00081 } 00082 return e_neutral; 00083 }
const vector<Particle>& particles | ( | ) | const [inline] |
Get the Rivet::Particles (full information) in this jet (const version).
Definition at line 60 of file Jet.hh.
References Jet::_fullParticles.
00060 { 00061 return _fullParticles; 00062 }
vector<Particle>& particles | ( | ) | [inline] |
Get the Rivet::Particles (full information) in this jet.
Definition at line 55 of file Jet.hh.
References Jet::_fullParticles.
Referenced by ATLAS_2010_CONF_2010_049::analyze(), JetShape::calc(), Jet::containsBottom(), Jet::containsCharm(), Jet::containsParticle(), Jet::containsParticleId(), Jet::hadronicEnergy(), and Jet::neutralEnergy().
00055 { 00056 return _fullParticles; 00057 }
double phi | ( | ) | const |
Get the unweighted average for this jet. (caches).
Definition at line 160 of file Jet.cc.
References Jet::momentum(), and FourVector::phi().
Referenced by Jet::_calcPtAvgs().
00160 { 00161 return momentum().phi(); 00162 }
double ptSum | ( | ) | const |
Get the sum of the values of the constituent tracks. (caches).
Definition at line 177 of file Jet.cc.
References Jet::momentum(), and FourMomentum::pT().
Referenced by CDF_2001_S4751469::analyze(), Rivet::cmpJetsByAscPt(), and Rivet::cmpJetsByPt().
00177 { 00178 return momentum().pT(); 00179 }
double ptWeightedEta | ( | ) | const |
Get the average for this jet, with the average weighted by the
values of the constituent tracks. (caches)
Definition at line 140 of file Jet.cc.
References Jet::_calcPtAvgs(), Jet::_okPtWeightedEta, and Jet::_ptWeightedEta.
00140 { 00141 _calcPtAvgs(); 00142 assert(_okPtWeightedEta); 00143 return _ptWeightedEta; 00144 }
double ptWeightedPhi | ( | ) | const |
Get the average for this jet, with the average weighted by the
values of the constituent tracks. (caches)
Definition at line 147 of file Jet.cc.
References Jet::_calcPtAvgs(), Jet::_okPtWeightedPhi, and Jet::_ptWeightedPhi.
Referenced by CDF_2001_S4751469::analyze().
00147 { 00148 _calcPtAvgs(); 00149 assert(_okPtWeightedPhi); 00150 return _ptWeightedPhi; 00151 }
Jet & setParticles | ( | const vector< FourMomentum > & | particles | ) |
Set the particles/tracks collection.
Definition at line 17 of file Jet.cc.
References Jet::_particles, and Jet::_resetCaches().
00017 { 00018 _particles = particles; 00019 _resetCaches(); 00020 return *this; 00021 }
size_t size | ( | ) | const [inline] |
Number of particles (tracks) in this jet.
Definition at line 65 of file Jet.hh.
References Jet::_particles.
00065 { 00066 return _particles.size(); 00067 }
double totalEnergy | ( | ) | const |
Get the total energy of this jet.
Definition at line 69 of file Jet.cc.
References FourMomentum::E(), and Jet::momentum().
Referenced by STAR_2009_UE_HELEN::analyze().
00069 { 00070 return momentum().E(); 00071 }
ParticleVector _fullParticles [private] |
Full particle information including tracks, ID etc.
Definition at line 159 of file Jet.hh.
Referenced by Jet::addParticle(), Jet::clear(), and Jet::particles().
FourMomentum _momentum [mutable, private] |
Cached effective jet 4-vector.
Definition at line 166 of file Jet.hh.
Referenced by Jet::_calcMomVector(), and Jet::momentum().
bool _okMomentum [mutable, private] |
Definition at line 167 of file Jet.hh.
Referenced by Jet::_calcMomVector(), and Jet::_resetCaches().
bool _okPtWeightedEta [mutable, private] |
Definition at line 163 of file Jet.hh.
Referenced by Jet::_calcPtAvgs(), Jet::_resetCaches(), and Jet::ptWeightedEta().
bool _okPtWeightedPhi [mutable, private] |
Definition at line 163 of file Jet.hh.
Referenced by Jet::_calcPtAvgs(), Jet::_resetCaches(), and Jet::ptWeightedPhi().
std::vector<FourMomentum> _particles [private] |
The particle tracks.
Definition at line 156 of file Jet.hh.
Referenced by Jet::addParticle(), Jet::begin(), Jet::clear(), Jet::end(), Jet::momenta(), Jet::setParticles(), and Jet::size().
double _ptWeightedEta [mutable, private] |
Definition at line 162 of file Jet.hh.
Referenced by Jet::_calcPtAvgs(), and Jet::ptWeightedEta().
double _ptWeightedPhi [mutable, private] |
Cached values of the -weighted
and
.
Definition at line 162 of file Jet.hh.
Referenced by Jet::_calcPtAvgs(), and Jet::ptWeightedPhi().