rivet is hosted by Hepforge, IPPP Durham
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,
00045              double rparameter, double seed_threshold=1.0)
00046       : JetAlg(fsp), _adef(0) { _init1(alg, rparameter, seed_threshold); }
00047 
00048     /// Native argument constructor, using FastJet alg/scheme enums.
00049     FastJets(const FinalState& fsp, fastjet::JetAlgorithm type,
00050              fastjet::RecombinationScheme recom, double rparameter)
00051       : JetAlg(fsp), _adef(0) { _init2(type, recom, rparameter); }
00052 
00053     /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete)
00054     FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin* plugin)
00055       : JetAlg(fsp), _adef(0) { _init3(plugin); }
00056     /// Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete)
00057     FastJets(const FinalState& fsp, fastjet::JetDefinition::Plugin& plugin)
00058       : JetAlg(fsp), _adef(0) { _init3(&plugin); }
00059 
00060 
00061     /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method)
00062     FastJets(JetAlgName alg, double rparameter, double seed_threshold=1.0)
00063       : _adef(0) { _init1(alg, rparameter, seed_threshold); }
00064     /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method)
00065     FastJets(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter)
00066       : _adef(0) { _init2(type, recom, rparameter); }
00067     /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method)
00068     FastJets(fastjet::JetDefinition::Plugin* plugin)
00069       : _adef(0) { _init3(plugin); }
00070     /// Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method)
00071     FastJets(fastjet::JetDefinition::Plugin& plugin)
00072       : _adef(0) { _init3(&plugin); }
00073 
00074 
00075     /// Clone on the heap.
00076     virtual const Projection* clone() const {
00077       return new FastJets(*this);
00078     }
00079 
00080     //@}
00081 
00082 
00083   public:
00084 
00085     /// Reset the projection. Jet def, etc. are unchanged.
00086     void reset();
00087 
00088     /// @brief Use provided jet area definition
00089     ///
00090     /// @warning @a adef is NOT copied, the user must ensure that it remains valid!
00091     ///
00092     /// Provide an adef null pointer to re-disable jet area calculation
00093     void useJetArea(fastjet::AreaDefinition* adef) {
00094       _adef = adef;
00095     }
00096 
00097     /// @name Access to the jets
00098     //@{
00099 
00100     /// Number of jets above the \f$ p_\perp \f$ cut.
00101     size_t numJets(double ptmin=0.0) const;
00102 
00103     /// Number of jets.
00104     size_t size() const {
00105       return numJets();
00106     }
00107 
00108     /// Get the jets (unordered).
00109     Jets _jets(double ptmin=0.0) const;
00110 
00111     /// Get the pseudo jets (unordered).
00112     PseudoJets pseudoJets(double ptmin=0.0) const;
00113     /// Alias
00114     PseudoJets pseudojets(double ptmin=0.0) const { return pseudoJets(ptmin); }
00115 
00116     /// Get the pseudo jets, ordered by \f$ p_T \f$.
00117     PseudoJets pseudoJetsByPt(double ptmin=0.0) const {
00118       return sorted_by_pt(pseudoJets(ptmin));
00119     }
00120     /// Alias
00121     PseudoJets pseudojetsByPt(double ptmin=0.0) const { return pseudoJetsByPt(ptmin); }
00122 
00123     /// Get the pseudo jets, ordered by \f$ E \f$.
00124     PseudoJets pseudoJetsByE(double ptmin=0.0) const {
00125       return sorted_by_E(pseudoJets(ptmin));
00126     }
00127     /// Alias
00128     PseudoJets pseudojetsByE(double ptmin=0.0) const { return pseudoJetsByE(ptmin); }
00129 
00130     /// Get the pseudo jets, ordered by rapidity.
00131     PseudoJets pseudoJetsByRapidity(double ptmin=0.0) const {
00132       return sorted_by_rapidity(pseudoJets(ptmin));
00133     }
00134     /// Alias
00135     PseudoJets pseudojetsByRapidity(double ptmin=0.0) const { return pseudoJetsByRapidity(ptmin); }
00136 
00137     //@}
00138 
00139 
00140     /// @name Access to the cluster sequence and splitting info
00141     //@{
00142 
00143     /// Return the cluster sequence.
00144     const fastjet::ClusterSequence* clusterSeq() const {
00145       return _cseq.get();
00146     }
00147 
00148     /// Return the area-enabled cluster sequence (if an area defn exists, otherwise returns a null ptr).
00149     const fastjet::ClusterSequenceArea* clusterSeqArea() const {
00150       if (areaDef() == NULL) return (fastjet::ClusterSequenceArea*) 0;
00151       return dynamic_cast<fastjet::ClusterSequenceArea*>(_cseq.get());
00152     }
00153 
00154     /// Return the jet definition.
00155     const fastjet::JetDefinition& jetDef() const {
00156       return _jdef;
00157     }
00158 
00159     /// Return the area definition.. May be null.
00160     const fastjet::AreaDefinition* areaDef() const {
00161       return _adef;
00162     }
00163 
00164 
00165     /// Get the subjet splitting variables for the given jet.
00166     vector<double> ySubJet(const fastjet::PseudoJet& jet) const;
00167 
00168     /// @brief Split a jet a la PRL100,242001(2008).
00169     ///
00170     /// Based on code from G.Salam, A.Davison.
00171     fastjet::PseudoJet splitJet(fastjet::PseudoJet jet, double& last_R) const;
00172 
00173     /// @brief Filter a jet a la PRL100,242001(2008).
00174     ///
00175     /// Based on code from G.Salam, A.Davison.
00176     fastjet::PseudoJet filterJet(fastjet::PseudoJet jet, double& stingy_R, const double def_R) const;
00177 
00178     //@}
00179 
00180 
00181   private:
00182 
00183     /// Shared utility functions to implement constructor behaviour
00184     void _init1(JetAlgName alg, double rparameter, double seed_threshold);
00185     void _init2(fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter);
00186     void _init3(fastjet::JetDefinition::Plugin* plugin);
00187 
00188   protected:
00189 
00190     /// Perform the projection on the Event.
00191     void project(const Event& e);
00192 
00193     /// Compare projections.
00194     int compare(const Projection& p) const;
00195 
00196   public:
00197 
00198     /// Do the calculation locally (no caching).
00199     void calc(const Particles& fsparticles, const Particles& tagparticles=Particles());
00200 
00201 
00202   private:
00203 
00204     /// Jet definition
00205     fastjet::JetDefinition _jdef;
00206 
00207     /// Pointer to user-handled area definition
00208     fastjet::AreaDefinition* _adef;
00209 
00210     /// Cluster sequence
00211     shared_ptr<fastjet::ClusterSequence> _cseq;
00212 
00213     /// FastJet external plugin
00214     shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
00215 
00216     /// Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation.
00217     mutable map<int, vector<double> > _yscales;
00218 
00219     /// set of particles sorted by their PT2
00220     //set<Particle, ParticleBase::byPTAscending> _particles;
00221     map<int, Particle> _particles;
00222 
00223   };
00224 
00225 }
00226 
00227 #endif