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 /// @todo Add a constructor from PseudoJet 00027 00028 //@} 00029 00030 00031 /// @name Access jet constituents 00032 //@{ 00033 00034 /// Number of particles in this jet. 00035 size_t size() const { return _particles.size(); } 00036 00037 /// Get the particles in this jet. 00038 vector<Particle>& particles() { return _particles; } 00039 00040 /// Get the particles in this jet (const version) 00041 const vector<Particle>& particles() const { return _particles; } 00042 00043 /// Check whether this jet contains a particular particle. 00044 bool containsParticle(const Particle& particle) const; 00045 00046 /// Check whether this jet contains a certain particle type. 00047 bool containsParticleId(PdgId pid) const; 00048 00049 /// Check whether this jet contains at least one of certain particle types. 00050 bool containsParticleId(const vector<PdgId>& pids) const; 00051 00052 /// Check whether this jet contains a charm-flavoured hadron (or decay products from one). 00053 bool containsCharm() const; 00054 00055 /// Check whether this jet contains a bottom-flavoured hadron (or decay products from one). 00056 bool containsBottom() const; 00057 00058 //@} 00059 00060 00061 /// @name Access additional effective jet 4-vector properties 00062 //@{ 00063 00064 /// Get equivalent single momentum four-vector. 00065 const FourMomentum& momentum() const { return _momentum; } 00066 00067 /// Get the total energy of this jet. 00068 double totalEnergy() const { return momentum().E(); } 00069 00070 /// Get the energy carried in this jet by neutral particles. 00071 double neutralEnergy() const; 00072 00073 /// Get the energy carried in this jet by hadrons. 00074 double hadronicEnergy() const; 00075 00076 //@} 00077 00078 00079 /// @name Set the jet constituents and properties 00080 //@{ 00081 00082 /// Set all the jet data, with full particle information. 00083 Jet& setState(const vector<Particle>& particles, const FourMomentum& pjet); 00084 00085 /// Set the effective 4-momentum of the jet. 00086 Jet& setMomentum(const FourMomentum& momentum); 00087 00088 /// Set the particles collection with full particle information. 00089 Jet& setParticles(const vector<Particle>& particles); 00090 00091 /// Reset this jet as empty. 00092 Jet& clear(); 00093 00094 //@} 00095 00096 00097 private: 00098 00099 /// Full particle information including tracks, ID etc 00100 Particles _particles; 00101 00102 /// Effective jet 4-vector 00103 FourMomentum _momentum; 00104 00105 /// @todo Add a FJ3 PseudoJet member to unify PseudoJet and Jet 00106 00107 }; 00108 00109 00110 /// Typedef for a collection of Jet objects. 00111 typedef vector<Jet> Jets; 00112 00113 00114 /// @name deltaR, deltaEta, deltaPhi functions specifically for Jet arguments 00115 //@{ 00116 00117 inline double deltaR(const Jet& j1, const Jet& j2, 00118 RapScheme scheme = PSEUDORAPIDITY) { 00119 return deltaR(j1.momentum(), j2.momentum(), scheme); 00120 } 00121 00122 inline double deltaR(const Jet& j, const Particle& p, 00123 RapScheme scheme = PSEUDORAPIDITY) { 00124 return deltaR(j.momentum(), p.momentum(), scheme); 00125 } 00126 00127 inline double deltaR(const Particle& p, const Jet& j, 00128 RapScheme scheme = PSEUDORAPIDITY) { 00129 return deltaR(p.momentum(), j.momentum(), scheme); 00130 } 00131 00132 inline double deltaR(const Jet& j, const FourMomentum& v, 00133 RapScheme scheme = PSEUDORAPIDITY) { 00134 return deltaR(j.momentum(), v, scheme); 00135 } 00136 00137 inline double deltaR(const Jet& j, const FourVector& v, 00138 RapScheme scheme = PSEUDORAPIDITY) { 00139 return deltaR(j.momentum(), v, scheme); 00140 } 00141 00142 inline double deltaR(const Jet& j, const Vector3& v) { 00143 return deltaR(j.momentum(), v); 00144 } 00145 00146 inline double deltaR(const Jet& j, double eta, double phi) { 00147 return deltaR(j.momentum(), eta, phi); 00148 } 00149 00150 inline double deltaR(const FourMomentum& v, const Jet& j, 00151 RapScheme scheme = PSEUDORAPIDITY) { 00152 return deltaR(v, j.momentum(), scheme); 00153 } 00154 00155 inline double deltaR(const FourVector& v, const Jet& j, 00156 RapScheme scheme = PSEUDORAPIDITY) { 00157 return deltaR(v, j.momentum(), scheme); 00158 } 00159 00160 inline double deltaR(const Vector3& v, const Jet& j) { 00161 return deltaR(v, j.momentum()); 00162 } 00163 00164 inline double deltaR(double eta, double phi, const Jet& j) { 00165 return deltaR(eta, phi, j.momentum()); 00166 } 00167 00168 00169 inline double deltaPhi(const Jet& j1, const Jet& j2) { 00170 return deltaPhi(j1.momentum(), j2.momentum()); 00171 } 00172 00173 inline double deltaPhi(const Jet& j, const Particle& p) { 00174 return deltaPhi(j.momentum(), p.momentum()); 00175 } 00176 00177 inline double deltaPhi(const Particle& p, const Jet& j) { 00178 return deltaPhi(p.momentum(), j.momentum()); 00179 } 00180 00181 inline double deltaPhi(const Jet& j, const FourMomentum& v) { 00182 return deltaPhi(j.momentum(), v); 00183 } 00184 00185 inline double deltaPhi(const Jet& j, const FourVector& v) { 00186 return deltaPhi(j.momentum(), v); 00187 } 00188 00189 inline double deltaPhi(const Jet& j, const Vector3& v) { 00190 return deltaPhi(j.momentum(), v); 00191 } 00192 00193 inline double deltaPhi(const Jet& j, double phi) { 00194 return deltaPhi(j.momentum(), phi); 00195 } 00196 00197 inline double deltaPhi(const FourMomentum& v, const Jet& j) { 00198 return deltaPhi(v, j.momentum()); 00199 } 00200 00201 inline double deltaPhi(const FourVector& v, const Jet& j) { 00202 return deltaPhi(v, j.momentum()); 00203 } 00204 00205 inline double deltaPhi(const Vector3& v, const Jet& j) { 00206 return deltaPhi(v, j.momentum()); 00207 } 00208 00209 inline double deltaPhi(double phi, const Jet& j) { 00210 return deltaPhi(phi, j.momentum()); 00211 } 00212 00213 00214 inline double deltaEta(const Jet& j1, const Jet& j2) { 00215 return deltaEta(j1.momentum(), j2.momentum()); 00216 } 00217 00218 inline double deltaEta(const Jet& j, const Particle& p) { 00219 return deltaEta(j.momentum(), p.momentum()); 00220 } 00221 00222 inline double deltaEta(const Particle& p, const Jet& j) { 00223 return deltaEta(p.momentum(), j.momentum()); 00224 } 00225 00226 inline double deltaEta(const Jet& j, const FourMomentum& v) { 00227 return deltaEta(j.momentum(), v); 00228 } 00229 00230 inline double deltaEta(const Jet& j, const FourVector& v) { 00231 return deltaEta(j.momentum(), v); 00232 } 00233 00234 inline double deltaEta(const Jet& j, const Vector3& v) { 00235 return deltaEta(j.momentum(), v); 00236 } 00237 00238 inline double deltaEta(const Jet& j, double eta) { 00239 return deltaEta(j.momentum(), eta); 00240 } 00241 00242 inline double deltaEta(const FourMomentum& v, const Jet& j) { 00243 return deltaEta(v, j.momentum()); 00244 } 00245 00246 inline double deltaEta(const FourVector& v, const Jet& j) { 00247 return deltaEta(v, j.momentum()); 00248 } 00249 00250 inline double deltaEta(const Vector3& v, const Jet& j) { 00251 return deltaEta(v, j.momentum()); 00252 } 00253 00254 inline double deltaEta(double eta, const Jet& j) { 00255 return deltaEta(eta, j.momentum()); 00256 } 00257 00258 //@} 00259 00260 00261 } 00262 00263 #endif Generated on Thu Feb 6 2014 17:38:45 for The Rivet MC analysis system by 1.7.6.1 |