Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-30 10:06:44

0001 //  VariableR Package
0002 //  Questions/Comments?  jthaler@jthaler.net
0003 //
0004 //  Copyright (c) 2009-2016
0005 //  David Krohn, Gregory Soyez, Jesse Thaler, and Lian-Tao Wang
0006 //
0007 //  $Id$
0008 //----------------------------------------------------------------------
0009 // This file is part of FastJet contrib.
0010 //
0011 // It is free software; you can redistribute it and/or modify it under
0012 // the terms of the GNU General Public License as published by the
0013 // Free Software Foundation; either version 2 of the License, or (at
0014 // your option) any later version.
0015 //
0016 // It is distributed in the hope that it will be useful, but WITHOUT
0017 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
0018 // or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public
0019 // License for more details.
0020 //
0021 // You should have received a copy of the GNU General Public License
0022 // along with this code. If not, see <http://www.gnu.org/licenses/>.
0023 //----------------------------------------------------------------------
0024 
0025 #ifndef __FASTJET_CONTRIB_VARIABLERPLUGIN_HH__
0026 #define __FASTJET_CONTRIB_VARIABLERPLUGIN_HH__
0027 
0028 #include <fastjet/internal/base.hh>
0029 #include <fastjet/config.h>
0030 
0031 #include "fastjet/JetDefinition.hh"
0032 #include "fastjet/PseudoJet.hh"
0033 #include "fastjet/ClusterSequence.hh"
0034 #include <fastjet/LimitedWarning.hh>
0035 
0036 #include <map>
0037 #include <queue>
0038 
0039 using namespace std;
0040 
0041 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0042 
0043 namespace contrib {
0044   
0045   ////////
0046   //
0047   //  Core VR Code
0048   //
0049   ////////
0050    
0051   // In version 1.1, this replaces CoreJetAlgorithm.
0052   // This acts like any fastjet plugin since it implements run_clustering
0053   class VariableRPlugin : public JetDefinition::Plugin {
0054 
0055   public:
0056     /// Type of clustering
0057     ///
0058     /// Since version 1.2.0 of VariableR, the clustering is treated as
0059     /// a generalised-kt algorithm and the previous "ClusterType"
0060     /// argument is replaced by the "p" parameter of the
0061     /// generalised-kt algorithm. The definitions below are shorthand
0062     /// for the antikt, C/A and kt algorithm which also allow for
0063     /// backwards compatibility.
0064     /// (These are initialized in VariableRPlugin.cc.)
0065     static const double CALIKE;  //  =  0.0;
0066     static const double KTLIKE;  //  =  1.0;
0067     static const double AKTLIKE; //  = -1.0;
0068 
0069     /// for backwards compatibility reasons, we also define ClusterType
0070     /// as "double"
0071     typedef double ClusterType;
0072 
0073     /// The strategy to be used with the clustering
0074     enum Strategy{
0075       Best,      ///< currently N2Tiled or N2Plain for FJ>3.2.0, Native for FastJet<3.2.0
0076       N2Tiled,   ///< the default (faster in most cases) [requires FastJet>=3.2.0]
0077       N2Plain,   ///< [requires FastJet>=3.2.0]
0078       NNH,       ///< slower but already available for FastJet<3.2.0
0079       Native     ///< original local implemtation of the clustering [the default for FastJet<3.2.0]
0080     };
0081 
0082     
0083     /// Constructor that sets VR algorithm parameters
0084     ///
0085     ///  \param rho         mass scale for effective radius (i.e. R ~ rho/pT)
0086     ///  \param min_r       minimum jet radius
0087     ///  \param max_r       maximum jet radius
0088     ///  \param clust_type  whether to use CA-like, kT-like, or anti-kT-like distance measure
0089     ///                     (this value is the same as the p exponent in generalized-kt, with
0090     ///                     anti-kt = -1.0, CA = 0.0, and kT = 1.0)
0091     ///  \param precluster  whether to use optional kT subjets (of size min_r) for preclustering
0092     ///                     (true is much faster, default=false). At the moment, the only option
0093     ///                     for preclustering is kT (use fastjet::NestedDefsPjugin otherwise)
0094     ///  \param strategy    decodes which algorithm to apply for the clustering
0095     ///
0096     /// Note that pre-clustering is deprecated and will likely be
0097     /// removed in a future releasse of this contrib. You can get the
0098     /// same behaviour using the NestedDefs fastjet plugin:
0099     ///
0100     /// \code
0101     ///    std::list<JetDefinition> jet_defs;
0102     ///    jet_defs.push_back(JetDefinition(kt_algorithm, min_r));
0103     ///    contrib::VariableRPlugin vr_plugin(rho,min_r,max_r,clust_type);
0104     ///    jet_defs.push_back(JetDefinition(&vr_plugin));
0105     ///
0106     ///    NestedDefsPlugin nd_plugin(jet_defs);
0107     ///    JetDefinition jet_def(&nd_plugin);
0108     /// \endcode
0109     ///
0110     /// For FastJet>=3.2.0, we will use by default the N2Tiled
0111     /// strategy unless pre-clustering is requested. Otherwise, we use
0112     /// the "old" native strategy.
0113     VariableRPlugin(double rho, double min_r, double max_r, ClusterType clust_type, bool precluster = false,
0114                     Strategy requested_strategy = Best);
0115       
0116     // virtual function from JetDefinition::Plugin that implements the actual VR algorithm
0117     void run_clustering(ClusterSequence & cs) const;
0118       
0119     // information string
0120     virtual string description() const;
0121       
0122     // TODO:  have this return a non-trivial answer.
0123     virtual double R() const{ return _max_r; }
0124       
0125   private:
0126       
0127     // parameters to define VR algorithm
0128     double _rho2, _min_r2, _max_r, _max_r2;
0129     ClusterType _clust_type;
0130     Strategy _requested_strategy;
0131       
0132     // For preclustering, can use kT algorithm to make subclusters of size min_r
0133     bool _precluster;
0134     JetDefinition _pre_jet_def;
0135 
0136     // warn about pre-clustering being deprecated
0137     static LimitedWarning _preclustering_deprecated_warning;
0138     
0139     // helper function to decide what strategy is best
0140     //
0141     // The optimal strategy will depend on the multiplicity and _max_r
0142     Strategy _best_strategy(unsigned int N) const;
0143     
0144     // helper function to apply kT preclustering if desired
0145     void _preclustering(ClusterSequence & cs, set<int>& unmerged_jets) const;
0146 
0147     // native implementation of the clustering
0148     void _native_clustering(ClusterSequence &cs) const;
0149 
0150     // implementation of the clustering using FastJet NN*** classes
0151     template<typename NN>
0152     void _NN_clustering(ClusterSequence &cs, NN &nn) const;
0153 
0154     // Helper struct to store two jets and a distance measure
0155     struct JetDistancePair{
0156       int j1,j2;
0157       double distance;
0158     };
0159       
0160     // Helper comparitor class for comparing JetDistancePairs
0161     class CompareJetDistancePair {
0162     public:
0163       CompareJetDistancePair(){};
0164       bool operator() (const JetDistancePair & lhs, const JetDistancePair &rhs) const {
0165         return (lhs.distance > rhs.distance);
0166       }
0167     };
0168       
0169     // helper function to merge jet with beam
0170     void _merge_jet_with_beam(ClusterSequence & clust_seq, JetDistancePair & jdp, set<int>& unmerged_jets) const;
0171       
0172     // helper function to merge two jets into a new pseudojet
0173     void _merge_jets(ClusterSequence & clust_seq,
0174                      JetDistancePair & jdp,
0175                      priority_queue< JetDistancePair, vector<JetDistancePair>, CompareJetDistancePair > &jet_queue,
0176                      set<int>& unmerged_jets) const;
0177       
0178     // use ClusterMode to determine jet-jet and jet-beam distance
0179     inline double _get_JJ_distance_measure(const PseudoJet& j1, const PseudoJet& j2) const;
0180     inline double _get_JB_distance_measure(const PseudoJet& jet) const;
0181       
0182     // helper function to establish measures
0183     void _setup_distance_measures(ClusterSequence & clust_seq,
0184                                   vector<JetDistancePair> &jet_vec,
0185                                   set<int>& unmerged_jets) const;
0186 
0187   };
0188    
0189 } // namespace contrib
0190 
0191 FASTJET_END_NAMESPACE
0192 
0193 #endif  // __FASTJET_CONTRIB_VARIABLERPLUGIN_HH__