Jet.cc
Go to the documentation of this file.
00001 #include "Rivet/Jet.hh" 00002 #include "Rivet/Tools/Logging.hh" 00003 #include "Rivet/Tools/ParticleIdUtils.hh" 00004 #include "Rivet/Tools/RivetBoost.hh" 00005 #include "Rivet/ParticleName.hh" 00006 00007 namespace Rivet { 00008 00009 00010 Jet& Jet::clear() { 00011 _momentum = FourMomentum(); 00012 _pseudojet.reset(0,0,0,0); 00013 _particles.clear(); 00014 return *this; 00015 } 00016 00017 00018 Jet& Jet::setState(const FourMomentum& mom, const Particles& particles, const Particles& tags) { 00019 clear(); 00020 _momentum = mom; 00021 _pseudojet = fastjet::PseudoJet(mom.px(), mom.py(), mom.pz(), mom.E()); 00022 _particles = particles; 00023 _tags = tags; 00024 return *this; 00025 } 00026 00027 00028 Jet& Jet::setState(const fastjet::PseudoJet& pj, const vector<Particle>& particles, const Particles& tags) { 00029 clear(); 00030 _pseudojet = pj; 00031 _momentum = FourMomentum(pj.e(), pj.px(), pj.py(), pj.pz()); 00032 _particles = particles; 00033 _tags = tags; 00034 // if (_particles.empty()) { 00035 // foreach (const fastjet::PseudoJet pjc, _pseudojet.constituents()) { 00036 // // If there is no attached user info, we can't create a meaningful particle, so skip 00037 // if (!pjc.has_user_info<RivetFJInfo>()) continue; 00038 // const RivetFJInfo& fjinfo = pjc.user_info<RivetFJInfo>(); 00039 // // Don't add ghosts to the particles list 00040 // if (fjinfo.isGhost) continue; 00041 // // Otherwise construct a Particle from the PseudoJet, preferably from an associated GenParticle 00042 // ?if (fjinfo.genParticle != NULL) { 00043 // _particles.push_back(Particle(fjinfo.genParticle)); 00044 // } else { 00045 // if (fjinfo.pid == 0) continue; // skip if there is a null PID entry in the FJ info 00046 // const FourMomentum pjcmom(pjc.e(), pjc.px(), pjc.py(), pjc.pz()); 00047 // _particles.push_back(Particle(fjinfo.pid, pjcmom)); 00048 // } 00049 // } 00050 // } 00051 return *this; 00052 } 00053 00054 00055 Jet& Jet::setParticles(const vector<Particle>& particles) { 00056 _particles = particles; 00057 return *this; 00058 } 00059 00060 00061 bool Jet::containsParticle(const Particle& particle) const { 00062 const int barcode = particle.genParticle()->barcode(); 00063 foreach (const Particle& p, particles()) { 00064 if (p.genParticle()->barcode() == barcode) return true; 00065 } 00066 return false; 00067 } 00068 00069 00070 bool Jet::containsParticleId(PdgId pid) const { 00071 foreach (const Particle& p, particles()) { 00072 if (p.pid() == pid) return true; 00073 } 00074 return false; 00075 } 00076 00077 00078 bool Jet::containsParticleId(const vector<PdgId>& pids) const { 00079 foreach (const Particle& p, particles()) { 00080 foreach (PdgId pid, pids) { 00081 if (p.pid() == pid) return true; 00082 } 00083 } 00084 return false; 00085 } 00086 00087 00088 /// @todo Jet::containsMatch(Matcher m) { ... if m(pid) return true; ... } 00089 00090 00091 double Jet::neutralEnergy() const { 00092 double e_neutral = 0.0; 00093 foreach (const Particle& p, particles()) { 00094 const PdgId pid = p.pid(); 00095 if (PID::threeCharge(pid) == 0) { 00096 e_neutral += p.E(); 00097 } 00098 } 00099 return e_neutral; 00100 } 00101 00102 00103 double Jet::hadronicEnergy() const { 00104 double e_hadr = 0.0; 00105 foreach (const Particle& p, particles()) { 00106 const PdgId pid = p.pid(); 00107 if (PID::isHadron(pid)) { 00108 e_hadr += p.E(); 00109 } 00110 } 00111 return e_hadr; 00112 } 00113 00114 00115 bool Jet::containsCharm(bool include_decay_products) const { 00116 foreach (const Particle& p, particles()) { 00117 const PdgId pid = p.pid(); 00118 if (abs(pid) == PID::CQUARK) return true; 00119 if (PID::isHadron(pid) && PID::hasCharm(pid)) return true; 00120 if (include_decay_products) { 00121 HepMC::GenVertex* gv = p.genParticle()->production_vertex(); 00122 if (gv) { 00123 foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { 00124 const PdgId pid2 = pi->pdg_id(); 00125 if (PID::isHadron(pid2) && PID::hasCharm(pid2)) return true; 00126 } 00127 } 00128 } 00129 } 00130 return false; 00131 } 00132 00133 00134 bool Jet::containsBottom(bool include_decay_products) const { 00135 foreach (const Particle& p, particles()) { 00136 const PdgId pid = p.pid(); 00137 if (abs(pid) == PID::BQUARK) return true; 00138 if (PID::isHadron(pid) && PID::hasBottom(pid)) return true; 00139 if (include_decay_products) { 00140 HepMC::GenVertex* gv = p.genParticle()->production_vertex(); 00141 if (gv) { 00142 foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) { 00143 const PdgId pid2 = pi->pdg_id(); 00144 if (PID::isHadron(pid2) && PID::hasBottom(pid2)) return true; 00145 } 00146 } 00147 } 00148 } 00149 return false; 00150 } 00151 00152 00153 } Generated on Tue Sep 30 2014 19:45:45 for The Rivet MC analysis system by ![]() |