Jet Class Reference Representation of a clustered jet of particles. More...
Inheritance diagram for Jet:
![]()
Collaboration diagram for Jet:
![]()
Detailed DescriptionConstructor & Destructor Documentation
Set all the jet data, with full particle information. Definition at line 21 of file Jet.hh. References Jet::setState(). : ParticleBase() { setState(particles, pjet); } Member Function DocumentationReset this jet as empty. Definition at line 146 of file Jet.cc. References Jet::_momentum, and Jet::_particles. Referenced by Jet::Jet(). { //_momenta.clear(); _particles.clear(); _momentum = FourMomentum(); return *this; }
Check whether this jet contains a bottom-flavoured hadron (or decay products from one). Definition at line 129 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 CDF_2008_S7782535::analyze(), EXAMPLE::analyze(), ATLAS_2012_I1126136::analyze(), MC_TTBAR::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_I1095236::analyze(), and MC_VH2BB::analyze(). { foreach (const Particle& p, particles()) { const PdgId pid = p.pdgId(); if (abs(pid) == BQUARK) return true; if (PID::isHadron(pid) && PID::hasBottom(pid)) return true; HepMC::GenVertex* gv = p.genParticle().production_vertex(); if (gv) { foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { const PdgId pid2 = pi->pdg_id(); if (PID::isHadron(pid2) && PID::hasBottom(pid2)) return true; } } } return false; }
Check whether this jet contains a charm-flavoured hadron (or decay products from one). Definition at line 112 of file Jet.cc. References Rivet::CQUARK, Particle::genParticle(), Rivet::PID::hasCharm(), Rivet::PID::isHadron(), Rivet::particles(), Jet::particles(), Particle::pdgId(), and Rivet::pi. { foreach (const Particle& p, particles()) { const PdgId pid = p.pdgId(); if (abs(pid) == CQUARK) return true; if (PID::isHadron(pid) && PID::hasCharm(pid)) return true; HepMC::GenVertex* gv = p.genParticle().production_vertex(); if (gv) { foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { const PdgId pid2 = pi->pdg_id(); if (PID::isHadron(pid2) && PID::hasCharm(pid2)) return true; } } } return false; }
Check whether this jet contains a particular particle. Definition at line 58 of file Jet.cc. References Particle::genParticle(), and Jet::particles(). { const int barcode = particle.genParticle().barcode(); foreach (const Particle& p, particles()) { if (p.genParticle().barcode() == barcode) return true; } return false; }
Check whether this jet contains a certain particle type. Definition at line 67 of file Jet.cc. References Jet::particles(), and Particle::pdgId(). { foreach (const Particle& p, particles()) { if (p.pdgId() == pid) return true; } return false; }
Check whether this jet contains at least one of certain particle types. Definition at line 75 of file Jet.cc. References Jet::particles(), and Particle::pdgId().
Get the unweighted average Definition at line 108 of file Jet.hh. References FourVector::eta(), and Jet::momentum(). Referenced by STAR_2006_S6870392::analyze(), ATLAS_2012_I1082009::analyze(), and ATLAS_2010_S8919674::analyze().
Get the sum of the Definition at line 126 of file Jet.hh. References FourMomentum::Et(), and Jet::momentum(). Referenced by CDF_2004_S5839831::analyze(), Rivet::cmpJetsByEt(), and Rivet::cmpJetsByEtDesc().
Get the energy carried in this jet by hadrons. Definition at line 100 of file Jet.cc. References FourMomentum::E(), Rivet::PID::isHadron(), Particle::momentum(), Jet::particles(), and Particle::pdgId(). { double e_hadr = 0.0; foreach (const Particle& p, particles()) { const PdgId pid = p.pdgId(); if (PID::isHadron(pid)) { e_hadr += p.momentum().E(); } } return e_hadr; }
Get equivalent single momentum four-vector. Implements ParticleBase. Definition at line 105 of file Jet.hh. References Jet::_momentum. Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), MC_JetAnalysis::analyze(), CMS_2012_I1087342::analyze(), CMS_2011_S9088458::analyze(), CDF_2006_S6450792::analyze(), CMS_2011_S9086218::analyze(), CMS_2011_S8957746::analyze(), STAR_2006_S6870392::analyze(), CDF_2007_S7057202::analyze(), D0_1996_S3324664::analyze(), CDF_1998_S3618439::analyze(), CDF_2001_S4563131::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CDF_2008_S7828950::analyze(), D0_2008_S6879055::analyze(), ZEUS_2001_S4815815::analyze(), CDF_1993_S2742446::analyze(), CDF_2008_S7540469::analyze(), MC_DIJET::analyze(), CDF_2008_S7782535::analyze(), ATLAS_2011_S8971293::analyze(), CDF_1996_S3108457::analyze(), CDF_1997_S3541940::analyze(), D0_2008_S7863608::analyze(), ATLAS_2010_S8817804::analyze(), D0_2008_S7662670::analyze(), ATLAS_2012_I1082009::analyze(), D0_2009_S8349509::analyze(), STAR_2009_UE_HELEN::analyze(), D0_2009_S8202443::analyze(), ATLAS_2012_I1082936::analyze(), MC_PHOTONJETUE::analyze(), CDF_2008_S7541902::analyze(), D0_1996_S3214044::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1186556::analyze(), MC_TTBAR::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), CDF_1996_S3349578::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2011_S9128077::analyze(), ATLAS_2012_I1125961::analyze(), MC_ZZJETS::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1117704::analyze(), MC_WWJETS::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2012_CONF_2012_153::analyze(), MC_SUSY::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2011_S9041966::analyze(), ATLAS_2012_I1094568::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2011_S9126244::analyze(), ATLAS_2011_I919017::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(), ATLAS_2011_I945498::selectJets(), Jet::setMomentum(), and Jet::totalEnergy(). { return _momentum; }
Get the energy carried in this jet by neutral particles.
Definition at line 88 of file Jet.cc. References FourMomentum::E(), Particle::momentum(), Jet::particles(), Particle::pdgId(), and Rivet::PID::threeCharge(). Referenced by STAR_2009_UE_HELEN::analyze(). { double e_neutral = 0.0; foreach (const Particle& p, particles()) { const PdgId pid = p.pdgId(); if (PID::threeCharge(pid) == 0) { e_neutral += p.momentum().E(); } } return e_neutral; } Get the particles in this jet. Definition at line 78 of file Jet.hh. References Jet::_particles. Referenced by ATLAS_2010_CONF_2010_049::analyze(), ATLAS_2011_I919017::analyze(), JetShape::calc(), Jet::containsBottom(), Jet::containsCharm(), Jet::containsParticle(), Jet::containsParticleId(), Jet::hadronicEnergy(), Jet::neutralEnergy(), and Jet::setParticles(). { return _particles; } Get the particles in this jet (const version) Definition at line 81 of file Jet.hh. References Jet::_particles. { return _particles; }
Get the unweighted average Definition at line 111 of file Jet.hh. References Jet::momentum(), and FourVector::phi(). Referenced by CDF_2001_S4751469::analyze().
Get the sum of the Definition at line 123 of file Jet.hh. References Jet::momentum(), and FourMomentum::pT(). Referenced by CDF_2001_S4751469::analyze(), Rivet::cmpJetsByAscPt(), and Rivet::cmpJetsByPt().
Set the effective 4-momentum of the jet. Definition at line 24 of file Jet.cc. References Jet::_momentum, and Jet::momentum(). Referenced by Jet::setState().
Set the particles collection with full particle information. Definition at line 30 of file Jet.cc. References Jet::_particles, and Jet::particles(). Referenced by Jet::setState(). { _particles = particles; // foreach (const Particle& p, particles) { // _momenta.push_back(p.momentum()); // } return *this; }
Set all the jet data, with full particle information. Definition at line 10 of file Jet.cc. References Jet::setMomentum(), and Jet::setParticles(). Referenced by Jet::Jet(). { setParticles(particles); setMomentum(pjet); return *this; }
Number of particles in this jet. Definition at line 39 of file Jet.hh. References Jet::_particles. Referenced by MC_TTBAR::analyze(). { return _particles.size(); }
Get the total energy of this jet. Definition at line 114 of file Jet.hh. References FourMomentum::E(), and Jet::momentum(). Referenced by STAR_2009_UE_HELEN::analyze(). Member Data Documentation
Effective jet 4-vector. Definition at line 171 of file Jet.hh. Referenced by Jet::clear(), Jet::momentum(), and Jet::setMomentum().
Full particle information including tracks, ID etc. Definition at line 164 of file Jet.hh. Referenced by Jet::clear(), Jet::particles(), Jet::setParticles(), and Jet::size(). The documentation for this class was generated from the following files: Generated on Fri Dec 21 2012 15:03:51 for The Rivet MC analysis system by ![]() |