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/Jet.hh" 00006 #include "Rivet/Particle.hh" 00007 #include "Rivet/Projection.hh" 00008 #include "Rivet/Projections/JetAlg.hh" 00009 #include "Rivet/Projections/FinalState.hh" 00010 #include "Rivet/Tools/RivetFastJet.hh" 00011 00012 #include "fastjet/SISConePlugin.hh" 00013 #include "fastjet/ATLASConePlugin.hh" 00014 #include "fastjet/CMSIterativeConePlugin.hh" 00015 #include "fastjet/CDFJetCluPlugin.hh" 00016 #include "fastjet/CDFMidPointPlugin.hh" 00017 #include "fastjet/D0RunIIConePlugin.hh" 00018 #include "fastjet/TrackJetPlugin.hh" 00019 #include "fastjet/JadePlugin.hh" 00020 //#include "fastjet/PxConePlugin.hh" 00021 00022 namespace Rivet { 00023 00024 00025 /// Project out jets found using the FastJet package jet algorithms. 00026 class FastJets : public JetAlg { 00027 public: 00028 00029 /// Wrapper enum for selected Fastjet jet algorithms. 00030 enum JetAlgName { KT, CAM, SISCONE, ANTIKT, 00031 // PXCONE, 00032 ATLASCONE, CMSCONE, 00033 CDFJETCLU, CDFMIDPOINT, D0ILCONE, 00034 JADE, DURHAM, TRACKJET }; 00035 00036 00037 /// @name Constructors etc. 00038 //@{ 00039 00040 /// "Wrapped" argument constructor using Rivet enums for most common 00041 /// jet alg choices (including some plugins). For the built-in algs, 00042 /// E-scheme recombination is used. For full control of 00043 /// FastJet built-in jet algs, use the native arg constructor. 00044 FastJets(const FinalState& fsp, JetAlgName alg, double rparameter, 00045 JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, 00046 JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES, 00047 double seed_threshold=1.0) 00048 : JetAlg(fsp, usemuons, useinvis), _adef(0) { _init1(alg, rparameter, seed_threshold); } 00049 00050 /// Native argument constructor, using FastJet alg/scheme enums. 00051 FastJets(const FinalState& fsp, fastjet::JetAlgorithm type, 00052 fastjet::RecombinationScheme recom, double rparameter, 00053 JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, 00054 JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES) 00055 : JetAlg(fsp, usemuons, useinvis), _adef(0) { _init2(type, recom, rparameter); } 00056 00057 /// Explicitly pass in a JetDefinition 00058 FastJets(const FinalState& fsp, const fastjet::JetDefinition& jdef, 00059 JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, 00060 JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES) 00061 : JetAlg(fsp, usemuons, useinvis), _adef(0) { _init3(jdef); } 00062 00063 /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) 00064 FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin* plugin, 00065 JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, 00066 JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES) 00067 : JetAlg(fsp, usemuons, useinvis), _adef(0) { _init4(plugin); } 00068 00069 /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) 00070 /// @deprecated Use the pointer version -- it makes the lifetime & ownership more obvious 00071 FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin& plugin, 00072 JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, 00073 JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES) 00074 : JetAlg(fsp, usemuons, useinvis), _adef(0) { _init4(&plugin); } 00075 00076 00077 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00078 FastJets(JetAlgName alg, double rparameter, double seed_threshold=1.0) 00079 : _adef(0) { _init1(alg, rparameter, seed_threshold); } 00080 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00081 FastJets(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter) 00082 : _adef(0) { _init2(type, recom, rparameter); } 00083 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00084 FastJets(fastjet::JetDefinition::Plugin* plugin) 00085 : _adef(0) { _init3(plugin); } 00086 /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) 00087 FastJets(fastjet::JetDefinition::Plugin& plugin) 00088 : _adef(0) { _init3(&plugin); } 00089 00090 00091 /// Clone on the heap. 00092 virtual const Projection* clone() const { 00093 return new FastJets(*this); 00094 } 00095 00096 //@} 00097 00098 00099 public: 00100 00101 /// Reset the projection. Jet def, etc. are unchanged. 00102 void reset(); 00103 00104 /// @brief Use provided jet area definition 00105 /// 00106 /// @warning @a adef is NOT copied, the user must ensure that it remains valid! 00107 /// 00108 /// Provide an adef null pointer to re-disable jet area calculation 00109 void useJetArea(fastjet::AreaDefinition* adef) { 00110 _adef = adef; 00111 } 00112 00113 /// @name Access to the jets 00114 //@{ 00115 00116 /// Number of jets above the \f$ p_\perp \f$ cut. 00117 size_t numJets(double ptmin=0.0) const; 00118 00119 /// Number of jets. 00120 size_t size() const { 00121 return numJets(); 00122 } 00123 00124 /// Get the jets (unordered). 00125 Jets _jets(double ptmin=0.0) const; 00126 00127 /// Get the pseudo jets (unordered). 00128 PseudoJets pseudoJets(double ptmin=0.0) const; 00129 /// Alias 00130 PseudoJets pseudojets(double ptmin=0.0) const { return pseudoJets(ptmin); } 00131 00132 /// Get the pseudo jets, ordered by \f$ p_T \f$. 00133 PseudoJets pseudoJetsByPt(double ptmin=0.0) const { 00134 return sorted_by_pt(pseudoJets(ptmin)); 00135 } 00136 /// Alias 00137 PseudoJets pseudojetsByPt(double ptmin=0.0) const { return pseudoJetsByPt(ptmin); } 00138 00139 /// Get the pseudo jets, ordered by \f$ E \f$. 00140 PseudoJets pseudoJetsByE(double ptmin=0.0) const { 00141 return sorted_by_E(pseudoJets(ptmin)); 00142 } 00143 /// Alias 00144 PseudoJets pseudojetsByE(double ptmin=0.0) const { return pseudoJetsByE(ptmin); } 00145 00146 /// Get the pseudo jets, ordered by rapidity. 00147 PseudoJets pseudoJetsByRapidity(double ptmin=0.0) const { 00148 return sorted_by_rapidity(pseudoJets(ptmin)); 00149 } 00150 /// Alias 00151 PseudoJets pseudojetsByRapidity(double ptmin=0.0) const { return pseudoJetsByRapidity(ptmin); } 00152 00153 //@} 00154 00155 00156 /// @name Access to the FastJet clustering objects such as jet def, area def, and cluster 00157 //@{ 00158 00159 /// Return the cluster sequence. 00160 const fastjet::ClusterSequence* clusterSeq() const { 00161 return _cseq.get(); 00162 } 00163 00164 /// Return the area-enabled cluster sequence (if an area defn exists, otherwise returns a null ptr). 00165 const fastjet::ClusterSequenceArea* clusterSeqArea() const { 00166 if (areaDef() == NULL) return (fastjet::ClusterSequenceArea*) 0; 00167 return dynamic_cast<fastjet::ClusterSequenceArea*>(_cseq.get()); 00168 } 00169 00170 /// Return the jet definition. 00171 const fastjet::JetDefinition& jetDef() const { 00172 return _jdef; 00173 } 00174 00175 /// Return the area definition.. May be null. 00176 const fastjet::AreaDefinition* areaDef() const { 00177 return _adef; 00178 } 00179 00180 //@} 00181 00182 00183 /// @name Access to jet splitting variables: DISABLED FROM 2.3.0, USE FASTJET OBJECTS DIRECTLY INSTEAD 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 /// @todo Replace with calls between constructors when C++11 available? 00206 void _initBase(); 00207 void _init1(JetAlgName alg, double rparameter, double seed_threshold); 00208 void _init2(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter); 00209 void _init3(const fastjet::JetDefinition& plugin); 00210 void _init4(fastjet::JetDefinition::Plugin* plugin); 00211 00212 protected: 00213 00214 /// Perform the projection on the Event. 00215 void project(const Event& e); 00216 00217 /// Compare projections. 00218 int compare(const Projection& p) const; 00219 00220 public: 00221 00222 /// Do the calculation locally (no caching). 00223 void calc(const Particles& fsparticles, const Particles& tagparticles=Particles()); 00224 00225 00226 private: 00227 00228 /// Jet definition 00229 fastjet::JetDefinition _jdef; 00230 00231 /// Pointer to user-handled area definition 00232 fastjet::AreaDefinition* _adef; 00233 00234 /// Cluster sequence 00235 shared_ptr<fastjet::ClusterSequence> _cseq; 00236 00237 /// FastJet external plugin 00238 shared_ptr<fastjet::JetDefinition::Plugin> _plugin; 00239 00240 /// Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation. 00241 mutable map<int, vector<double> > _yscales; 00242 00243 /// set of particles sorted by their PT2 00244 //set<Particle, ParticleBase::byPTAscending> _particles; 00245 map<int, Particle> _particles; 00246 00247 }; 00248 00249 } 00250 00251 #endif Generated on Thu Mar 10 2016 08:29:50 for The Rivet MC analysis system by ![]() |