00001
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
00012 class Jet : public ParticleBase {
00013 public:
00014
00015
00016 Jet();
00017
00018
00019 typedef vector<FourMomentum>::iterator iterator;
00020
00021
00022 typedef vector<FourMomentum>::const_iterator const_iterator;
00023
00024
00025 iterator begin() {
00026 return _particles.begin();
00027 }
00028
00029
00030 iterator end() {
00031 return _particles.end();
00032 }
00033
00034
00035 const_iterator begin() const {
00036 return _particles.begin();
00037 }
00038
00039
00040 const_iterator end() const {
00041 return _particles.end();
00042 }
00043
00044
00045 vector<FourMomentum>& momenta() {
00046 return _particles;
00047 }
00048
00049
00050 const vector<FourMomentum>& momenta() const {
00051 return _particles;
00052 }
00053
00054
00055 vector<Particle>& particles() {
00056 return _fullParticles;
00057 }
00058
00059
00060 const vector<Particle>& particles() const {
00061 return _fullParticles;
00062 }
00063
00064
00065 size_t size() const {
00066 return _particles.size();
00067 }
00068
00069
00070 Jet& setParticles(const vector<FourMomentum>& particles);
00071
00072
00073 Jet& addParticle(const FourMomentum& particle);
00074
00075
00076 Jet& addParticle(const Particle& particle);
00077
00078
00079 bool containsParticle(const Particle& particle) const;
00080
00081
00082 bool containsParticleId(PdgId pid) const;
00083
00084
00085 bool containsParticleId(const vector<PdgId>& pids) const;
00086
00087
00088 bool containsCharm() const;
00089
00090
00091 bool containsBottom() const;
00092
00093
00094 Jet& clear();
00095
00096
00097
00098 double ptWeightedEta() const;
00099
00100
00101
00102 double ptWeightedPhi() const;
00103
00104
00105 double eta() const;
00106
00107
00108 double phi() const;
00109
00110
00111 const FourMomentum& momentum() const;
00112
00113
00114
00115
00116
00117
00118 public:
00119
00120
00121 double totalEnergy() const;
00122
00123
00124 double neutralEnergy() const;
00125
00126
00127 double hadronicEnergy() const;
00128
00129
00130 double ptSum() const;
00131
00132
00133 double EtSum() const;
00134
00135
00136 private:
00137
00138
00139
00140 void _resetCaches() const;
00141
00142
00143
00144 void _calcMomVector() const;
00145
00146
00147
00148
00149
00150 void _calcPtAvgs() const;
00151
00152
00153 private:
00154
00155
00156 std::vector<FourMomentum> _particles;
00157
00158
00159 ParticleVector _fullParticles;
00160
00161
00162 mutable double _ptWeightedPhi, _ptWeightedEta;
00163 mutable bool _okPtWeightedPhi, _okPtWeightedEta;
00164
00165
00166 mutable FourMomentum _momentum;
00167 mutable bool _okMomentum;
00168
00169 };
00170
00171
00172
00173 typedef std::vector<Jet> Jets;
00174
00175
00176
00177
00178
00179
00180
00181 inline bool cmpJetsByPt(const Jet& a, const Jet& b) {
00182 return a.ptSum() > b.ptSum();
00183 }
00184
00185
00186 inline bool cmpJetsByAscPt(const Jet& a, const Jet& b) {
00187 return a.ptSum() < b.ptSum();
00188 }
00189
00190
00191 inline bool cmpJetsByP(const Jet& a, const Jet& b) {
00192 return a.momentum().vector3().mod() > b.momentum().vector3().mod();
00193 }
00194
00195 inline bool cmpJetsByAscP(const Jet& a, const Jet& b) {
00196 return a.momentum().vector3().mod() < b.momentum().vector3().mod();
00197 }
00198
00199
00200
00201 inline bool cmpJetsByEt(const Jet& a, const Jet& b) {
00202 return a.EtSum() > b.EtSum();
00203 }
00204
00205
00206 inline bool cmpJetsByEtDesc(const Jet& a, const Jet& b) {
00207 return a.EtSum() < b.EtSum();
00208 }
00209
00210
00211
00212 inline bool cmpJetsByE(const Jet& a, const Jet& b) {
00213 return a.momentum().E() > b.momentum().E();
00214 }
00215
00216
00217 inline bool cmpJetsByAscE(const Jet& a, const Jet& b) {
00218 return a.momentum().E() < b.momentum().E();
00219 }
00220
00221
00222
00223 inline bool cmpJetsByDescPseudorapidity(const Jet& a, const Jet& b) {
00224 return a.momentum().pseudorapidity() > b.momentum().pseudorapidity();
00225 }
00226
00227
00228 inline bool cmpJetsByAscPseudorapidity(const Jet& a, const Jet& b) {
00229 return a.momentum().pseudorapidity() < b.momentum().pseudorapidity();
00230 }
00231
00232
00233
00234 inline bool cmpJetsByDescAbsPseudorapidity(const Jet& a, const Jet& b) {
00235 return fabs(a.momentum().pseudorapidity()) > fabs(b.momentum().pseudorapidity());
00236 }
00237
00238
00239 inline bool cmpJetsByAscAbsPseudorapidity(const Jet& a, const Jet& b) {
00240 return fabs(a.momentum().pseudorapidity()) < fabs(b.momentum().pseudorapidity());
00241 }
00242
00243
00244
00245 inline bool cmpJetsByDescRapidity(const Jet& a, const Jet& b) {
00246 return a.momentum().rapidity() > b.momentum().rapidity();
00247 }
00248
00249
00250 inline bool cmpJetsByAscRapidity(const Jet& a, const Jet& b) {
00251 return a.momentum().rapidity() < b.momentum().rapidity();
00252 }
00253
00254
00255
00256 inline bool cmpJetsByDescAbsRapidity(const Jet& a, const Jet& b) {
00257 return fabs(a.momentum().rapidity()) > fabs(b.momentum().rapidity());
00258 }
00259
00260
00261 inline bool cmpJetsByAscAbsRapidity(const Jet& a, const Jet& b) {
00262 return fabs(a.momentum().rapidity()) < fabs(b.momentum().rapidity());
00263 }
00264
00265
00266
00267 inline double deltaR(const Jet& j1, const Jet& j2,
00268 RapScheme scheme = PSEUDORAPIDITY) {
00269 return deltaR(j1.momentum(), j2.momentum(), scheme);
00270 }
00271
00272 inline double deltaR(const Jet& j, const Particle& p,
00273 RapScheme scheme = PSEUDORAPIDITY) {
00274 return deltaR(j.momentum(), p.momentum(), scheme);
00275 }
00276
00277 inline double deltaR(const Particle& p, const Jet& j,
00278 RapScheme scheme = PSEUDORAPIDITY) {
00279 return deltaR(p.momentum(), j.momentum(), scheme);
00280 }
00281
00282 inline double deltaR(const Jet& j, const FourMomentum& v,
00283 RapScheme scheme = PSEUDORAPIDITY) {
00284 return deltaR(j.momentum(), v, scheme);
00285 }
00286
00287 inline double deltaR(const Jet& j, const FourVector& v,
00288 RapScheme scheme = PSEUDORAPIDITY) {
00289 return deltaR(j.momentum(), v, scheme);
00290 }
00291
00292 inline double deltaR(const Jet& j, const Vector3& v) {
00293 return deltaR(j.momentum(), v);
00294 }
00295
00296 inline double deltaR(const Jet& j, double eta, double phi) {
00297 return deltaR(j.momentum(), eta, phi);
00298 }
00299
00300 inline double deltaR(const FourMomentum& v, const Jet& j,
00301 RapScheme scheme = PSEUDORAPIDITY) {
00302 return deltaR(v, j.momentum(), scheme);
00303 }
00304
00305 inline double deltaR(const FourVector& v, const Jet& j,
00306 RapScheme scheme = PSEUDORAPIDITY) {
00307 return deltaR(v, j.momentum(), scheme);
00308 }
00309
00310 inline double deltaR(const Vector3& v, const Jet& j) {
00311 return deltaR(v, j.momentum());
00312 }
00313
00314 inline double deltaR(double eta, double phi, const Jet& j) {
00315 return deltaR(eta, phi, j.momentum());
00316 }
00317
00318
00319 inline double deltaPhi(const Jet& j1, const Jet& j2) {
00320 return deltaPhi(j1.momentum(), j2.momentum());
00321 }
00322
00323 inline double deltaPhi(const Jet& j, const Particle& p) {
00324 return deltaPhi(j.momentum(), p.momentum());
00325 }
00326
00327 inline double deltaPhi(const Particle& p, const Jet& j) {
00328 return deltaPhi(p.momentum(), j.momentum());
00329 }
00330
00331 inline double deltaPhi(const Jet& j, const FourMomentum& v) {
00332 return deltaPhi(j.momentum(), v);
00333 }
00334
00335 inline double deltaPhi(const Jet& j, const FourVector& v) {
00336 return deltaPhi(j.momentum(), v);
00337 }
00338
00339 inline double deltaPhi(const Jet& j, const Vector3& v) {
00340 return deltaPhi(j.momentum(), v);
00341 }
00342
00343 inline double deltaPhi(const Jet& j, double phi) {
00344 return deltaPhi(j.momentum(), phi);
00345 }
00346
00347 inline double deltaPhi(const FourMomentum& v, const Jet& j) {
00348 return deltaPhi(v, j.momentum());
00349 }
00350
00351 inline double deltaPhi(const FourVector& v, const Jet& j) {
00352 return deltaPhi(v, j.momentum());
00353 }
00354
00355 inline double deltaPhi(const Vector3& v, const Jet& j) {
00356 return deltaPhi(v, j.momentum());
00357 }
00358
00359 inline double deltaPhi(double phi, const Jet& j) {
00360 return deltaPhi(phi, j.momentum());
00361 }
00362
00363
00364 inline double deltaEta(const Jet& j1, const Jet& j2) {
00365 return deltaEta(j1.momentum(), j2.momentum());
00366 }
00367
00368 inline double deltaEta(const Jet& j, const Particle& p) {
00369 return deltaEta(j.momentum(), p.momentum());
00370 }
00371
00372 inline double deltaEta(const Particle& p, const Jet& j) {
00373 return deltaEta(p.momentum(), j.momentum());
00374 }
00375
00376 inline double deltaEta(const Jet& j, const FourMomentum& v) {
00377 return deltaEta(j.momentum(), v);
00378 }
00379
00380 inline double deltaEta(const Jet& j, const FourVector& v) {
00381 return deltaEta(j.momentum(), v);
00382 }
00383
00384 inline double deltaEta(const Jet& j, const Vector3& v) {
00385 return deltaEta(j.momentum(), v);
00386 }
00387
00388 inline double deltaEta(const Jet& j, double eta) {
00389 return deltaEta(j.momentum(), eta);
00390 }
00391
00392 inline double deltaEta(const FourMomentum& v, const Jet& j) {
00393 return deltaEta(v, j.momentum());
00394 }
00395
00396 inline double deltaEta(const FourVector& v, const Jet& j) {
00397 return deltaEta(v, j.momentum());
00398 }
00399
00400 inline double deltaEta(const Vector3& v, const Jet& j) {
00401 return deltaEta(v, j.momentum());
00402 }
00403
00404 inline double deltaEta(double eta, const Jet& j) {
00405 return deltaEta(eta, j.momentum());
00406 }
00407
00408 }
00409
00410 #endif