Jet.hh
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #ifndef RIVET_Jet_HH 00003 #define RIVET_Jet_HH 00004 00005 #include "Rivet/Rivet.hh" 00006 #include <numeric> 00007 00008 namespace Rivet { 00009 00010 00011 /// @brief Representation of a clustered jet of particles. 00012 class Jet : public ParticleBase { 00013 public: 00014 00015 /// @name Constructors 00016 //@{ 00017 00018 Jet() : ParticleBase() { clear(); } 00019 00020 /// Set all the jet data, with full particle information. 00021 Jet(const vector<Particle>& particles, const FourMomentum& pjet) 00022 : ParticleBase() { 00023 setState(particles, pjet); 00024 } 00025 00026 //@} 00027 00028 00029 /// @name Access jet constituents 00030 //@{ 00031 00032 /// Number of particles in this jet. 00033 size_t size() const { return _particles.size(); } 00034 00035 /// Get the particles in this jet. 00036 vector<Particle>& particles() { return _particles; } 00037 00038 /// Get the particles in this jet (const version) 00039 const vector<Particle>& particles() const { return _particles; } 00040 00041 /// Check whether this jet contains a particular particle. 00042 bool containsParticle(const Particle& particle) const; 00043 00044 /// Check whether this jet contains a certain particle type. 00045 bool containsParticleId(PdgId pid) const; 00046 00047 /// Check whether this jet contains at least one of certain particle types. 00048 bool containsParticleId(const vector<PdgId>& pids) const; 00049 00050 /// Check whether this jet contains a charm-flavoured hadron (or decay products from one). 00051 bool containsCharm() const; 00052 00053 /// Check whether this jet contains a bottom-flavoured hadron (or decay products from one). 00054 bool containsBottom() const; 00055 00056 //@} 00057 00058 00059 /// @name Access additional effective jet 4-vector properties 00060 //@{ 00061 00062 /// Get equivalent single momentum four-vector. 00063 const FourMomentum& momentum() const { return _momentum; } 00064 00065 /// Get the total energy of this jet. 00066 double totalEnergy() const { return momentum().E(); } 00067 00068 /// Get the energy carried in this jet by neutral particles. 00069 double neutralEnergy() const; 00070 00071 /// Get the energy carried in this jet by hadrons. 00072 double hadronicEnergy() const; 00073 00074 //@} 00075 00076 00077 /// @name Set the jet constituents and properties 00078 //@{ 00079 00080 /// Set all the jet data, with full particle information. 00081 Jet& setState(const vector<Particle>& particles, const FourMomentum& pjet); 00082 00083 /// Set the effective 4-momentum of the jet. 00084 Jet& setMomentum(const FourMomentum& momentum); 00085 00086 /// Set the particles collection with full particle information. 00087 Jet& setParticles(const vector<Particle>& particles); 00088 00089 /// Reset this jet as empty. 00090 Jet& clear(); 00091 00092 //@} 00093 00094 00095 private: 00096 00097 /// Full particle information including tracks, ID etc 00098 Particles _particles; 00099 00100 /// Effective jet 4-vector 00101 FourMomentum _momentum; 00102 00103 }; 00104 00105 00106 /// Typedef for a collection of Jet objects. 00107 typedef std::vector<Jet> Jets; 00108 00109 00110 /// @name Jet comparison functions for STL sorting 00111 //@{ 00112 00113 /// @brief Compare jets by \f$ p_\perp \f$ (descending - usual sorting for HEP) 00114 /// Use this so that highest \f$ p_\perp \f$ is at the front of the list 00115 inline bool cmpJetsByPt(const Jet& a, const Jet& b) { 00116 return a.pT() > b.pT(); 00117 } 00118 /// @brief Compare jets by \f$ p_\perp \f$ (ascending) 00119 /// Use this so that lowest \f$ p_\perp \f$ is at the front of the list 00120 inline bool cmpJetsByAscPt(const Jet& a, const Jet& b) { 00121 return a.pT() < b.pT(); 00122 } 00123 00124 /// @brief Compare jets by descending momentum, \f$ p \f$ 00125 inline bool cmpJetsByP(const Jet& a, const Jet& b) { 00126 return a.momentum().vector3().mod() > b.momentum().vector3().mod(); 00127 } 00128 /// @brief Compare jets by ascending momentum, \f$ p \f$ 00129 inline bool cmpJetsByAscP(const Jet& a, const Jet& b) { 00130 return a.momentum().vector3().mod() < b.momentum().vector3().mod(); 00131 } 00132 00133 // @brief Compare jets by \f$ E_\perp \f$ (descending - usual sorting for HEP) 00134 /// Use this so that highest \f$ E_\perp \f$ is at the front of the list 00135 inline bool cmpJetsByEt(const Jet& a, const Jet& b) { 00136 return a.Et() > b.Et(); 00137 } 00138 // @brief Compare jets by \f$ E_\perp \f$ (ascending) 00139 /// Use this so that lowest \f$ E_\perp \f$ is at the front of the list 00140 inline bool cmpJetsByEtDesc(const Jet& a, const Jet& b) { 00141 return a.Et() < b.Et(); 00142 } 00143 00144 /// @brief Compare jets by \f$ E \f$ (descending - usual sorting for HEP) 00145 /// Use this so that highest \f$ E \f$ is at the front of the list 00146 inline bool cmpJetsByE(const Jet& a, const Jet& b) { 00147 return a.E() > b.E(); 00148 } 00149 /// @brief Compare jets by \f$ E \f$ (ascending) 00150 /// Use this so that lowest \f$ E \f$ is at the front of the list 00151 inline bool cmpJetsByAscE(const Jet& a, const Jet& b) { 00152 return a.E() < b.E(); 00153 } 00154 00155 /// @brief Compare jets by \f$ \eta \f$ (descending) 00156 /// Use this so that highest \f$ \eta \f$ is at the front of the list 00157 inline bool cmpJetsByDescPseudorapidity(const Jet& a, const Jet& b) { 00158 return a.eta() > b.eta(); 00159 } 00160 /// @brief Compare jets by \f$ \eta \f$ (ascending) 00161 /// Use this so that lowest \f$ \eta \f$ is at the front of the list 00162 inline bool cmpJetsByAscPseudorapidity(const Jet& a, const Jet& b) { 00163 return a.eta() < b.eta(); 00164 } 00165 00166 /// @brief Compare jets by \f$ |\eta| \f$ (descending) 00167 /// Use this so that highest \f$ |\eta| \f$ is at the front of the list 00168 inline bool cmpJetsByDescAbsPseudorapidity(const Jet& a, const Jet& b) { 00169 return fabs(a.eta()) > fabs(b.eta()); 00170 } 00171 /// @brief Compare jets by \f$ |\eta| \f$ (ascending) 00172 /// Use this so that lowest \f$ |\eta| \f$ is at the front of the list 00173 inline bool cmpJetsByAscAbsPseudorapidity(const Jet& a, const Jet& b) { 00174 return fabs(a.eta()) < fabs(b.eta()); 00175 } 00176 00177 /// @brief Compare jets by \f$ y \f$ (descending) 00178 /// Use this so that highest \f$ y \f$ is at the front of the list 00179 inline bool cmpJetsByDescRapidity(const Jet& a, const Jet& b) { 00180 return a.rapidity() > b.rapidity(); 00181 } 00182 /// @brief Compare jets by \f$ y \f$ (ascending) 00183 /// Use this so that lowest \f$ y \f$ is at the front of the list 00184 inline bool cmpJetsByAscRapidity(const Jet& a, const Jet& b) { 00185 return a.rapidity() < b.rapidity(); 00186 } 00187 00188 /// @brief Compare jets by \f$ |y| \f$ (descending) 00189 /// Use this so that highest \f$ |y| \f$ is at the front of the list 00190 inline bool cmpJetsByDescAbsRapidity(const Jet& a, const Jet& b) { 00191 return fabs(a.rapidity()) > fabs(b.rapidity()); 00192 } 00193 /// @brief Compare jets by \f$ |y| \f$ (ascending) 00194 /// Use this so that lowest \f$ |y| \f$ is at the front of the list 00195 inline bool cmpJetsByAscAbsRapidity(const Jet& a, const Jet& b) { 00196 return fabs(a.rapidity()) < fabs(b.rapidity()); 00197 } 00198 00199 //@} 00200 00201 inline double deltaR(const Jet& j1, const Jet& j2, 00202 RapScheme scheme = PSEUDORAPIDITY) { 00203 return deltaR(j1.momentum(), j2.momentum(), scheme); 00204 } 00205 00206 inline double deltaR(const Jet& j, const Particle& p, 00207 RapScheme scheme = PSEUDORAPIDITY) { 00208 return deltaR(j.momentum(), p.momentum(), scheme); 00209 } 00210 00211 inline double deltaR(const Particle& p, const Jet& j, 00212 RapScheme scheme = PSEUDORAPIDITY) { 00213 return deltaR(p.momentum(), j.momentum(), scheme); 00214 } 00215 00216 inline double deltaR(const Jet& j, const FourMomentum& v, 00217 RapScheme scheme = PSEUDORAPIDITY) { 00218 return deltaR(j.momentum(), v, scheme); 00219 } 00220 00221 inline double deltaR(const Jet& j, const FourVector& v, 00222 RapScheme scheme = PSEUDORAPIDITY) { 00223 return deltaR(j.momentum(), v, scheme); 00224 } 00225 00226 inline double deltaR(const Jet& j, const Vector3& v) { 00227 return deltaR(j.momentum(), v); 00228 } 00229 00230 inline double deltaR(const Jet& j, double eta, double phi) { 00231 return deltaR(j.momentum(), eta, phi); 00232 } 00233 00234 inline double deltaR(const FourMomentum& v, const Jet& j, 00235 RapScheme scheme = PSEUDORAPIDITY) { 00236 return deltaR(v, j.momentum(), scheme); 00237 } 00238 00239 inline double deltaR(const FourVector& v, const Jet& j, 00240 RapScheme scheme = PSEUDORAPIDITY) { 00241 return deltaR(v, j.momentum(), scheme); 00242 } 00243 00244 inline double deltaR(const Vector3& v, const Jet& j) { 00245 return deltaR(v, j.momentum()); 00246 } 00247 00248 inline double deltaR(double eta, double phi, const Jet& j) { 00249 return deltaR(eta, phi, j.momentum()); 00250 } 00251 00252 00253 inline double deltaPhi(const Jet& j1, const Jet& j2) { 00254 return deltaPhi(j1.momentum(), j2.momentum()); 00255 } 00256 00257 inline double deltaPhi(const Jet& j, const Particle& p) { 00258 return deltaPhi(j.momentum(), p.momentum()); 00259 } 00260 00261 inline double deltaPhi(const Particle& p, const Jet& j) { 00262 return deltaPhi(p.momentum(), j.momentum()); 00263 } 00264 00265 inline double deltaPhi(const Jet& j, const FourMomentum& v) { 00266 return deltaPhi(j.momentum(), v); 00267 } 00268 00269 inline double deltaPhi(const Jet& j, const FourVector& v) { 00270 return deltaPhi(j.momentum(), v); 00271 } 00272 00273 inline double deltaPhi(const Jet& j, const Vector3& v) { 00274 return deltaPhi(j.momentum(), v); 00275 } 00276 00277 inline double deltaPhi(const Jet& j, double phi) { 00278 return deltaPhi(j.momentum(), phi); 00279 } 00280 00281 inline double deltaPhi(const FourMomentum& v, const Jet& j) { 00282 return deltaPhi(v, j.momentum()); 00283 } 00284 00285 inline double deltaPhi(const FourVector& v, const Jet& j) { 00286 return deltaPhi(v, j.momentum()); 00287 } 00288 00289 inline double deltaPhi(const Vector3& v, const Jet& j) { 00290 return deltaPhi(v, j.momentum()); 00291 } 00292 00293 inline double deltaPhi(double phi, const Jet& j) { 00294 return deltaPhi(phi, j.momentum()); 00295 } 00296 00297 00298 inline double deltaEta(const Jet& j1, const Jet& j2) { 00299 return deltaEta(j1.momentum(), j2.momentum()); 00300 } 00301 00302 inline double deltaEta(const Jet& j, const Particle& p) { 00303 return deltaEta(j.momentum(), p.momentum()); 00304 } 00305 00306 inline double deltaEta(const Particle& p, const Jet& j) { 00307 return deltaEta(p.momentum(), j.momentum()); 00308 } 00309 00310 inline double deltaEta(const Jet& j, const FourMomentum& v) { 00311 return deltaEta(j.momentum(), v); 00312 } 00313 00314 inline double deltaEta(const Jet& j, const FourVector& v) { 00315 return deltaEta(j.momentum(), v); 00316 } 00317 00318 inline double deltaEta(const Jet& j, const Vector3& v) { 00319 return deltaEta(j.momentum(), v); 00320 } 00321 00322 inline double deltaEta(const Jet& j, double eta) { 00323 return deltaEta(j.momentum(), eta); 00324 } 00325 00326 inline double deltaEta(const FourMomentum& v, const Jet& j) { 00327 return deltaEta(v, j.momentum()); 00328 } 00329 00330 inline double deltaEta(const FourVector& v, const Jet& j) { 00331 return deltaEta(v, j.momentum()); 00332 } 00333 00334 inline double deltaEta(const Vector3& v, const Jet& j) { 00335 return deltaEta(v, j.momentum()); 00336 } 00337 00338 inline double deltaEta(double eta, const Jet& j) { 00339 return deltaEta(eta, j.momentum()); 00340 } 00341 00342 } 00343 00344 #endif Generated on Fri Oct 25 2013 12:41:45 for The Rivet MC analysis system by ![]() |