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
00192 inline bool cmpJetsByEt(const Jet& a, const Jet& b) {
00193 return a.EtSum() > b.EtSum();
00194 }
00195
00196
00197 inline bool cmpJetsByEtDesc(const Jet& a, const Jet& b) {
00198 return a.EtSum() < b.EtSum();
00199 }
00200
00201
00202
00203 inline bool cmpJetsByE(const Jet& a, const Jet& b) {
00204 return a.momentum().E() > b.momentum().E();
00205 }
00206
00207
00208 inline bool cmpJetsByAscE(const Jet& a, const Jet& b) {
00209 return a.momentum().E() < b.momentum().E();
00210 }
00211
00212
00213
00214 inline bool cmpJetsByDescPseudorapidity(const Jet& a, const Jet& b) {
00215 return a.momentum().pseudorapidity() > b.momentum().pseudorapidity();
00216 }
00217
00218
00219 inline bool cmpJetsByAscPseudorapidity(const Jet& a, const Jet& b) {
00220 return a.momentum().pseudorapidity() < b.momentum().pseudorapidity();
00221 }
00222
00223
00224
00225 inline bool cmpJetsByDescAbsPseudorapidity(const Jet& a, const Jet& b) {
00226 return fabs(a.momentum().pseudorapidity()) > fabs(b.momentum().pseudorapidity());
00227 }
00228
00229
00230 inline bool cmpJetsByAscAbsPseudorapidity(const Jet& a, const Jet& b) {
00231 return fabs(a.momentum().pseudorapidity()) < fabs(b.momentum().pseudorapidity());
00232 }
00233
00234
00235
00236 inline bool cmpJetsByDescRapidity(const Jet& a, const Jet& b) {
00237 return a.momentum().rapidity() > b.momentum().rapidity();
00238 }
00239
00240
00241 inline bool cmpJetsByAscRapidity(const Jet& a, const Jet& b) {
00242 return a.momentum().rapidity() < b.momentum().rapidity();
00243 }
00244
00245
00246
00247 inline bool cmpJetsByDescAbsRapidity(const Jet& a, const Jet& b) {
00248 return fabs(a.momentum().rapidity()) > fabs(b.momentum().rapidity());
00249 }
00250
00251
00252 inline bool cmpJetsByAscAbsRapidity(const Jet& a, const Jet& b) {
00253 return fabs(a.momentum().rapidity()) < fabs(b.momentum().rapidity());
00254 }
00255
00256
00257
00258
00259 }
00260
00261 #endif