00001
00002 #ifndef RIVET_FastJets_HH
00003 #define RIVET_FastJets_HH
00004
00005 #include "Rivet/Projection.hh"
00006 #include "Rivet/Projections/JetAlg.hh"
00007 #include "Rivet/Projections/FinalState.hh"
00008 #include "Rivet/Particle.hh"
00009 #include "Rivet/Jet.hh"
00010
00011 #include "fastjet/JetDefinition.hh"
00012 #include "fastjet/AreaDefinition.hh"
00013 #include "fastjet/ClusterSequence.hh"
00014 #include "fastjet/ClusterSequenceArea.hh"
00015 #include "fastjet/PseudoJet.hh"
00016
00017 namespace Rivet {
00018
00019
00020
00021 inline Vector3 momentum3(const fastjet::PseudoJet& pj) {
00022 return Vector3(pj.px(), pj.py(), pj.pz());
00023 }
00024
00025
00026 inline FourMomentum momentum(const fastjet::PseudoJet& pj) {
00027 return FourMomentum(pj.E(), pj.px(), pj.py(), pj.pz());
00028 }
00029
00030
00031
00032 typedef vector<fastjet::PseudoJet> PseudoJets;
00033
00034
00035
00036
00037
00038
00039
00040
00041 class FastJets : public JetAlg {
00042
00043 public:
00044
00045 enum JetAlgName { KT, CAM, SISCONE, ANTIKT, PXCONE,
00046 CDFJETCLU, CDFMIDPOINT, D0ILCONE,
00047 JADE, DURHAM, TRACKJET };
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057 FastJets(const FinalState& fsp, JetAlgName alg,
00058 double rparameter, double seed_threshold=1.0);
00059
00060
00061 FastJets(const FinalState& fsp, fastjet::JetAlgorithm type,
00062 fastjet::RecombinationScheme recom, double rparameter);
00063
00064
00065 FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin& plugin);
00066
00067
00068
00069
00070
00071 virtual const Projection* clone() const {
00072 return new FastJets(*this);
00073 }
00074
00075
00076
00077
00078 public:
00079
00080
00081 void reset();
00082
00083
00084
00085
00086 void useJetArea(fastjet::AreaDefinition* adef) {
00087 _adef = adef;
00088 }
00089
00090
00091 size_t numJets(double ptmin = 0.0) const;
00092
00093
00094 size_t size() const {
00095 return numJets();
00096 }
00097
00098
00099 Jets jets(double ptmin = 0.0) const;
00100
00101
00102 Jets jetsByPt(double ptmin = 0.0) const;
00103
00104
00105 Jets jetsByE(double ptmin = 0.0) const;
00106
00107
00108 Jets jetsByRapidity(double ptmin = 0.0) const;
00109
00110
00111 PseudoJets pseudoJets(double ptmin = 0.0) const;
00112
00113
00114 PseudoJets pseudoJetsByPt(double ptmin = 0.0) const {
00115 return sorted_by_pt(pseudoJets(ptmin));
00116 }
00117
00118
00119 PseudoJets pseudoJetsByE(double ptmin = 0.0) const {
00120 return sorted_by_E(pseudoJets(ptmin));
00121 }
00122
00123
00124 PseudoJets pseudoJetsByRapidity(double ptmin = 0.0) const {
00125 return sorted_by_rapidity(pseudoJets(ptmin));
00126 }
00127
00128
00129 const fastjet::ClusterSequence* clusterSeq() const {
00130 return _cseq.get();
00131 }
00132
00133
00134 const fastjet::ClusterSequenceArea* clusterSeqArea() const {
00135
00136 if (_adef == 0) return (fastjet::ClusterSequenceArea*) 0;
00137 return dynamic_cast<fastjet::ClusterSequenceArea*>(_cseq.get());
00138 }
00139
00140
00141 const fastjet::JetDefinition& jetDef() const {
00142 return _jdef;
00143 }
00144
00145
00146 const fastjet::AreaDefinition* areaDef() const {
00147 return _adef;
00148 }
00149
00150
00151 vector<double> ySubJet(const fastjet::PseudoJet& jet) const;
00152
00153
00154
00155 fastjet::PseudoJet splitJet(fastjet::PseudoJet jet, double& last_R) const;
00156
00157
00158
00159 fastjet::PseudoJet filterJet(fastjet::PseudoJet jet, double& stingy_R, const double def_R) const;
00160
00161 private:
00162
00163 Jets _pseudojetsToJets(const PseudoJets& pjets) const;
00164
00165 protected:
00166
00167
00168 void project(const Event& e);
00169
00170
00171 int compare(const Projection& p) const;
00172
00173 public:
00174
00175
00176 void calc(const ParticleVector& ps);
00177
00178
00179 private:
00180
00181
00182 fastjet::JetDefinition _jdef;
00183
00184
00185 fastjet::AreaDefinition* _adef;
00186
00187
00188 shared_ptr<fastjet::ClusterSequence> _cseq;
00189
00190
00191 shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
00192
00193
00194 mutable map<int, vector<double> > _yscales;
00195
00196
00197
00198 map<int, Particle> _particles;
00199
00200 };
00201
00202 }
00203
00204 #endif