FastJets.hh
Go to the documentation of this file.
00001 // -*- C++ -*- 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 #include "fastjet/SISConePlugin.hh" 00018 #include "fastjet/ATLASConePlugin.hh" 00019 #include "fastjet/CMSIterativeConePlugin.hh" 00020 #include "fastjet/CDFJetCluPlugin.hh" 00021 #include "fastjet/CDFMidPointPlugin.hh" 00022 #include "fastjet/D0RunIIConePlugin.hh" 00023 #include "fastjet/TrackJetPlugin.hh" 00024 #include "fastjet/JadePlugin.hh" 00025 //#include "fastjet/PxConePlugin.hh" 00026 00027 namespace Rivet { 00028 00029 00030 /// Make a 3-momentum vector from a FastJet pseudojet 00031 inline Vector3 momentum3(const fastjet::PseudoJet& pj) { 00032 return Vector3(pj.px(), pj.py(), pj.pz()); 00033 } 00034 00035 /// Make a 4-momentum vector from a FastJet pseudojet 00036 inline FourMomentum momentum(const fastjet::PseudoJet& pj) { 00037 return FourMomentum(pj.E(), pj.px(), pj.py(), pj.pz()); 00038 } 00039 00040 00041 00042 ///////////////////////// 00043 00044 00045 00046 /// Project out jets found using the FastJet package jet algorithms. 00047 class FastJets : public JetAlg { 00048 public: 00049 00050 /// Wrapper enum for selected Fastjet jet algorithms. 00051 enum JetAlgName { KT, CAM, SISCONE, ANTIKT, 00052 PXCONE, 00053 ATLASCONE, CMSCONE, 00054 CDFJETCLU, CDFMIDPOINT, D0ILCONE, 00055 JADE, DURHAM, TRACKJET }; 00056 00057 00058 /// @name Constructors etc. 00059 //@{ 00060 00061 /// "Wrapped" argument constructor using Rivet enums for most common 00062 /// jet alg choices (including some plugins). For the built-in algs, 00063 /// E-scheme recombination is used. For full control of 00064 /// FastJet built-in jet algs, use the native arg constructor. 00065 FastJets(const FinalState& fsp, JetAlgName alg, 00066 double rparameter, double seed_threshold=1.0) 00067 : JetAlg(fsp), _adef(0) { _init1(alg, rparameter, seed_threshold); } 00068 00069 /// Native argument constructor, using FastJet alg/scheme enums. 00070 FastJets(const FinalState& fsp, fastjet::JetAlgorithm type, 00071 fastjet::RecombinationScheme recom, double rparameter) 00072 : JetAlg(fsp), _adef(0) { _init2(type, recom, rparameter); } 00073 00074 /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) 00075 FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin* plugin) 00076 : JetAlg(fsp), _adef(0) { _init3(plugin); } 00077 /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) 00078 FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin& plugin) 00079 : JetAlg(fsp), _adef(0) { _init3(&plugin); } 00080 00081 00082 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00083 FastJets(JetAlgName alg, double rparameter, double seed_threshold=1.0) 00084 : _adef(0) { _init1(alg, rparameter, seed_threshold); } 00085 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00086 FastJets(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter) 00087 : _adef(0) { _init2(type, recom, rparameter); } 00088 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00089 FastJets(fastjet::JetDefinition::Plugin* plugin) 00090 : _adef(0) { _init3(plugin); } 00091 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00092 FastJets(fastjet::JetDefinition::Plugin& plugin) 00093 : _adef(0) { _init3(&plugin); } 00094 00095 00096 /// Clone on the heap. 00097 virtual const Projection* clone() const { 00098 return new FastJets(*this); 00099 } 00100 00101 //@} 00102 00103 00104 public: 00105 00106 /// Reset the projection. Jet def, etc. are unchanged. 00107 void reset(); 00108 00109 /// @brief Use provided jet area definition 00110 /// 00111 /// @warning @a adef is NOT copied, the user must ensure that it remains valid! 00112 /// 00113 /// Provide an adef null pointer to re-disable jet area calculation 00114 void useJetArea(fastjet::AreaDefinition* adef) { 00115 _adef = adef; 00116 } 00117 00118 /// @name Access to the jets 00119 //@{ 00120 00121 /// Number of jets above the \f$ p_\perp \f$ cut. 00122 size_t numJets(double ptmin = 0.0) const; 00123 00124 /// Number of jets. 00125 size_t size() const { 00126 return numJets(); 00127 } 00128 00129 /// Get the jets (unordered). 00130 Jets _jets(double ptmin = 0.0) const; 00131 00132 /// Get the pseudo jets (unordered). 00133 PseudoJets pseudoJets(double ptmin = 0.0) const; 00134 /// Alias 00135 PseudoJets pseudojets(double ptmin = 0.0) const { return pseudoJets(ptmin); } 00136 00137 /// Get the pseudo jets, ordered by \f$ p_T \f$. 00138 PseudoJets pseudoJetsByPt(double ptmin = 0.0) const { 00139 return sorted_by_pt(pseudoJets(ptmin)); 00140 } 00141 /// Alias 00142 PseudoJets pseudojetsByPt(double ptmin = 0.0) const { return pseudoJetsByPt(ptmin); } 00143 00144 /// Get the pseudo jets, ordered by \f$ E \f$. 00145 PseudoJets pseudoJetsByE(double ptmin = 0.0) const { 00146 return sorted_by_E(pseudoJets(ptmin)); 00147 } 00148 /// Alias 00149 PseudoJets pseudojetsByE(double ptmin = 0.0) const { return pseudoJetsByE(ptmin); } 00150 00151 /// Get the pseudo jets, ordered by rapidity. 00152 PseudoJets pseudoJetsByRapidity(double ptmin = 0.0) const { 00153 return sorted_by_rapidity(pseudoJets(ptmin)); 00154 } 00155 /// Alias 00156 PseudoJets pseudojetsByRapidity(double ptmin = 0.0) const { return pseudoJetsByRapidity(ptmin); } 00157 00158 //@} 00159 00160 00161 /// @name Access to the cluster sequence and splitting info 00162 //@{ 00163 00164 /// Return the cluster sequence. 00165 const fastjet::ClusterSequence* clusterSeq() const { 00166 return _cseq.get(); 00167 } 00168 00169 /// Return the area-enabled cluster sequence (if an area defn exists, otherwise returns a null ptr). 00170 const fastjet::ClusterSequenceArea* clusterSeqArea() const { 00171 if (areaDef() == NULL) return (fastjet::ClusterSequenceArea*) 0; 00172 return dynamic_cast<fastjet::ClusterSequenceArea*>(_cseq.get()); 00173 } 00174 00175 /// Return the jet definition. 00176 const fastjet::JetDefinition& jetDef() const { 00177 return _jdef; 00178 } 00179 00180 /// Return the area definition.. May be null. 00181 const fastjet::AreaDefinition* areaDef() const { 00182 return _adef; 00183 } 00184 00185 00186 /// Get the subjet splitting variables for the given jet. 00187 vector<double> ySubJet(const fastjet::PseudoJet& jet) const; 00188 00189 /// @brief Split a jet a la PRL100,242001(2008). 00190 /// 00191 /// Based on code from G.Salam, A.Davison. 00192 fastjet::PseudoJet splitJet(fastjet::PseudoJet jet, double& last_R) const; 00193 00194 /// @brief Filter a jet a la PRL100,242001(2008). 00195 /// 00196 /// Based on code from G.Salam, A.Davison. 00197 fastjet::PseudoJet filterJet(fastjet::PseudoJet jet, double& stingy_R, const double def_R) const; 00198 00199 //@} 00200 00201 00202 private: 00203 00204 /// Shared utility functions to implement constructor behaviour 00205 void _init1(JetAlgName alg, double rparameter, double seed_threshold); 00206 void _init2(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter); 00207 void _init3(fastjet::JetDefinition::Plugin* plugin); 00208 00209 protected: 00210 00211 /// Perform the projection on the Event. 00212 void project(const Event& e); 00213 00214 /// Compare projections. 00215 int compare(const Projection& p) const; 00216 00217 public: 00218 00219 /// Do the calculation locally (no caching). 00220 void calc(const Particles& fsparticles, const Particles& tagparticles=Particles()); 00221 00222 00223 private: 00224 00225 /// Jet definition 00226 fastjet::JetDefinition _jdef; 00227 00228 /// Pointer to user-handled area definition 00229 fastjet::AreaDefinition* _adef; 00230 00231 /// Cluster sequence 00232 shared_ptr<fastjet::ClusterSequence> _cseq; 00233 00234 /// FastJet external plugin 00235 shared_ptr<fastjet::JetDefinition::Plugin> _plugin; 00236 00237 /// Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation. 00238 mutable map<int, vector<double> > _yscales; 00239 00240 /// set of particles sorted by their PT2 00241 //set<Particle, ParticleBase::byPTAscending> _particles; 00242 map<int, Particle> _particles; 00243 00244 }; 00245 00246 } 00247 00248 #endif Generated on Tue Sep 30 2014 19:45:44 for The Rivet MC analysis system by ![]() |