Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // $Id$
0002 //
0003 // Copyright (c)-, 
0004 //
0005 //----------------------------------------------------------------------
0006 // This file is part of FastJet contrib.
0007 //
0008 // It is free software; you can redistribute it and/or modify it under
0009 // the terms of the GNU General Public License as published by the
0010 // Free Software Foundation; either version 2 of the License, or (at
0011 // your option) any later version.
0012 //
0013 // It is distributed in the hope that it will be useful, but WITHOUT
0014 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
0015 // or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public
0016 // License for more details.
0017 //
0018 // You should have received a copy of the GNU General Public License
0019 // along with this code. If not, see <http://www.gnu.org/licenses/>.
0020 //----------------------------------------------------------------------
0021 
0022 #ifndef __FASTJET_CONTRIB_SUBJETCOUNTING_HH__
0023 #define __FASTJET_CONTRIB_SUBJETCOUNTING_HH__
0024 
0025 #include <fastjet/internal/base.hh>
0026 #include "fastjet/FunctionOfPseudoJet.hh"
0027 #include "fastjet/PseudoJet.hh"
0028 #include "fastjet/ClusterSequence.hh"
0029 #include "fastjet/JetDefinition.hh"
0030 #include <vector>
0031 
0032 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0033 
0034 namespace contrib{
0035 
0036 //------------------------------------------------------------------------
0037 //
0038 //   Implementation of the two subjet counting algorithms used in 
0039 //   "Learning How to Count: A High Multiplicity Search for the LHC"
0040 //   Sonia El Hedri, Anson Hook, Martin Jankowiak, Jay G. Wacker
0041 //   JHEP 1308:136,2013
0042 //   http://arxiv.org/abs/1302.1870
0043 //
0044 //   for usage example see example.cc 
0045 //
0046 //   for questions contact jankowiak@gmail.com
0047 //
0048 //   tested with Fastjet 3.0.3
0049 //
0050 //------------------------------------------------------------------------
0051 
0052 class SubjetCountingKt : public FunctionOfPseudoJet<unsigned int> {
0053 
0054  // **************************************************************************************
0055  //
0056  // constructor for subjet counting with the exclusive Kt algorithm
0057  //
0058  // ~ input:  f_Kt (or rather f_Kt*total_jet_pt) defines the scale at which the exclusive Kt algorithm
0059  //           stops recombination; typical values might be f_Kt ~ 0.04 - 0.10 
0060  // ~ input:  pt_cut is the minimum pt required of the identified subjets included in n_Kt;
0061  //           typical values might be 20-70 GeV
0062  //
0063  // ~ return: a non-negative integer (accessed through the parent class's () operator as 
0064  //           stipulated by the FuntionOfPseudoJet paradigm: see example.cc);
0065  //           alternatively the subjets can be accessed directly via getSubjets
0066  // 
0067  //   see original reference for details and discussion
0068  //
0069  // **************************************************************************************
0070 
0071 public:
0072   SubjetCountingKt(double f_Kt, double pt_cut) : _f_Kt(f_Kt), _pt_cut(pt_cut) {}
0073 
0074   /// default dtor
0075   virtual ~SubjetCountingKt() {}
0076 
0077   /// returns the value of the subjet count for a jet's
0078   /// constituents. (Normally accessed by the parent class's operator()).
0079   unsigned int result(const PseudoJet& jet) const;
0080 
0081   /// get the actual subjets identified by the algorithm
0082   std::vector<fastjet::PseudoJet> getSubjets(const fastjet::PseudoJet &jet) const;
0083 
0084   /// returns the description of the class
0085   std::string description() const; 
0086 
0087 private: 
0088 double _f_Kt, _pt_cut; // parameters defining n_Kt
0089 
0090 // internal function that computes the result
0091 unsigned int _n_Kt(const fastjet::PseudoJet &jet) const;
0092 
0093 }; //end class SubjectCountingKt
0094 
0095 class SubjetCountingCA : public FunctionOfPseudoJet<unsigned int> {
0096  // *************************************************************************************
0097  //
0098  // constructor for subjet counting with the CA algorithm
0099  //
0100  // ~ input:  mass_cut_off is the mass scale down to which the fat jet is declustered;
0101  //           typical values might be 20-70 GeV 
0102  // ~ input:  ycut is the asymmetry parameter that determines how hard to cut on asymmetric
0103  //           splittings; typical values might be ycut ~ 0.10 - 0.15 
0104  // ~ input:  R_min defines a lower angular cutoff to how far the declustering proceeds;
0105  //           this parameter may not make sense or be useful depending on the source of the
0106  //           constituent objects making up the fat jet (topoclusters, etc.);
0107  //           can be made inoperative by being set to zero.
0108  //           typical values might be R_min ~ 0.10 - 0.20 
0109  // ~ input:  pt_cut is the minimum pt required of the identified subjets included in n_CA;
0110  //           typical values might be 20-70 GeV
0111  //
0112  // ~ return: a non-negative integer (accessed through the parent class's () operator as 
0113  //           stipulated by the FuntionOfPseudoJet paradigm: see example.cc);
0114  //           alternatively the subjets can be accessed directly via getSubjets
0115  //
0116  //   see original reference for details and discussion
0117  //
0118  // *************************************************************************************
0119 
0120 public:
0121   SubjetCountingCA(double mass_cut_off, double ycut, double R_min, double pt_cut) : 
0122                   _mass_cut_off(mass_cut_off), _ycut(ycut), _R_min(R_min), _pt_cut(pt_cut) {}
0123 
0124   /// default dtor
0125   virtual ~SubjetCountingCA(){}
0126 
0127   /// returns the value of the subjet count for a jet's
0128   /// constituents. (Normally accessed by the parent class's operator()).
0129   unsigned int result(const PseudoJet& jet) const;
0130 
0131   /// get the actual subjets identified by the algorithm
0132   std::vector<fastjet::PseudoJet> getSubjets(const fastjet::PseudoJet &jet) const;
0133 
0134   /// returns the description of the class
0135   std::string description() const; 
0136 
0137 private: 
0138 double _mass_cut_off, _ycut, _R_min,  _pt_cut; // parameters defining n_CA
0139 
0140 // internal function that calculates the result
0141 unsigned int _n_CA(const fastjet::PseudoJet &jet) const;
0142 
0143 // internal function called recursively by _n_CA
0144 void _FindHardSubst(const fastjet::PseudoJet & this_jet, 
0145                    std::vector<fastjet::PseudoJet> & t_parts) const;
0146 }; // end class SubjectCountingCA
0147 
0148 } // namespace contrib
0149 
0150 FASTJET_END_NAMESPACE
0151 
0152 #endif  // __FASTJET_CONTRIB_SUBJETCOUNTING_HH__