Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:13

0001 //  ClusteringVetoPlugin Package
0002 //  Questions/Comments? liew@hep-th.phys.s.u-tokyo.ac.jp
0003 //                      stoll@hep-th.phys.s.u-tokyo.ac.jp
0004 //
0005 //  Copyright (c) 2014-2015
0006 //  Seng Pei Liew, Martin Stoll
0007 //
0008 //  $Id: ClusteringVetoPlugin.hh 792 2015-05-04 03:42:26Z martinstoll $
0009 //----------------------------------------------------------------------
0010 // This file is part of FastJet contrib.
0011 //
0012 // It is free software; you can redistribute it and/or modify it under
0013 // the terms of the GNU General Public License as published by the
0014 // Free Software Foundation; either version 2 of the License, or (at
0015 // your option) any later version.
0016 //
0017 // It is distributed in the hope that it will be useful, but WITHOUT
0018 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
0019 // or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public
0020 // License for more details.
0021 //
0022 // You should have received a copy of the GNU General Public License
0023 // along with this code. If not, see <http://www.gnu.org/licenses/>.
0024 //----------------------------------------------------------------------
0025 
0026 #ifndef __FASTJET_CONTRIB_CLUSTERINGVETOPLUGIN_HH__
0027 #define __FASTJET_CONTRIB_CLUSTERINGVETOPLUGIN_HH__
0028 
0029 #include <fastjet/internal/base.hh>
0030 
0031 #include "fastjet/JetDefinition.hh"
0032 #include "fastjet/PseudoJet.hh"
0033 #include "fastjet/ClusterSequence.hh"
0034 #include <fastjet/LimitedWarning.hh>
0035 
0036 #include <queue>
0037 
0038 using namespace std;
0039 
0040 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0041 
0042 namespace contrib {
0043 
0044   //------------------------------------------------------------------------
0045   /// \class ClusteringVetoPlugin
0046   /// This class implements a terminating clustering veto in the sense that
0047   /// if a recombination step is vetoed, both jets do not further participate
0048   /// in jet clustering ("turn passive").
0049   ///
0050   /// The default veto function is the mass-jump veto:
0051   ///
0052   ///  - 2 veto parameters: mu, theta
0053   ///  - Maximum jet radius parameter: max_r
0054   ///  - Metric d_ab, d_aB: CA-like, kT-like, anti-kT-like
0055   ///
0056   /// Mass-jump veto function in a recombination step j_a+j_b -> j:
0057   ///
0058   ///  - If mass m(j) < mu, cluster.
0059   ///  - If mass m(j) > mu and theta*m(j) > max(m(j_a),m(j_b)), veto.
0060   ///  - If mass m(j) > mu but theta*m(j) < max(m(j_a),m(j_b)), no veto.
0061   ///    Check active-passive veto, i.e.
0062   ///     * find passive jet j_n with smallest d_an and d_an < d_nB
0063   ///     * if d_an < d_ab, check veto between j_a and j_n
0064   ///     * if a veto is called, j_a turns passive
0065   ///    Do the same for j_b.
0066   ///  - If no active-passive veto is called, cluster.
0067   ///
0068   /// A user-defined veto function can be passed, which takes the return values
0069   /// { CLUSTER, VETO, NOVETO }
0070   /// corresponding to the first three bullet points in the mass-jump function.
0071   class ClusteringVetoPlugin : public JetDefinition::Plugin {
0072 
0073   public:
0074     // Type of clustering
0075     enum ClusterType {
0076       CALIKE,
0077       KTLIKE,
0078       AKTLIKE
0079     };
0080 
0081     // Result of veto function
0082     enum VetoResult {
0083       CLUSTER,
0084       VETO,
0085       NOVETO
0086     };
0087       
0088     // Constructor
0089     // mu         : veto condition 1, m(i+j) > mu
0090     // theta      : veto condition 2, theta * m(i+j) > max( m_i, m_j )
0091     // max_r      : maximum jet radius
0092     // clust_type : whether to use CA-like, kT-like, or anti-kT-like distance measure
0093     ClusteringVetoPlugin(double mu, double theta, double max_r, ClusterType clust_type);
0094       
0095     // Virtual function from JetDefinition::Plugin that implements the algorithm
0096     void run_clustering(fastjet::ClusterSequence & cs) const;
0097       
0098     // Information string
0099     virtual string description() const;
0100       
0101     // NOTE: Required by JetDefinition::Plugin
0102     double R() const { return sqrt(_max_r2); }
0103 
0104     // Set user-defined mass-jump veto
0105     void set_veto_function(VetoResult (*f)(const PseudoJet& j1,
0106                        const PseudoJet& j2)) {
0107       _veto_function = f; }
0108       
0109   private:
0110 
0111     // Parameters of MJ clustering
0112     double _max_r2, _mu, _theta;
0113     ClusterType _clust_type;
0114 
0115     // pointer to CheckVeto_xxx function
0116     VetoResult (*_veto_function) (const PseudoJet& j1, const PseudoJet& j2);
0117 
0118   private:
0119       
0120     // Use ClusterMode to determine jet-jet and jet-beam distance
0121     inline double GetJJDistanceMeasure(const PseudoJet& j1, const PseudoJet& j2) const;
0122     inline double GetJBDistanceMeasure(const PseudoJet& jet) const;
0123 
0124     // Check veto condition
0125     VetoResult CheckVeto(const PseudoJet& j1, const PseudoJet& j2) const;
0126     // Pre-defined mass-jump veto
0127     VetoResult CheckVeto_MJ(const PseudoJet& j1, const PseudoJet& j2) const;
0128 
0129   };
0130 
0131 
0132   /// helper class to store some clustering parameters
0133   class ClusteringVetoJetInfo {
0134   public:
0135     ClusteringVetoPlugin::ClusterType clust_type;
0136     double max_r2;
0137   };
0138   /// helper class for NNH
0139   class ClusteringVetoJet {
0140   public:
0141     void init(const PseudoJet & jet,
0142           const ClusteringVetoJetInfo * info) {
0143       ph = jet.phi();
0144       rp = jet.rap();
0145       max_r2 = info->max_r2;
0146       switch(info->clust_type) {
0147       case ClusteringVetoPlugin::AKTLIKE: perpfactor = 1./jet.perp2();
0148     break;
0149       case ClusteringVetoPlugin::CALIKE: perpfactor = 1.;
0150     break;
0151       case ClusteringVetoPlugin::KTLIKE: perpfactor = jet.perp2();
0152     break;
0153       default: assert(false);
0154       }
0155     }
0156 
0157     double distance(const ClusteringVetoJet * jet ) const {
0158       double dij = min ( perpfactor, jet->perpfactor );
0159       double dphi = fabs(ph-jet->ph);
0160       if (dphi > pi) {dphi = twopi - dphi;}
0161       dij *= (dphi*dphi+((rp-jet->rp)*(rp-jet->rp))) / max_r2;
0162       return dij;
0163     }
0164 
0165     double beam_distance() const {
0166       return perpfactor;
0167     }
0168 
0169   private:
0170     double ph, rp, perpfactor, max_r2;
0171   };
0172    
0173 } // namespace contrib
0174 
0175 FASTJET_END_NAMESPACE
0176 
0177 #endif  // __FASTJET_CONTRIB_CLUSTERINGVETOPLUGIN_HH__