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/ClusterSequence.hh"
00012 #include "fastjet/PseudoJet.hh"
00013 
00014 namespace Rivet {
00015 
00016 
00017   /// Make a 3-momentum vector from a FastJet pseudo-jet
00018   inline Vector3 momentum3(const fastjet::PseudoJet& pj) {
00019     return Vector3(pj.px(), pj.py(), pj.pz());
00020   }
00021 
00022   /// Make a 4-momentum vector from a FastJet pseudo-jet
00023   inline FourMomentum momentum(const fastjet::PseudoJet& pj) {
00024     return FourMomentum(pj.E(), pj.px(), pj.py(), pj.pz());
00025   }
00026 
00027 
00028   /// Typedef for a collection of PseudoJets.
00029   typedef vector<fastjet::PseudoJet> PseudoJets;
00030 
00031 
00032   /// Project out jets found using fastJet package.
00033   class FastJets : public JetAlg {
00034 
00035   public:
00036     /// Wrapper enum for selected Fastjet jet algorithms.
00037     enum JetAlgName { KT, CAM, SISCONE, ANTIKT, PXCONE,
00038                   CDFJETCLU, CDFMIDPOINT, D0ILCONE,
00039                   JADE, DURHAM, TRACKJET };
00040 
00041 
00042     /// @name Constructors etc.
00043     //@{
00044 
00045     /// "Wrapped" argument constructor using Rivet enums for most common
00046     /// jet alg choices (including some plugins). For the built-in algs,
00047     /// E-scheme recombination is used. For full control of
00048     /// FastJet built-in jet algs, use the native arg constructor.
00049     FastJets(const FinalState& fsp, JetAlgName alg,
00050              double rparameter, double seed_threshold=1.0);
00051 
00052     /// Native argument constructor, using FastJet alg/scheme enums.
00053     FastJets(const FinalState& fsp, fastjet::JetAlgorithm type,
00054              fastjet::RecombinationScheme recom, double rparameter);
00055 
00056     /// Explicitly pass in an externally-constructed plugin
00057     FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin& plugin);
00058 
00059     // /// Explicit copy constructor.
00060     // FastJets(const FastJets& other);
00061 
00062     /// Clone on the heap.
00063     virtual const Projection* clone() const {
00064       return new FastJets(*this);
00065     }
00066 
00067     //@}
00068 
00069 
00070   public:
00071 
00072     /// Reset the projection. Jet def etc are unchanged.
00073     void reset();
00074 
00075     /// Number of jets above the \f$ p_\perp \f$ cut.
00076     size_t numJets(double ptmin = 0.0) const;
00077 
00078     /// Number of jets.
00079     size_t size() const {
00080       return numJets();
00081     }
00082 
00083     /// Get the jets (unordered).
00084     Jets jets(double ptmin = 0.0) const;
00085 
00086     /// Get the jets, ordered by \f$ p_T \f$.
00087     Jets jetsByPt(double ptmin = 0.0) const;
00088 
00089     /// Get the jets, ordered by \f$ E \f$.
00090     Jets jetsByE(double ptmin = 0.0) const;
00091 
00092     /// Get the jets, ordered by rapidity.
00093     Jets jetsByRapidity(double ptmin = 0.0) const;
00094 
00095     /// Get the pseudo jets (unordered).
00096     PseudoJets pseudoJets(double ptmin = 0.0) const;
00097 
00098     /// Get the pseudo jets, ordered by \f$ p_T \f$.
00099     PseudoJets pseudoJetsByPt(double ptmin = 0.0) const {
00100       return sorted_by_pt(pseudoJets(ptmin));
00101     }
00102 
00103     /// Get the pseudo jets, ordered by \f$ E \f$.
00104     PseudoJets pseudoJetsByE(double ptmin = 0.0) const {
00105       return sorted_by_E(pseudoJets(ptmin));
00106     }
00107 
00108     /// Get the pseudo jets, ordered by rapidity.
00109     PseudoJets pseudoJetsByRapidity(double ptmin = 0.0) const {
00110       return sorted_by_rapidity(pseudoJets(ptmin));
00111     }
00112 
00113     /// Return the cluster sequence (FastJet-specific).
00114     const fastjet::ClusterSequence* clusterSeq() const {
00115       return _cseq.get();
00116     }
00117 
00118     /// Return the jet definition (FastJet-specific).
00119     const fastjet::JetDefinition& jetDef() const {
00120       return _jdef;
00121     }
00122 
00123     /// Get the subjet splitting variables for the given jet.
00124     vector<double> ySubJet(const fastjet::PseudoJet& jet) const;
00125 
00126     /// @brief Split a jet a la PRL100,242001(2008).
00127     /// Based on code from G.Salam, A.Davison.
00128     fastjet::PseudoJet splitJet(fastjet::PseudoJet jet, double& last_R) const;
00129 
00130     /// @brief Filter a jet a la PRL100,242001(2008).
00131     /// Based on code from G.Salam, A.Davison.
00132     fastjet::PseudoJet filterJet(fastjet::PseudoJet jet, double& stingy_R, const double def_R) const;
00133 
00134   private:
00135 
00136     Jets _pseudojetsToJets(const PseudoJets& pjets) const;
00137 
00138   protected:
00139 
00140     /// Perform the projection on the Event.
00141     void project(const Event& e);
00142 
00143     /// Compare projections.
00144     int compare(const Projection& p) const;
00145 
00146   public:
00147 
00148     /// Do the calculation locally (no caching).
00149     void calc(const ParticleVector& ps);
00150 
00151 
00152   private:
00153 
00154     /// Jet definition
00155     fastjet::JetDefinition _jdef;
00156 
00157     /// Cluster sequence
00158     shared_ptr<fastjet::ClusterSequence> _cseq;
00159 
00160     /// FastJet external plugin
00161     shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
00162 
00163     /// Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation.
00164     mutable map<int, vector<double> > _yscales;
00165 
00166     /// set of particles sorted by their PT2
00167     //set<Particle, ParticleBase::byPTAscending> _particles;
00168     map<int, Particle> _particles;
00169 
00170   };
00171 
00172 }
00173 
00174 #endif